1 #include <symengine/matrix.h>
6 #include <symengine/serialize-cereal.h>
7 #include <symengine/subs.h>
8 #include <symengine/symengine_exception.h>
9 #include <symengine/polys/uexprpoly.h>
11 #include <symengine/test_visitors.h>
17 DenseMatrix::DenseMatrix() {}
19 DenseMatrix::DenseMatrix(
unsigned row,
unsigned col) : row_(row), col_(col)
24 DenseMatrix::DenseMatrix(
unsigned row,
unsigned col,
const vec_basic &l)
25 : m_{l}, row_(row), col_(col){SYMENGINE_ASSERT(m_.size() == row * col)}
27 DenseMatrix::DenseMatrix(
const vec_basic &column_elements)
28 : m_(column_elements), row_(static_cast<unsigned>(column_elements.size())),
34 void DenseMatrix::resize(
unsigned row,
unsigned col)
42 RCP<const Basic> DenseMatrix::get(
unsigned i,
unsigned j)
const
44 SYMENGINE_ASSERT(i < row_ and j < col_);
45 return m_[i * col_ + j];
48 void DenseMatrix::set(
unsigned i,
unsigned j,
const RCP<const Basic> &e)
50 SYMENGINE_ASSERT(i < row_ and j < col_);
55 vec_basic DenseMatrix::as_vec_basic()
const
60 RCP<const Basic> DenseMatrix::trace()
const
62 SYMENGINE_ASSERT(row_ == col_);
65 for (
unsigned i = 0; i < row_; i++) {
66 diag.push_back(m_[offset]);
73 unsigned DenseMatrix::rank()
const
75 throw NotImplementedError(
"Not Implemented");
78 bool DenseMatrix::is_lower()
const
81 unsigned n = A.nrows();
82 for (
unsigned i = 1; i < n; ++i) {
83 for (
unsigned j = 0; j < i; ++j) {
92 bool DenseMatrix::is_upper()
const
95 unsigned n = A.nrows();
96 for (
unsigned i = 0; i < n - 1; ++i) {
97 for (
unsigned j = i + 1; j < n; ++j) {
109 tribool cur = tribool::tritrue;
111 cur = and_tribool(cur, SymEngine::is_zero(*e));
119 tribool DenseMatrix::is_diagonal()
const
122 if (not A.is_square()) {
123 return tribool::trifalse;
125 unsigned ncols = A.ncols();
127 tribool cur = tribool::tritrue;
128 for (
unsigned i = 0; i < ncols; i++) {
130 for (
unsigned j = 0; j < ncols; j++) {
132 auto &e = m_[offset];
133 cur = and_tribool(cur, SymEngine::is_zero(*e));
144 tribool DenseMatrix::is_real(
const Assumptions *assumptions)
const
146 RealVisitor visitor(assumptions);
147 tribool cur = tribool::tritrue;
149 cur = and_tribool(cur, visitor.apply(*e));
157 tribool DenseMatrix::is_symmetric()
const
160 if (not A.is_square()) {
161 return tribool::trifalse;
163 unsigned ncols = A.ncols();
164 tribool cur = tribool::tritrue;
165 for (
unsigned i = 0; i < ncols; i++) {
166 for (
unsigned j = 0; j <= i; j++) {
168 auto &e = m_[i * ncols + j];
169 auto &e2 = m_[j * ncols + i];
170 cur = and_tribool(cur, SymEngine::is_zero(*
sub(e, e2)));
180 tribool DenseMatrix::is_hermitian()
const
183 if (not A.is_square()) {
184 return tribool::trifalse;
186 unsigned ncols = A.ncols();
187 tribool cur = tribool::tritrue;
188 for (
unsigned i = 0; i < ncols; i++) {
189 for (
unsigned j = 0; j <= i; j++) {
190 auto &e = m_[i * ncols + j];
192 auto &e2 = m_[j * ncols + i];
196 cur = and_tribool(cur, SymEngine::is_real(*e));
206 tribool DenseMatrix::is_weakly_diagonally_dominant()
const
209 if (not A.is_square()) {
210 return tribool::trifalse;
213 unsigned ncols = A.ncols();
214 RCP<const Basic> diag;
215 RCP<const Basic> sum;
216 tribool diagdom = tribool::tritrue;
217 for (
unsigned i = 0; i < ncols; i++) {
219 for (
unsigned j = 0; j < ncols; j++) {
220 auto &e = m_[i * ncols + j];
227 diagdom = and_tribool(diagdom, is_nonnegative(*
sub(diag, sum)));
228 if (is_false(diagdom)) {
235 tribool DenseMatrix::is_strictly_diagonally_dominant()
const
238 if (not A.is_square()) {
239 return tribool::trifalse;
242 unsigned ncols = A.ncols();
243 RCP<const Basic> diag;
244 RCP<const Basic> sum;
245 tribool diagdom = tribool::tritrue;
246 for (
unsigned i = 0; i < ncols; i++) {
248 for (
unsigned j = 0; j < ncols; j++) {
249 auto &e = m_[i * ncols + j];
256 diagdom = and_tribool(diagdom, is_positive(*
sub(diag, sum)));
257 if (is_false(diagdom)) {
264 tribool DenseMatrix::shortcut_to_posdef()
const
266 tribool is_diagonal_positive = tribool::tritrue;
268 for (
unsigned i = 0; i < row_; i++) {
270 = and_tribool(is_diagonal_positive, is_positive(*m_[offset]));
271 if (is_false(is_diagonal_positive))
272 return is_diagonal_positive;
275 if (is_true(and_tribool(is_diagonal_positive,
276 this->is_strictly_diagonally_dominant())))
277 return tribool::tritrue;
278 return tribool::indeterminate;
281 tribool DenseMatrix::is_positive_definite_GE()
284 for (
unsigned i = 0; i < size; i++) {
285 auto ispos = is_positive(*m_[i * size + i]);
288 for (
unsigned j = i + 1; j < size; j++) {
289 for (
unsigned k = i + 1; k < size; k++) {
290 m_[j * size + k] =
sub(
mul(m_[i * size + i], m_[j * size + k]),
291 mul(m_[j * size + i], m_[i * size + k]));
295 return tribool::tritrue;
298 tribool DenseMatrix::is_positive_definite()
const
303 const DenseMatrix *H;
304 if (!is_true(A.is_hermitian())) {
306 return tribool::trifalse;
307 DenseMatrix tmp1 = DenseMatrix(A.row_, A.col_);
309 A.conjugate_transpose(tmp1);
310 add_dense_dense(A, tmp1, *B.
get());
316 tribool shortcut = H->shortcut_to_posdef();
317 if (!is_indeterminate(shortcut))
323 return B->is_positive_definite_GE();
326 tribool DenseMatrix::is_negative_definite()
const
328 auto res = DenseMatrix(row_, col_);
329 mul_dense_scalar(*
this,
integer(-1), res);
330 return res.is_positive_definite();
333 RCP<const Basic> DenseMatrix::det()
const
335 return det_bareis(*
this);
338 void DenseMatrix::inv(MatrixBase &result)
const
340 if (is_a<DenseMatrix>(result)) {
341 DenseMatrix &r = down_cast<DenseMatrix &>(result);
342 inverse_pivoted_LU(*
this, r);
346 void DenseMatrix::add_matrix(
const MatrixBase &other, MatrixBase &result)
const
348 SYMENGINE_ASSERT(row_ == result.nrows() and col_ == result.ncols());
350 if (is_a<DenseMatrix>(other) and is_a<DenseMatrix>(result)) {
351 const DenseMatrix &o = down_cast<const DenseMatrix &>(other);
352 DenseMatrix &r = down_cast<DenseMatrix &>(result);
353 add_dense_dense(*
this, o, r);
357 void DenseMatrix::mul_matrix(
const MatrixBase &other, MatrixBase &result)
const
359 SYMENGINE_ASSERT(row_ == result.nrows()
360 and other.ncols() == result.ncols());
362 if (is_a<DenseMatrix>(other) and is_a<DenseMatrix>(result)) {
363 const DenseMatrix &o = down_cast<const DenseMatrix &>(other);
364 DenseMatrix &r = down_cast<DenseMatrix &>(result);
365 mul_dense_dense(*
this, o, r);
369 void DenseMatrix::elementwise_mul_matrix(
const MatrixBase &other,
370 MatrixBase &result)
const
372 SYMENGINE_ASSERT(row_ == result.nrows() and col_ == result.ncols()
373 and row_ == other.nrows() and col_ == other.ncols());
375 if (is_a<DenseMatrix>(other) and is_a<DenseMatrix>(result)) {
376 const DenseMatrix &o = down_cast<const DenseMatrix &>(other);
377 DenseMatrix &r = down_cast<DenseMatrix &>(result);
378 elementwise_mul_dense_dense(*
this, o, r);
383 void DenseMatrix::add_scalar(
const RCP<const Basic> &k,
384 MatrixBase &result)
const
386 if (is_a<DenseMatrix>(result)) {
387 DenseMatrix &r = down_cast<DenseMatrix &>(result);
388 add_dense_scalar(*
this, k, r);
393 void DenseMatrix::mul_scalar(
const RCP<const Basic> &k,
394 MatrixBase &result)
const
396 if (is_a<DenseMatrix>(result)) {
397 DenseMatrix &r = down_cast<DenseMatrix &>(result);
398 mul_dense_scalar(*
this, k, r);
403 void DenseMatrix::conjugate(MatrixBase &result)
const
405 if (is_a<DenseMatrix>(result)) {
406 DenseMatrix &r = down_cast<DenseMatrix &>(result);
407 conjugate_dense(*
this, r);
412 void DenseMatrix::transpose(MatrixBase &result)
const
414 if (is_a<DenseMatrix>(result)) {
415 DenseMatrix &r = down_cast<DenseMatrix &>(result);
416 transpose_dense(*
this, r);
421 void DenseMatrix::conjugate_transpose(MatrixBase &result)
const
423 if (is_a<DenseMatrix>(result)) {
424 DenseMatrix &r = down_cast<DenseMatrix &>(result);
425 conjugate_transpose_dense(*
this, r);
430 void DenseMatrix::submatrix(MatrixBase &result,
unsigned row_start,
431 unsigned col_start,
unsigned row_end,
432 unsigned col_end,
unsigned row_step,
433 unsigned col_step)
const
435 if (is_a<DenseMatrix>(result)) {
436 DenseMatrix &r = down_cast<DenseMatrix &>(result);
437 submatrix_dense(*
this, r, row_start, col_start, row_end, col_end,
443 void DenseMatrix::LU(MatrixBase &L, MatrixBase &U)
const
445 if (is_a<DenseMatrix>(L) and is_a<DenseMatrix>(U)) {
446 DenseMatrix &L_ = down_cast<DenseMatrix &>(L);
447 DenseMatrix &U_ = down_cast<DenseMatrix &>(U);
448 SymEngine::LU(*
this, L_, U_);
453 void DenseMatrix::LDL(MatrixBase &L, MatrixBase &D)
const
455 if (is_a<DenseMatrix>(L) and is_a<DenseMatrix>(D)) {
456 DenseMatrix &L_ = down_cast<DenseMatrix &>(L);
457 DenseMatrix &D_ = down_cast<DenseMatrix &>(D);
458 SymEngine::LDL(*
this, L_, D_);
463 void DenseMatrix::LU_solve(
const MatrixBase &b, MatrixBase &x)
const
465 if (is_a<DenseMatrix>(b) and is_a<DenseMatrix>(x)) {
466 const DenseMatrix &b_ = down_cast<const DenseMatrix &>(b);
467 DenseMatrix &x_ = down_cast<DenseMatrix &>(x);
468 SymEngine::LU_solve(*
this, b_, x_);
473 void DenseMatrix::FFLU(MatrixBase &LU)
const
475 if (is_a<DenseMatrix>(LU)) {
476 DenseMatrix &LU_ = down_cast<DenseMatrix &>(LU);
477 fraction_free_LU(*
this, LU_);
482 void DenseMatrix::FFLDU(MatrixBase &L, MatrixBase &D, MatrixBase &U)
const
484 if (is_a<DenseMatrix>(L) and is_a<DenseMatrix>(D)
485 and is_a<DenseMatrix>(U)) {
486 DenseMatrix &L_ = down_cast<DenseMatrix &>(L);
487 DenseMatrix &D_ = down_cast<DenseMatrix &>(D);
488 DenseMatrix &U_ = down_cast<DenseMatrix &>(U);
489 fraction_free_LDU(*
this, L_, D_, U_);
494 void DenseMatrix::QR(MatrixBase &Q, MatrixBase &R)
const
496 if (is_a<DenseMatrix>(Q) and is_a<DenseMatrix>(R)) {
497 DenseMatrix &Q_ = down_cast<DenseMatrix &>(Q);
498 DenseMatrix &R_ = down_cast<DenseMatrix &>(R);
499 SymEngine::QR(*
this, Q_, R_);
504 void DenseMatrix::cholesky(MatrixBase &L)
const
506 if (is_a<DenseMatrix>(L)) {
507 DenseMatrix &L_ = down_cast<DenseMatrix &>(L);
508 SymEngine::cholesky(*
this, L_);
514 void jacobian(
const DenseMatrix &A,
const DenseMatrix &x, DenseMatrix &result,
517 SYMENGINE_ASSERT(A.col_ == 1);
518 SYMENGINE_ASSERT(x.col_ == 1);
519 SYMENGINE_ASSERT(A.row_ == result.nrows() and x.row_ == result.ncols());
521 #pragma omp parallel for
522 for (
unsigned i = 0; i < result.row_; i++) {
523 for (
unsigned j = 0; j < result.col_; j++) {
524 if (is_a<Symbol>(*(x.m_[j]))) {
525 const RCP<const Symbol> x_
526 = rcp_static_cast<const Symbol>(x.m_[j]);
527 result.m_[i * result.col_ + j] = A.m_[i]->diff(x_, diff_cache);
535 throw SymEngineException(
536 "'x' must contain Symbols only. "
537 "Use sjacobian for SymPy style differentiation");
541 void sjacobian(
const DenseMatrix &A,
const DenseMatrix &x, DenseMatrix &result,
544 SYMENGINE_ASSERT(A.col_ == 1);
545 SYMENGINE_ASSERT(x.col_ == 1);
546 SYMENGINE_ASSERT(A.row_ == result.nrows() and x.row_ == result.ncols());
547 #pragma omp parallel for
548 for (
unsigned i = 0; i < result.row_; i++) {
549 for (
unsigned j = 0; j < result.col_; j++) {
550 if (is_a<Symbol>(*(x.m_[j]))) {
551 const RCP<const Symbol> x_
552 = rcp_static_cast<const Symbol>(x.m_[j]);
553 result.m_[i * result.col_ + j] = A.m_[i]->diff(x_, diff_cache);
556 const RCP<const Symbol> x_ =
symbol(
"x_");
557 result.m_[i * result.col_ + j] =
ssubs(
558 ssubs(A.m_[i], {{x.m_[j], x_}})->diff(x_, diff_cache),
567 void diff(
const DenseMatrix &A,
const RCP<const Symbol> &x, DenseMatrix &result,
570 SYMENGINE_ASSERT(A.row_ == result.nrows() and A.col_ == result.ncols());
571 #pragma omp parallel for
572 for (
unsigned i = 0; i < result.row_; i++) {
573 for (
unsigned j = 0; j < result.col_; j++) {
574 result.m_[i * result.col_ + j]
575 = A.m_[i * result.col_ + j]->diff(x, diff_cache);
580 void sdiff(
const DenseMatrix &A,
const RCP<const Basic> &x, DenseMatrix &result,
583 SYMENGINE_ASSERT(A.row_ == result.nrows() and A.col_ == result.ncols());
584 #pragma omp parallel for
585 for (
unsigned i = 0; i < result.row_; i++) {
586 for (
unsigned j = 0; j < result.col_; j++) {
587 if (is_a<Symbol>(*x)) {
588 const RCP<const Symbol> x_ = rcp_static_cast<const Symbol>(x);
589 result.m_[i * result.col_ + j]
590 = A.m_[i * result.col_ + j]->diff(x_, diff_cache);
593 const RCP<const Symbol> x_ =
symbol(
"_x");
594 result.m_[i * result.col_ + j]
595 =
ssubs(
ssubs(A.m_[i * result.col_ + j], {{x, x_}})
596 ->diff(x_, diff_cache),
604 void conjugate_dense(
const DenseMatrix &A, DenseMatrix &B)
606 SYMENGINE_ASSERT(B.col_ == A.col_ and B.row_ == A.row_);
608 for (
unsigned i = 0; i < A.row_; i++)
609 for (
unsigned j = 0; j < A.col_; j++)
610 B.m_[i * B.col_ + j] =
conjugate(A.m_[i * A.col_ + j]);
614 void transpose_dense(
const DenseMatrix &A, DenseMatrix &B)
616 SYMENGINE_ASSERT(B.row_ == A.col_ and B.col_ == A.row_);
618 for (
unsigned i = 0; i < A.row_; i++)
619 for (
unsigned j = 0; j < A.col_; j++)
620 B.m_[j * B.col_ + i] = A.m_[i * A.col_ + j];
624 void conjugate_transpose_dense(
const DenseMatrix &A, DenseMatrix &B)
626 SYMENGINE_ASSERT(B.row_ == A.col_ and B.col_ == A.row_);
628 for (
unsigned i = 0; i < A.row_; i++)
629 for (
unsigned j = 0; j < A.col_; j++)
630 B.m_[j * B.col_ + i] =
conjugate(A.m_[i * A.col_ + j]);
634 void submatrix_dense(
const DenseMatrix &A, DenseMatrix &B,
unsigned row_start,
635 unsigned col_start,
unsigned row_end,
unsigned col_end,
636 unsigned row_step,
unsigned col_step)
638 SYMENGINE_ASSERT(row_end >= row_start and col_end >= col_start);
639 SYMENGINE_ASSERT(row_end < A.row_);
640 SYMENGINE_ASSERT(col_end < A.col_);
641 SYMENGINE_ASSERT(B.row_ == row_end - row_start + 1
642 and B.col_ == col_end - col_start + 1);
644 unsigned row = B.row_, col = B.col_;
646 for (
unsigned i = 0; i < row; i += row_step)
647 for (
unsigned j = 0; j < col; j += col_step)
648 B.m_[i * col + j] = A.m_[(row_start + i) * A.col_ + col_start + j];
652 void add_dense_dense(
const DenseMatrix &A,
const DenseMatrix &B, DenseMatrix &C)
654 SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_ and A.row_ == C.row_
655 and A.col_ == C.col_);
657 unsigned row = A.row_, col = A.col_;
659 for (
unsigned i = 0; i < row; i++) {
660 for (
unsigned j = 0; j < col; j++) {
661 C.m_[i * col + j] =
add(A.m_[i * col + j], B.m_[i * col + j]);
666 void add_dense_scalar(
const DenseMatrix &A,
const RCP<const Basic> &k,
669 SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_);
671 unsigned row = A.row_, col = A.col_;
673 for (
unsigned i = 0; i < row; i++) {
674 for (
unsigned j = 0; j < col; j++) {
675 B.m_[i * col + j] =
add(A.m_[i * col + j], k);
681 void mul_dense_dense(
const DenseMatrix &A,
const DenseMatrix &B, DenseMatrix &C)
683 SYMENGINE_ASSERT(A.col_ == B.row_ and C.row_ == A.row_
684 and C.col_ == B.col_);
686 unsigned row = A.row_, col = B.col_;
688 if (&A != &C and &B != &C) {
689 for (
unsigned r = 0; r < row; r++) {
690 for (
unsigned c = 0; c < col; c++) {
691 C.m_[r * col + c] = zero;
692 for (
unsigned k = 0; k < A.col_; k++)
694 =
add(C.m_[r * col + c],
695 mul(A.m_[r * A.col_ + k], B.m_[k * col + c]));
699 DenseMatrix tmp = DenseMatrix(A.row_, B.col_);
700 mul_dense_dense(A, B, tmp);
705 void elementwise_mul_dense_dense(
const DenseMatrix &A,
const DenseMatrix &B,
708 SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_ and A.row_ == C.row_
709 and A.col_ == C.col_);
711 unsigned row = A.row_, col = A.col_;
713 for (
unsigned i = 0; i < row; i++) {
714 for (
unsigned j = 0; j < col; j++) {
715 C.m_[i * col + j] =
mul(A.m_[i * col + j], B.m_[i * col + j]);
720 void mul_dense_scalar(
const DenseMatrix &A,
const RCP<const Basic> &k,
723 SYMENGINE_ASSERT(A.col_ == B.col_ and A.row_ == B.row_);
725 unsigned row = A.row_, col = A.col_;
727 for (
unsigned i = 0; i < row; i++) {
728 for (
unsigned j = 0; j < col; j++) {
729 B.m_[i * col + j] =
mul(A.m_[i * col + j], k);
735 void DenseMatrix::row_join(
const DenseMatrix &B)
737 this->col_insert(B, col_);
740 void DenseMatrix::col_join(
const DenseMatrix &B)
742 this->row_insert(B, row_);
745 void DenseMatrix::row_insert(
const DenseMatrix &B,
unsigned pos)
747 SYMENGINE_ASSERT(col_ == B.col_ and pos <= row_)
749 unsigned row = row_, col = col_;
750 this->resize(row_ + B.row_, col_);
752 for (
unsigned i = row; i-- > pos;) {
753 for (
unsigned j = col; j-- > 0;) {
754 this->m_[(i + B.row_) * col + j] = this->m_[i * col + j];
758 for (
unsigned i = 0; i < B.row_; i++) {
759 for (
unsigned j = 0; j < col; j++) {
760 this->m_[(i + pos) * col + j] = B.m_[i * col + j];
765 void DenseMatrix::col_insert(
const DenseMatrix &B,
unsigned pos)
767 SYMENGINE_ASSERT(row_ == B.row_ and pos <= col_)
769 unsigned row = row_, col = col_;
770 this->resize(row_, col_ + B.col_);
772 for (
unsigned i = row; i-- > 0;) {
773 for (
unsigned j = col; j-- > 0;) {
775 this->m_[i * (col + B.col_) + j + B.col_]
776 = this->m_[i * col + j];
778 this->m_[i * (col + B.col_) + j] = this->m_[i * col + j];
783 for (
unsigned i = 0; i < row; i++) {
784 for (
unsigned j = 0; j < B.col_; j++) {
785 this->m_[i * (col + B.col_) + j + pos] = B.m_[i * B.col_ + j];
790 void DenseMatrix::row_del(
unsigned k)
792 SYMENGINE_ASSERT(k < row_)
797 for (
unsigned i = k; i < row_ - 1; i++) {
798 row_exchange_dense(*
this, i, i + 1);
800 this->resize(row_ - 1, col_);
804 void DenseMatrix::col_del(
unsigned k)
806 SYMENGINE_ASSERT(k < col_)
811 unsigned row = row_, col = col_, m = 0;
813 for (
unsigned i = 0; i < row; i++) {
814 for (
unsigned j = 0; j < col; j++) {
816 this->m_[m] = this->m_[i * col + j];
821 this->resize(row_, col_ - 1);
826 void row_exchange_dense(DenseMatrix &A,
unsigned i,
unsigned j)
828 SYMENGINE_ASSERT(i != j and i < A.row_ and j < A.row_);
830 unsigned col = A.col_;
832 for (
unsigned k = 0; k < A.col_; k++)
833 std::swap(A.m_[i * col + k], A.m_[j * col + k]);
836 void row_mul_scalar_dense(DenseMatrix &A,
unsigned i, RCP<const Basic> &c)
838 SYMENGINE_ASSERT(i < A.row_);
840 unsigned col = A.col_;
842 for (
unsigned j = 0; j < A.col_; j++)
843 A.m_[i * col + j] =
mul(c, A.m_[i * col + j]);
846 void row_add_row_dense(DenseMatrix &A,
unsigned i,
unsigned j,
849 SYMENGINE_ASSERT(i != j and i < A.row_ and j < A.row_);
851 unsigned col = A.col_;
853 for (
unsigned k = 0; k < A.col_; k++)
854 A.m_[i * col + k] =
add(A.m_[i * col + k],
mul(c, A.m_[j * col + k]));
857 void permuteFwd(DenseMatrix &A, permutelist &pl)
860 row_exchange_dense(A, p.first, p.second);
865 void column_exchange_dense(DenseMatrix &A,
unsigned i,
unsigned j)
867 SYMENGINE_ASSERT(i != j and i < A.col_ and j < A.col_);
869 unsigned col = A.col_;
871 for (
unsigned k = 0; k < A.row_; k++)
872 std::swap(A.m_[k * col + i], A.m_[k * col + j]);
876 void pivoted_gaussian_elimination(
const DenseMatrix &A, DenseMatrix &B,
879 SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_);
881 unsigned row = A.row_, col = A.col_;
882 unsigned index = 0, i, j, k;
885 RCP<const Basic> scale;
887 for (i = 0; i < col - 1; i++) {
891 k = pivot(B, index, i);
895 row_exchange_dense(B, k, index);
896 pl.push_back({k, index});
899 scale =
div(one, B.m_[index * col + i]);
900 row_mul_scalar_dense(B, index, scale);
902 for (j = i + 1; j < row; j++) {
903 for (k = i + 1; k < col; k++)
905 =
sub(B.m_[j * col + k],
906 mul(B.m_[j * col + i], B.m_[i * col + k]));
907 B.m_[j * col + i] = zero;
917 void fraction_free_gaussian_elimination(
const DenseMatrix &A, DenseMatrix &B)
919 SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_);
921 unsigned col = A.col_;
924 for (
unsigned i = 0; i < col - 1; i++)
925 for (
unsigned j = i + 1; j < A.row_; j++) {
926 for (
unsigned k = i + 1; k < col; k++) {
928 =
sub(
mul(B.m_[i * col + i], B.m_[j * col + k]),
929 mul(B.m_[j * col + i], B.m_[i * col + k]));
932 =
div(B.m_[j * col + k], B.m_[i * col - col + i - 1]);
934 B.m_[j * col + i] = zero;
939 void pivoted_fraction_free_gaussian_elimination(
const DenseMatrix &A,
940 DenseMatrix &B, permutelist &pl)
942 SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_);
944 unsigned col = A.col_, row = A.row_;
945 unsigned index = 0, i, k, j;
948 for (i = 0; i < col - 1; i++) {
952 k = pivot(B, index, i);
956 row_exchange_dense(B, k, index);
957 pl.push_back({k, index});
960 for (j = i + 1; j < row; j++) {
961 for (k = i + 1; k < col; k++) {
963 =
sub(
mul(B.m_[i * col + i], B.m_[j * col + k]),
964 mul(B.m_[j * col + i], B.m_[i * col + k]));
967 =
div(B.m_[j * col + k], B.m_[i * col - col + i - 1]);
969 B.m_[j * col + i] = zero;
976 void pivoted_gauss_jordan_elimination(
const DenseMatrix &A, DenseMatrix &B,
979 SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_);
981 unsigned row = A.row_, col = A.col_;
982 unsigned index = 0, i, j, k;
983 RCP<const Basic> scale;
986 for (i = 0; i < col; i++) {
990 k = pivot(B, index, i);
994 row_exchange_dense(B, k, index);
995 pl.push_back({k, index});
998 scale =
div(one, B.m_[index * col + i]);
999 row_mul_scalar_dense(B, index, scale);
1001 for (j = 0; j < row; j++) {
1005 scale =
mul(minus_one, B.m_[j * col + i]);
1006 row_add_row_dense(B, j, index, scale);
1016 void fraction_free_gauss_jordan_elimination(
const DenseMatrix &A,
1019 SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_);
1021 unsigned row = A.row_, col = A.col_;
1027 for (i = 0; i < col; i++) {
1029 d = B.m_[i * col - col + i - 1];
1030 for (j = 0; j < row; j++)
1032 for (k = 0; k < col; k++) {
1035 =
sub(
mul(B.m_[i * col + i], B.m_[j * col + k]),
1036 mul(B.m_[j * col + i], B.m_[i * col + k]));
1038 B.m_[j * col + k] =
div(B.m_[j * col + k], d);
1041 for (j = 0; j < row; j++)
1043 B.m_[j * col + i] = zero;
1047 void pivoted_fraction_free_gauss_jordan_elimination(
const DenseMatrix &A,
1051 SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_);
1053 unsigned row = A.row_, col = A.col_;
1054 unsigned index = 0, i, k, j;
1059 for (i = 0; i < col; i++) {
1063 k = pivot(B, index, i);
1067 row_exchange_dense(B, k, index);
1068 pl.push_back({k, index});
1071 for (j = 0; j < row; j++) {
1074 for (k = 0; k < col; k++) {
1078 =
sub(
mul(B.m_[index * col + i], B.m_[j * col + k]),
1079 mul(B.m_[j * col + i], B.m_[index * col + k]));
1082 B.m_[j * col + k] =
div(B.m_[j * col + k], d);
1086 d = B.m_[index * col + i];
1088 for (j = 0; j < row; j++) {
1091 B.m_[j * col + i] = zero;
1097 unsigned pivot(DenseMatrix &B,
unsigned r,
unsigned c)
1099 for (
unsigned k = r; k < B.row_; k++) {
1100 if (!is_true(
is_zero(*(B.m_[k * B.col_ + c])))) {
1107 void reduced_row_echelon_form(
const DenseMatrix &A, DenseMatrix &b,
1108 vec_uint &pivot_cols,
bool normalize_last)
1111 if (normalize_last) {
1112 pivoted_fraction_free_gauss_jordan_elimination(A, b, pl);
1114 pivoted_gauss_jordan_elimination(A, b, pl);
1117 for (
unsigned col = 0; col < b.col_ && row < b.row_; col++) {
1118 if (is_true(
is_zero(*b.get(row, col))))
1120 pivot_cols.push_back(col);
1121 if (row == 0 and normalize_last) {
1122 RCP<const Basic> m =
div(one, b.get(row, col));
1131 void diagonal_solve(
const DenseMatrix &A,
const DenseMatrix &b, DenseMatrix &x)
1133 SYMENGINE_ASSERT(A.row_ == A.col_);
1134 SYMENGINE_ASSERT(x.row_ == A.col_ and x.col_ == b.col_);
1136 const unsigned sys = b.col_;
1139 for (
unsigned k = 0; k < sys; k++) {
1140 for (
unsigned i = 0; i < A.col_; i++) {
1141 x.m_[i * sys + k] =
div(b.m_[i * sys + k], A.m_[i * A.col_ + i]);
1146 void back_substitution(
const DenseMatrix &U,
const DenseMatrix &b,
1149 SYMENGINE_ASSERT(U.row_ == U.col_);
1150 SYMENGINE_ASSERT(b.row_ == U.row_);
1151 SYMENGINE_ASSERT(x.row_ == U.col_ and x.col_ == b.col_);
1153 const unsigned col = U.col_;
1154 const unsigned sys = b.col_;
1157 for (
unsigned k = 0; k < sys; k++) {
1158 for (
int i = col - 1; i >= 0; i--) {
1159 for (
unsigned j = i + 1; j < col; j++)
1161 =
sub(x.m_[i * sys + k],
1162 mul(U.m_[i * col + j], x.m_[j * sys + k]));
1163 x.m_[i * sys + k] =
div(x.m_[i * sys + k], U.m_[i * col + i]);
1168 void forward_substitution(
const DenseMatrix &A,
const DenseMatrix &b,
1171 SYMENGINE_ASSERT(A.row_ == A.col_);
1172 SYMENGINE_ASSERT(b.row_ == A.row_);
1173 SYMENGINE_ASSERT(x.row_ == A.col_ and x.col_ == b.col_);
1175 unsigned col = A.col_;
1176 const unsigned sys = b.col_;
1179 for (
unsigned k = 0; k < b.col_; k++) {
1180 for (
unsigned i = 0; i < col - 1; i++) {
1181 for (
unsigned j = i + 1; j < col; j++) {
1183 =
sub(
mul(A.m_[i * col + i], x.m_[j * sys + k]),
1184 mul(A.m_[j * col + i], x.m_[i * sys + k]));
1187 =
div(x.m_[j * sys + k], A.m_[(i - 1) * col + i - 1]);
1193 void fraction_free_gaussian_elimination_solve(
const DenseMatrix &A,
1194 const DenseMatrix &b,
1197 SYMENGINE_ASSERT(A.row_ == A.col_);
1198 SYMENGINE_ASSERT(b.row_ == A.row_ and x.row_ == A.row_);
1199 SYMENGINE_ASSERT(x.col_ == b.col_);
1201 int i, j, k, col = A.col_, bcol = b.col_;
1202 DenseMatrix A_ = DenseMatrix(A.row_, A.col_, A.m_);
1203 DenseMatrix b_ = DenseMatrix(b.row_, b.col_, b.m_);
1205 for (i = 0; i < col - 1; i++)
1206 for (j = i + 1; j < col; j++) {
1207 for (k = 0; k < bcol; k++) {
1209 =
sub(
mul(A_.m_[i * col + i], b_.m_[j * bcol + k]),
1210 mul(A_.m_[j * col + i], b_.m_[i * bcol + k]));
1212 b_.m_[j * bcol + k] =
div(b_.m_[j * bcol + k],
1213 A_.m_[i * col - col + i - 1]);
1216 for (k = i + 1; k < col; k++) {
1218 =
sub(
mul(A_.m_[i * col + i], A_.m_[j * col + k]),
1219 mul(A_.m_[j * col + i], A_.m_[i * col + k]));
1222 =
div(A_.m_[j * col + k], A_.m_[i * col - col + i - 1]);
1224 A_.m_[j * col + i] = zero;
1227 for (i = 0; i < col * bcol; i++)
1230 for (k = 0; k < bcol; k++) {
1231 for (i = col - 1; i >= 0; i--) {
1232 for (j = i + 1; j < col; j++)
1234 =
sub(b_.m_[i * bcol + k],
1235 mul(A_.m_[i * col + j], x.m_[j * bcol + k]));
1236 x.m_[i * bcol + k] =
div(b_.m_[i * bcol + k], A_.m_[i * col + i]);
1241 void fraction_free_gauss_jordan_solve(
const DenseMatrix &A,
1242 const DenseMatrix &b, DenseMatrix &x,
1245 SYMENGINE_ASSERT(A.row_ == A.col_);
1246 SYMENGINE_ASSERT(b.row_ == A.row_ and x.row_ == A.row_);
1247 SYMENGINE_ASSERT(x.col_ == b.col_);
1249 unsigned i, j, k, p, col = A.col_, bcol = b.col_;
1251 DenseMatrix A_ = DenseMatrix(A.row_, A.col_, A.m_);
1252 DenseMatrix b_ = DenseMatrix(b.row_, b.col_, b.m_);
1254 for (i = 0; i < col; i++) {
1256 d = A_.m_[i * col - col + i - 1];
1259 while (p < col and
eq(*A_.m_[p * col + i], *zero)) {
1262 SYMENGINE_ASSERT(p != col);
1265 for (k = i; k < col; k++) {
1266 std::swap(A_.m_[p * col + k], A_.m_[i * col + k]);
1268 for (k = 0; k < bcol; k++) {
1269 std::swap(b_.m_[p * bcol + k], b_.m_[i * bcol + k]);
1273 for (j = 0; j < col; j++)
1275 for (k = 0; k < bcol; k++) {
1277 =
sub(
mul(A_.m_[i * col + i], b_.m_[j * bcol + k]),
1278 mul(A_.m_[j * col + i], b_.m_[i * bcol + k]));
1280 b_.m_[j * bcol + k] =
div(b_.m_[j * bcol + k], d);
1283 for (k = 0; k < col; k++) {
1286 =
sub(
mul(A_.m_[i * col + i], A_.m_[j * col + k]),
1287 mul(A_.m_[j * col + i], A_.m_[i * col + k]));
1289 A_.m_[j * col + k] =
div(A_.m_[j * col + k], d);
1294 for (j = 0; j < col; j++)
1296 A_.m_[j * col + i] = zero;
1300 for (k = 0; k < bcol; k++)
1301 for (i = 0; i < col; i++)
1302 x.m_[i * bcol + k] =
div(b_.m_[i * bcol + k], A_.m_[i * col + i]);
1305 void fraction_free_LU_solve(
const DenseMatrix &A,
const DenseMatrix &b,
1308 DenseMatrix LU = DenseMatrix(A.nrows(), A.ncols());
1309 DenseMatrix x_ = DenseMatrix(b.nrows(), b.ncols());
1311 fraction_free_LU(A, LU);
1312 forward_substitution(LU, b, x_);
1313 back_substitution(LU, x_, x);
1316 void LU_solve(
const DenseMatrix &A,
const DenseMatrix &b, DenseMatrix &x)
1318 DenseMatrix L = DenseMatrix(A.nrows(), A.ncols());
1319 DenseMatrix U = DenseMatrix(A.nrows(), A.ncols());
1320 DenseMatrix x_ = DenseMatrix(b.nrows(), b.ncols());
1323 forward_substitution(L, b, x_);
1324 back_substitution(U, x_, x);
1327 void pivoted_LU_solve(
const DenseMatrix &A,
const DenseMatrix &b,
1330 DenseMatrix L = DenseMatrix(A.nrows(), A.ncols());
1331 DenseMatrix U = DenseMatrix(A.nrows(), A.ncols());
1332 DenseMatrix x_ = DenseMatrix(b);
1335 pivoted_LU(A, L, U, pl);
1337 forward_substitution(L, x_, x_);
1338 back_substitution(U, x_, x);
1341 void LDL_solve(
const DenseMatrix &A,
const DenseMatrix &b, DenseMatrix &x)
1343 DenseMatrix L = DenseMatrix(A.nrows(), A.ncols());
1344 DenseMatrix D = DenseMatrix(A.nrows(), A.ncols());
1345 DenseMatrix x_ = DenseMatrix(b.nrows(), b.ncols());
1347 if (not is_symmetric_dense(A))
1348 throw SymEngineException(
"Matrix must be symmetric");
1351 forward_substitution(L, b, x);
1352 diagonal_solve(D, x, x_);
1353 transpose_dense(L, D);
1354 back_substitution(D, x_, x);
1365 void fraction_free_LU(
const DenseMatrix &A, DenseMatrix &LU)
1367 SYMENGINE_ASSERT(A.row_ == A.col_ and LU.row_ == LU.col_
1368 and A.row_ == LU.row_);
1370 unsigned n = A.row_;
1375 for (i = 0; i < n - 1; i++)
1376 for (j = i + 1; j < n; j++)
1377 for (k = i + 1; k < n; k++) {
1378 LU.m_[j * n + k] =
sub(
mul(LU.m_[i * n + i], LU.m_[j * n + k]),
1379 mul(LU.m_[j * n + i], LU.m_[i * n + k]));
1382 =
div(LU.m_[j * n + k], LU.m_[i * n - n + i - 1]);
1389 void LU(
const DenseMatrix &A, DenseMatrix &L, DenseMatrix &U)
1391 SYMENGINE_ASSERT(A.row_ == A.col_ and L.row_ == L.col_
1392 and U.row_ == U.col_);
1393 SYMENGINE_ASSERT(A.row_ == L.row_ and A.row_ == U.row_);
1395 unsigned n = A.row_;
1397 RCP<const Basic> scale;
1401 for (j = 0; j < n; j++) {
1402 for (i = 0; i < j; i++)
1403 for (k = 0; k < i; k++)
1404 U.m_[i * n + j] =
sub(U.m_[i * n + j],
1405 mul(U.m_[i * n + k], U.m_[k * n + j]));
1407 for (i = j; i < n; i++) {
1408 for (k = 0; k < j; k++)
1409 U.m_[i * n + j] =
sub(U.m_[i * n + j],
1410 mul(U.m_[i * n + k], U.m_[k * n + j]));
1413 scale =
div(one, U.m_[j * n + j]);
1415 for (i = j + 1; i < n; i++)
1416 U.m_[i * n + j] =
mul(U.m_[i * n + j], scale);
1419 for (i = 0; i < n; i++) {
1420 for (j = 0; j < i; j++) {
1421 L.m_[i * n + j] = U.m_[i * n + j];
1422 U.m_[i * n + j] = zero;
1424 L.m_[i * n + i] = one;
1425 for (j = i + 1; j < n; j++)
1426 L.m_[i * n + j] = zero;
1433 void pivoted_LU(
const DenseMatrix &A, DenseMatrix &LU, permutelist &pl)
1435 SYMENGINE_ASSERT(A.row_ == A.col_ and LU.row_ == LU.col_);
1436 SYMENGINE_ASSERT(A.row_ == LU.row_);
1438 unsigned n = A.row_;
1440 RCP<const Basic> scale;
1445 for (j = 0; j < n; j++) {
1446 for (i = 0; i < j; i++)
1447 for (k = 0; k < i; k++)
1448 LU.m_[i * n + j] =
sub(LU.m_[i * n + j],
1449 mul(LU.m_[i * n + k], LU.m_[k * n + j]));
1451 for (i = j; i < n; i++) {
1452 for (k = 0; k < j; k++) {
1453 LU.m_[i * n + j] =
sub(LU.m_[i * n + j],
1454 mul(LU.m_[i * n + k], LU.m_[k * n + j]));
1456 if (pivot == -1 and !is_true(
is_zero(*LU.m_[i * n + j])))
1460 throw SymEngineException(
"Matrix is rank deficient");
1461 if (pivot - j != 0) {
1462 row_exchange_dense(LU, pivot, j);
1463 pl.push_back({pivot, j});
1465 scale =
div(one, LU.m_[j * n + j]);
1467 for (i = j + 1; i < n; i++)
1468 LU.m_[i * n + j] =
mul(LU.m_[i * n + j], scale);
1475 void pivoted_LU(
const DenseMatrix &A, DenseMatrix &L, DenseMatrix &U,
1478 SYMENGINE_ASSERT(A.row_ == A.col_ and L.row_ == L.col_
1479 and U.row_ == U.col_);
1480 SYMENGINE_ASSERT(A.row_ == L.row_ and A.row_ == U.row_);
1482 pivoted_LU(A, U, pl);
1483 unsigned n = A.col_;
1484 for (
unsigned i = 0; i < n; i++) {
1485 for (
unsigned j = 0; j < i; j++) {
1486 L.m_[i * n + j] = U.m_[i * n + j];
1487 U.m_[i * n + j] = zero;
1489 L.m_[i * n + i] = one;
1490 for (
unsigned j = i + 1; j < n; j++)
1491 L.m_[i * n + j] = zero;
1500 void fraction_free_LDU(
const DenseMatrix &A, DenseMatrix &L, DenseMatrix &D,
1503 SYMENGINE_ASSERT(A.row_ == L.row_ and A.row_ == U.row_);
1504 SYMENGINE_ASSERT(A.col_ == L.col_ and A.col_ == U.col_);
1506 unsigned row = A.row_, col = A.col_;
1508 RCP<const Basic> old =
integer(1);
1513 for (i = 0; i < row; i++)
1514 for (j = 0; j < row; j++)
1516 L.m_[i * col + j] = zero;
1518 L.m_[i * col + i] = one;
1521 for (i = 0; i < row * row; i++)
1524 for (k = 0; k < row - 1; k++) {
1525 L.m_[k * col + k] = U.m_[k * col + k];
1526 D.m_[k * col + k] =
mul(old, U.m_[k * col + k]);
1528 for (i = k + 1; i < row; i++) {
1529 L.m_[i * col + k] = U.m_[i * col + k];
1530 for (j = k + 1; j < col; j++)
1532 =
div(
sub(
mul(U.m_[k * col + k], U.m_[i * col + j]),
1533 mul(U.m_[k * col + j], U.m_[i * col + k])),
1535 U.m_[i * col + k] = zero;
1538 old = U.m_[k * col + k];
1541 D.m_[row * col - col + row - 1] = old;
1546 void QR(
const DenseMatrix &A, DenseMatrix &Q, DenseMatrix &R)
1548 unsigned row = A.row_;
1549 unsigned col = A.col_;
1551 SYMENGINE_ASSERT(Q.row_ == row and Q.col_ == col and R.row_ == col
1559 for (i = 0; i < row * col; i++)
1563 for (i = 0; i < col * col; i++)
1566 for (j = 0; j < col; j++) {
1568 for (k = 0; k < row; k++)
1569 tmp[k] = A.m_[k * col + j];
1571 for (i = 0; i < j; i++) {
1573 for (k = 0; k < row; k++)
1574 t =
add(t,
mul(A.m_[k * col + j], Q.m_[k * col + i]));
1575 for (k = 0; k < row; k++)
1576 tmp[k] =
expand(
sub(tmp[k],
mul(Q.m_[k * col + i], t)));
1581 for (k = 0; k < row; k++)
1586 R.m_[j * col + j] = t;
1587 for (k = 0; k < row; k++)
1588 Q.m_[k * col + j] =
div(tmp[k], t);
1590 for (i = 0; i < j; i++) {
1592 for (k = 0; k < row; k++)
1593 t =
add(t,
mul(Q.m_[k * col + i], A.m_[k * col + j]));
1594 R.m_[i * col + j] = t;
1601 void LDL(
const DenseMatrix &A, DenseMatrix &L, DenseMatrix &D)
1603 SYMENGINE_ASSERT(A.row_ == A.col_);
1604 SYMENGINE_ASSERT(L.row_ == A.row_ and L.col_ == A.row_);
1605 SYMENGINE_ASSERT(D.row_ == A.row_ and D.col_ == A.row_);
1607 unsigned col = A.col_;
1609 RCP<const Basic> sum;
1610 RCP<const Basic> i2 =
integer(2);
1613 for (i = 0; i < col; i++)
1614 for (j = 0; j < col; j++)
1615 D.m_[i * col + j] = zero;
1618 for (i = 0; i < col; i++)
1619 for (j = 0; j < col; j++)
1620 L.m_[i * col + j] = (i != j) ? zero : one;
1622 for (i = 0; i < col; i++) {
1623 for (j = 0; j < i; j++) {
1625 for (k = 0; k < j; k++)
1626 sum =
add(sum,
mul(
mul(L.m_[i * col + k], L.m_[j * col + k]),
1627 D.m_[k * col + k]));
1629 =
mul(
div(one, D.m_[j * col + j]),
sub(A.m_[i * col + j], sum));
1632 for (k = 0; k < i; k++)
1633 sum =
add(sum,
mul(
pow(L.m_[i * col + k], i2), D.m_[k * col + k]));
1634 D.m_[i * col + i] =
sub(A.m_[i * col + i], sum);
1639 void cholesky(
const DenseMatrix &A, DenseMatrix &L)
1641 SYMENGINE_ASSERT(A.row_ == A.col_);
1642 SYMENGINE_ASSERT(L.row_ == A.row_ and L.col_ == A.row_);
1644 unsigned col = A.col_;
1646 RCP<const Basic> sum;
1647 RCP<const Basic> i2 =
integer(2);
1648 RCP<const Basic> half =
div(one, i2);
1651 for (i = 0; i < col; i++)
1652 for (j = 0; j < col; j++)
1653 L.m_[i * col + j] = zero;
1655 for (i = 0; i < col; i++) {
1656 for (j = 0; j < i; j++) {
1658 for (k = 0; k < j; k++)
1659 sum =
add(sum,
mul(L.m_[i * col + k], L.m_[j * col + k]));
1662 =
mul(
div(one, L.m_[j * col + j]),
sub(A.m_[i * col + j], sum));
1665 for (k = 0; k < i; k++)
1666 sum =
add(sum,
pow(L.m_[i * col + k], i2));
1667 L.m_[i * col + i] =
pow(
sub(A.m_[i * col + i], sum), half);
1672 bool is_symmetric_dense(
const DenseMatrix &A)
1674 if (A.col_ != A.row_)
1677 unsigned col = A.col_;
1680 for (
unsigned i = 0; i < col; i++)
1681 for (
unsigned j = i + 1; j < col; j++)
1682 if (not
eq(*(A.m_[j * col + i]), *(A.m_[i * col + j]))) {
1691 RCP<const Basic> det_bareis(
const DenseMatrix &A)
1693 SYMENGINE_ASSERT(A.row_ == A.col_);
1695 unsigned n = A.row_;
1699 }
else if (n == 2) {
1701 return sub(
mul(A.m_[0], A.m_[3]),
mul(A.m_[1], A.m_[2]));
1702 }
else if (n == 3) {
1706 mul(
mul(A.m_[1], A.m_[5]), A.m_[6])),
1707 mul(
mul(A.m_[2], A.m_[3]), A.m_[7])),
1709 mul(
mul(A.m_[1], A.m_[3]), A.m_[8])),
1710 mul(
mul(A.m_[0], A.m_[5]), A.m_[7])));
1713 if (A.is_lower() or A.is_upper()) {
1714 RCP<const Basic> det = A.m_[0];
1715 for (
unsigned i = 1; i < n; ++i) {
1716 det =
mul(det, A.m_[i * n + i]);
1721 DenseMatrix B = DenseMatrix(n, n, A.m_);
1722 unsigned i,
sign = 1;
1725 for (
unsigned k = 0; k < n - 1; k++) {
1726 if (is_true(
is_zero(*B.m_[k * n + k]))) {
1727 for (i = k + 1; i < n; i++)
1728 if (!is_true(
is_zero(*B.m_[i * n + k]))) {
1729 row_exchange_dense(B, i, k);
1737 for (i = k + 1; i < n; i++) {
1738 for (
unsigned j = k + 1; j < n; j++) {
1739 d =
sub(
mul(B.m_[k * n + k], B.m_[i * n + j]),
1740 mul(B.m_[i * n + k], B.m_[k * n + j]));
1742 d =
div(d, B.m_[(k - 1) * n + k - 1]);
1743 B.m_[i * n + j] = d;
1748 return (sign == 1) ? B.m_[n * n - 1] :
mul(minus_one, B.m_[n * n - 1]);
1758 SYMENGINE_ASSERT(A.row_ == A.col_);
1760 unsigned col = A.col_;
1761 unsigned i, k, l, m;
1767 for (
unsigned n = col; n > 1; n--) {
1770 DenseMatrix T = DenseMatrix(n + 1, n);
1771 DenseMatrix C = DenseMatrix(k, 1);
1774 for (i = 0; i < n * (n + 1); i++)
1776 for (i = 0; i < k; i++)
1777 C.m_[i] = A.m_[i * col + k];
1780 for (i = 0; i < n - 2; i++) {
1781 DenseMatrix B = DenseMatrix(k, 1);
1782 for (l = 0; l < k; l++) {
1784 for (m = 0; m < k; m++)
1786 =
add(B.m_[l],
mul(A.m_[l * col + m], items[i].m_[m]));
1792 for (i = 0; i < n - 1; i++) {
1793 RCP<const Basic> element = zero;
1794 for (l = 0; l < k; l++)
1795 element =
add(element,
mul(A.m_[k * col + l], items[i].m_[l]));
1798 items_.
insert(items_.
begin(),
mul(minus_one, A.m_[k * col + k]));
1801 for (i = 0; i < n; i++) {
1802 for (l = 0; l < n - i + 1; l++)
1803 T.m_[(i + l) * n + i] = items_[l];
1809 polys.
push_back(DenseMatrix(2, 1, {one,
mul(A.m_[0], minus_one)}));
1811 for (i = 0; i < col - 1; i++) {
1812 unsigned t_row = transforms[col - 2 - i].nrows();
1813 unsigned t_col = transforms[col - 2 - i].ncols();
1814 DenseMatrix B = DenseMatrix(t_row, 1);
1816 for (l = 0; l < t_row; l++) {
1818 for (m = 0; m < t_col; m++) {
1819 B.m_[l] =
add(B.m_[l],
1820 mul(transforms[col - 2 - i].m_[l * t_col + m],
1822 B.m_[l] =
expand(B.m_[l]);
1829 RCP<const Basic> det_berkowitz(
const DenseMatrix &A)
1833 berkowitz(A, polys);
1834 DenseMatrix poly = polys[polys.
size() - 1];
1836 if (polys.
size() % 2 == 1)
1837 return mul(minus_one, poly.get(poly.nrows() - 1, 0));
1839 return poly.get(poly.nrows() - 1, 0);
1842 void char_poly(
const DenseMatrix &A, DenseMatrix &B)
1844 SYMENGINE_ASSERT(B.ncols() == 1 and B.nrows() == A.nrows() + 1);
1845 SYMENGINE_ASSERT(A.nrows() == A.ncols());
1849 berkowitz(A, polys);
1850 B = polys[polys.
size() - 1];
1853 void inverse_fraction_free_LU(
const DenseMatrix &A, DenseMatrix &B)
1855 SYMENGINE_ASSERT(A.row_ == A.col_ and B.row_ == B.col_
1856 and B.row_ == A.row_);
1858 unsigned n = A.row_, i;
1859 DenseMatrix LU = DenseMatrix(n, n);
1860 DenseMatrix e = DenseMatrix(n, 1);
1861 DenseMatrix x = DenseMatrix(n, 1);
1862 DenseMatrix x_ = DenseMatrix(n, 1);
1865 for (i = 0; i < n * n; i++) {
1869 for (i = 0; i < n; i++) {
1875 fraction_free_LU(A, LU);
1880 for (
unsigned j = 0; j < n; j++) {
1883 forward_substitution(LU, e, x_);
1884 back_substitution(LU, x_, x);
1886 for (i = 0; i < n; i++)
1887 B.m_[i * n + j] = x.m_[i];
1893 void inverse_LU(
const DenseMatrix &A, DenseMatrix &B)
1895 SYMENGINE_ASSERT(A.row_ == A.col_ and B.row_ == B.col_
1896 and B.row_ == A.row_);
1897 DenseMatrix e = DenseMatrix(A.row_, A.col_);
1902 void inverse_pivoted_LU(
const DenseMatrix &A, DenseMatrix &B)
1904 SYMENGINE_ASSERT(A.row_ == A.col_ and B.row_ == B.col_
1905 and B.row_ == A.row_);
1906 DenseMatrix e = DenseMatrix(A.row_, A.col_);
1908 pivoted_LU_solve(A, e, B);
1911 void inverse_gauss_jordan(
const DenseMatrix &A, DenseMatrix &B)
1913 SYMENGINE_ASSERT(A.row_ == A.col_ and B.row_ == B.col_
1914 and B.row_ == A.row_);
1916 unsigned n = A.row_;
1917 DenseMatrix e = DenseMatrix(n, n);
1920 for (
unsigned i = 0; i < n; i++)
1921 for (
unsigned j = 0; j < n; j++) {
1923 e.m_[i * n + j] = zero;
1925 e.m_[i * n + i] = one;
1927 B.m_[i * n + j] = zero;
1930 fraction_free_gauss_jordan_solve(A, e, B);
1935 void dot(
const DenseMatrix &A,
const DenseMatrix &B, DenseMatrix &C)
1937 if (A.col_ == B.row_) {
1939 DenseMatrix tmp1 = DenseMatrix(A.col_, A.row_);
1941 DenseMatrix tmp2 = DenseMatrix(B.col_, B.row_);
1943 C.resize(tmp1.row_, tmp2.col_);
1944 mul_dense_dense(tmp1, tmp2, C);
1946 C.resize(A.row_, B.col_);
1947 mul_dense_dense(A, B, C);
1949 C.resize(1, C.row_ * C.col_);
1950 }
else if (A.col_ == B.col_) {
1951 DenseMatrix tmp2 = DenseMatrix(B.col_, B.row_);
1954 }
else if (A.row_ == B.row_) {
1955 DenseMatrix tmp1 = DenseMatrix(A.col_, A.row_);
1959 throw SymEngineException(
"Dimensions incorrect for dot product");
1963 void cross(
const DenseMatrix &A,
const DenseMatrix &B, DenseMatrix &C)
1965 SYMENGINE_ASSERT((A.row_ * A.col_ == 3 and B.row_ * B.col_ == 3)
1966 and (A.row_ == C.row_ and A.col_ == C.col_));
1968 C.m_[0] =
sub(
mul(A.m_[1], B.m_[2]),
mul(A.m_[2], B.m_[1]));
1969 C.m_[1] =
sub(
mul(A.m_[2], B.m_[0]),
mul(A.m_[0], B.m_[2]));
1970 C.m_[2] =
sub(
mul(A.m_[0], B.m_[1]),
mul(A.m_[1], B.m_[0]));
1973 RCP<const Set> eigen_values(
const DenseMatrix &A)
1975 unsigned n = A.nrows();
1976 if (A.is_lower() or A.is_upper()) {
1977 RCP<const Set> eigenvals =
emptyset();
1979 for (
unsigned i = 0; i < n; ++i) {
1986 DenseMatrix B = DenseMatrix(A.nrows() + 1, 1);
1988 map_int_Expr coeffs;
1989 auto nr = A.nrows();
1990 for (
unsigned i = 0; i <= nr; i++)
1991 insert(coeffs, nr - i, B.get(i, 0));
1992 auto lambda =
symbol(
"lambda");
1993 return solve_poly(uexpr_poly(lambda, coeffs), lambda);
1999 void eye(DenseMatrix &A,
int k)
2001 if ((k >= 0 and (
unsigned) k >= A.col_) or k + A.row_ <= 0) {
2005 vec_basic v = vec_basic(k > 0 ? A.col_ - k : A.row_ + k, one);
2011 void diag(DenseMatrix &A, vec_basic &v,
int k)
2013 SYMENGINE_ASSERT(v.size() > 0);
2015 unsigned k_ = std::abs(k);
2018 for (
unsigned i = 0; i < A.row_; i++) {
2019 for (
unsigned j = 0; j < A.col_; j++) {
2020 if (j != (
unsigned)k) {
2021 A.m_[i * A.col_ + j] = zero;
2023 A.m_[i * A.col_ + j] = v[k - k_];
2031 for (
unsigned j = 0; j < A.col_; j++) {
2032 for (
unsigned i = 0; i < A.row_; i++) {
2033 if (i != (
unsigned)k) {
2034 A.m_[i * A.col_ + j] = zero;
2036 A.m_[i * A.col_ + j] = v[k - k_];
2045 void ones(DenseMatrix &A)
2047 for (
unsigned i = 0; i < A.row_ * A.col_; i++) {
2053 void zeros(DenseMatrix &A)
2055 for (
unsigned i = 0; i < A.row_ * A.col_; i++) {
2063 unsigned short major = SYMENGINE_MAJOR_VERSION;
2064 unsigned short minor = SYMENGINE_MINOR_VERSION;
2065 cereal::PortableBinaryOutputArchive{oss}(major, minor, row_, col_, m_);
2071 unsigned short major, minor;
2075 cereal::PortableBinaryInputArchive iarchive{iss};
2076 iarchive(major, minor);
2077 if (major != SYMENGINE_MAJOR_VERSION or minor != SYMENGINE_MINOR_VERSION) {
2078 throw SerializationError(StreamFmt()
2079 <<
"SymEngine-" << SYMENGINE_MAJOR_VERSION
2080 <<
"." << SYMENGINE_MINOR_VERSION
2081 <<
" was asked to deserialize an object "
2082 <<
"created using SymEngine-" << major <<
"."
2085 iarchive(row, col, obj);
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.