dense_matrix.cpp
1 #include <symengine/matrix.h>
2 #include <symengine/number.h>
3 #include <symengine/add.h>
4 #include <symengine/functions.h>
5 #include <symengine/pow.h>
6 #include <symengine/subs.h>
7 #include <symengine/symengine_exception.h>
8 #include <symengine/polys/uexprpoly.h>
9 #include <symengine/solve.h>
10 #include <symengine/test_visitors.h>
11 #if HAVE_SYMENGINE_RTTI
12 #include <symengine/serialize-cereal.h>
13 #endif
14 
15 namespace SymEngine
16 {
17 
18 // Constructors
19 DenseMatrix::DenseMatrix() {}
20 
21 DenseMatrix::DenseMatrix(unsigned row, unsigned col) : row_(row), col_(col)
22 {
23  m_ = std::vector<RCP<const Basic>>(row * col);
24 }
25 
26 DenseMatrix::DenseMatrix(unsigned row, unsigned col, const vec_basic &l)
27  : m_{l}, row_(row), col_(col){SYMENGINE_ASSERT(m_.size() == row * col)}
28 
29  DenseMatrix::DenseMatrix(const vec_basic &column_elements)
30  : m_(column_elements), row_(static_cast<unsigned>(column_elements.size())),
31  col_(1)
32 {
33 }
34 
35 // Resize the Matrix
36 void DenseMatrix::resize(unsigned row, unsigned col)
37 {
38  row_ = row;
39  col_ = col;
40  m_.resize(row * col);
41 }
42 
43 // Get and set elements
44 RCP<const Basic> DenseMatrix::get(unsigned i, unsigned j) const
45 {
46  SYMENGINE_ASSERT(i < row_ and j < col_);
47  return m_[i * col_ + j];
48 }
49 
50 void DenseMatrix::set(unsigned i, unsigned j, const RCP<const Basic> &e)
51 {
52  SYMENGINE_ASSERT(i < row_ and j < col_);
53  m_[i * col_ + j] = e;
54 }
55 
56 // Flat (row major) view of Matrix
57 vec_basic DenseMatrix::as_vec_basic() const
58 {
59  return m_;
60 }
61 
62 RCP<const Basic> DenseMatrix::trace() const
63 {
64  SYMENGINE_ASSERT(row_ == col_);
65  unsigned offset = 0;
66  vec_basic diag;
67  for (unsigned i = 0; i < row_; i++) {
68  diag.push_back(m_[offset]);
69  offset += row_ + 1;
70  }
71  auto sum = add(diag);
72  return sum;
73 }
74 
75 unsigned DenseMatrix::rank() const
76 {
77  throw NotImplementedError("Not Implemented");
78 }
79 
80 bool DenseMatrix::is_lower() const
81 {
82  auto A = *this;
83  unsigned n = A.nrows();
84  for (unsigned i = 1; i < n; ++i) {
85  for (unsigned j = 0; j < i; ++j) {
86  if (not is_number_and_zero(*A.get(i, j))) {
87  return false;
88  }
89  }
90  }
91  return true;
92 }
93 
94 bool DenseMatrix::is_upper() const
95 {
96  auto A = *this;
97  unsigned n = A.nrows();
98  for (unsigned i = 0; i < n - 1; ++i) {
99  for (unsigned j = i + 1; j < n; ++j) {
100  if (not is_number_and_zero(*A.get(i, j))) {
101  return false;
102  }
103  }
104  }
105  return true;
106 }
107 
108 tribool DenseMatrix::is_zero() const
109 {
110  auto A = *this;
111  tribool cur = tribool::tritrue;
112  for (auto &e : m_) {
113  cur = and_tribool(cur, SymEngine::is_zero(*e));
114  if (is_false(cur)) {
115  return cur;
116  }
117  }
118  return cur;
119 }
120 
121 tribool DenseMatrix::is_diagonal() const
122 {
123  auto A = *this;
124  if (not A.is_square()) {
125  return tribool::trifalse;
126  }
127  unsigned ncols = A.ncols();
128  unsigned offset;
129  tribool cur = tribool::tritrue;
130  for (unsigned i = 0; i < ncols; i++) {
131  offset = i * ncols;
132  for (unsigned j = 0; j < ncols; j++) {
133  if (j != i) {
134  auto &e = m_[offset];
135  cur = and_tribool(cur, SymEngine::is_zero(*e));
136  if (is_false(cur)) {
137  return cur;
138  }
139  }
140  offset++;
141  }
142  }
143  return cur;
144 }
145 
146 tribool DenseMatrix::is_real(const Assumptions *assumptions) const
147 {
148  RealVisitor visitor(assumptions);
149  tribool cur = tribool::tritrue;
150  for (auto &e : m_) {
151  cur = and_tribool(cur, visitor.apply(*e));
152  if (is_false(cur)) {
153  return cur;
154  }
155  }
156  return cur;
157 }
158 
159 tribool DenseMatrix::is_symmetric() const
160 {
161  auto A = *this;
162  if (not A.is_square()) {
163  return tribool::trifalse;
164  }
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++) {
169  if (j != i) {
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)));
173  }
174  if (is_false(cur)) {
175  return cur;
176  }
177  }
178  }
179  return cur;
180 }
181 
182 tribool DenseMatrix::is_hermitian() const
183 {
184  auto A = *this;
185  if (not A.is_square()) {
186  return tribool::trifalse;
187  }
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];
193  if (j != i) {
194  auto &e2 = m_[j * ncols + i];
195  cur = and_tribool(
196  cur, SymEngine::is_zero(*sub(e, SymEngine::conjugate(e2))));
197  } else {
198  cur = and_tribool(cur, SymEngine::is_real(*e));
199  }
200  if (is_false(cur)) {
201  return cur;
202  }
203  }
204  }
205  return cur;
206 }
207 
208 tribool DenseMatrix::is_weakly_diagonally_dominant() const
209 {
210  auto A = *this;
211  if (not A.is_square()) {
212  return tribool::trifalse;
213  }
214 
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++) {
220  sum = zero;
221  for (unsigned j = 0; j < ncols; j++) {
222  auto &e = m_[i * ncols + j];
223  if (i == j) {
224  diag = abs(e);
225  } else {
226  sum = add(sum, abs(e));
227  }
228  }
229  diagdom = and_tribool(diagdom, is_nonnegative(*sub(diag, sum)));
230  if (is_false(diagdom)) {
231  return diagdom;
232  }
233  }
234  return diagdom;
235 }
236 
237 tribool DenseMatrix::is_strictly_diagonally_dominant() const
238 {
239  auto A = *this;
240  if (not A.is_square()) {
241  return tribool::trifalse;
242  }
243 
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++) {
249  sum = zero;
250  for (unsigned j = 0; j < ncols; j++) {
251  auto &e = m_[i * ncols + j];
252  if (i == j) {
253  diag = abs(e);
254  } else {
255  sum = add(sum, abs(e));
256  }
257  }
258  diagdom = and_tribool(diagdom, is_positive(*sub(diag, sum)));
259  if (is_false(diagdom)) {
260  return diagdom;
261  }
262  }
263  return diagdom;
264 }
265 
266 tribool DenseMatrix::shortcut_to_posdef() const
267 {
268  tribool is_diagonal_positive = tribool::tritrue;
269  unsigned offset = 0;
270  for (unsigned i = 0; i < row_; i++) {
271  is_diagonal_positive
272  = and_tribool(is_diagonal_positive, is_positive(*m_[offset]));
273  if (is_false(is_diagonal_positive))
274  return is_diagonal_positive;
275  offset += row_ + 1;
276  }
277  if (is_true(and_tribool(is_diagonal_positive,
278  this->is_strictly_diagonally_dominant())))
279  return tribool::tritrue;
280  return tribool::indeterminate;
281 }
282 
283 tribool DenseMatrix::is_positive_definite_GE()
284 {
285  auto size = row_;
286  for (unsigned i = 0; i < size; i++) {
287  auto ispos = is_positive(*m_[i * size + i]);
288  if (!is_true(ispos))
289  return ispos;
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]));
294  }
295  }
296  }
297  return tribool::tritrue;
298 }
299 
300 tribool DenseMatrix::is_positive_definite() const
301 {
302 
303  auto A = *this;
305  const DenseMatrix *H;
306  if (!is_true(A.is_hermitian())) {
307  if (!A.is_square())
308  return tribool::trifalse;
309  DenseMatrix tmp1 = DenseMatrix(A.row_, A.col_);
310  B = std::unique_ptr<DenseMatrix>(new DenseMatrix(A.row_, A.col_));
311  A.conjugate_transpose(tmp1);
312  add_dense_dense(A, tmp1, *B.get());
313  H = B.get();
314  } else {
315  H = this;
316  }
317 
318  tribool shortcut = H->shortcut_to_posdef();
319  if (!is_indeterminate(shortcut))
320  return shortcut;
321 
322  if (!B) {
323  B = std::unique_ptr<DenseMatrix>(new DenseMatrix(A));
324  }
325  return B->is_positive_definite_GE();
326 }
327 
328 tribool DenseMatrix::is_negative_definite() const
329 {
330  auto res = DenseMatrix(row_, col_);
331  mul_dense_scalar(*this, integer(-1), res);
332  return res.is_positive_definite();
333 }
334 
335 RCP<const Basic> DenseMatrix::det() const
336 {
337  return det_bareis(*this);
338 }
339 
340 void DenseMatrix::inv(MatrixBase &result) const
341 {
342  if (is_a<DenseMatrix>(result)) {
343  DenseMatrix &r = down_cast<DenseMatrix &>(result);
344  inverse_pivoted_LU(*this, r);
345  }
346 }
347 
348 void DenseMatrix::add_matrix(const MatrixBase &other, MatrixBase &result) const
349 {
350  SYMENGINE_ASSERT(row_ == result.nrows() and col_ == result.ncols());
351 
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);
356  }
357 }
358 
359 void DenseMatrix::mul_matrix(const MatrixBase &other, MatrixBase &result) const
360 {
361  SYMENGINE_ASSERT(row_ == result.nrows()
362  and other.ncols() == result.ncols());
363 
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);
368  }
369 }
370 
371 void DenseMatrix::elementwise_mul_matrix(const MatrixBase &other,
372  MatrixBase &result) const
373 {
374  SYMENGINE_ASSERT(row_ == result.nrows() and col_ == result.ncols()
375  and row_ == other.nrows() and col_ == other.ncols());
376 
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);
381  }
382 }
383 
384 // Add a scalar
385 void DenseMatrix::add_scalar(const RCP<const Basic> &k,
386  MatrixBase &result) const
387 {
388  if (is_a<DenseMatrix>(result)) {
389  DenseMatrix &r = down_cast<DenseMatrix &>(result);
390  add_dense_scalar(*this, k, r);
391  }
392 }
393 
394 // Multiply by a scalar
395 void DenseMatrix::mul_scalar(const RCP<const Basic> &k,
396  MatrixBase &result) const
397 {
398  if (is_a<DenseMatrix>(result)) {
399  DenseMatrix &r = down_cast<DenseMatrix &>(result);
400  mul_dense_scalar(*this, k, r);
401  }
402 }
403 
404 // Matrix conjugate
405 void DenseMatrix::conjugate(MatrixBase &result) const
406 {
407  if (is_a<DenseMatrix>(result)) {
408  DenseMatrix &r = down_cast<DenseMatrix &>(result);
409  conjugate_dense(*this, r);
410  }
411 }
412 
413 // Matrix transpose
414 void DenseMatrix::transpose(MatrixBase &result) const
415 {
416  if (is_a<DenseMatrix>(result)) {
417  DenseMatrix &r = down_cast<DenseMatrix &>(result);
418  transpose_dense(*this, r);
419  }
420 }
421 
422 // Matrix conjugate transpose
423 void DenseMatrix::conjugate_transpose(MatrixBase &result) const
424 {
425  if (is_a<DenseMatrix>(result)) {
426  DenseMatrix &r = down_cast<DenseMatrix &>(result);
427  conjugate_transpose_dense(*this, r);
428  }
429 }
430 
431 // Extract out a submatrix
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
436 {
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,
440  row_step, col_step);
441  }
442 }
443 
444 // LU factorization
445 void DenseMatrix::LU(MatrixBase &L, MatrixBase &U) const
446 {
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_);
451  }
452 }
453 
454 // LDL factorization
455 void DenseMatrix::LDL(MatrixBase &L, MatrixBase &D) const
456 {
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_);
461  }
462 }
463 
464 // Solve Ax = b using LU factorization
465 void DenseMatrix::LU_solve(const MatrixBase &b, MatrixBase &x) const
466 {
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_);
471  }
472 }
473 
474 // Fraction free LU factorization
475 void DenseMatrix::FFLU(MatrixBase &LU) const
476 {
477  if (is_a<DenseMatrix>(LU)) {
478  DenseMatrix &LU_ = down_cast<DenseMatrix &>(LU);
479  fraction_free_LU(*this, LU_);
480  }
481 }
482 
483 // Fraction free LDU factorization
484 void DenseMatrix::FFLDU(MatrixBase &L, MatrixBase &D, MatrixBase &U) const
485 {
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_);
492  }
493 }
494 
495 // QR factorization
496 void DenseMatrix::QR(MatrixBase &Q, MatrixBase &R) const
497 {
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_);
502  }
503 }
504 
505 // Cholesky decomposition
506 void DenseMatrix::cholesky(MatrixBase &L) const
507 {
508  if (is_a<DenseMatrix>(L)) {
509  DenseMatrix &L_ = down_cast<DenseMatrix &>(L);
510  SymEngine::cholesky(*this, L_);
511  }
512 }
513 
514 // ---------------------------- Jacobian -------------------------------------//
515 
516 void jacobian(const DenseMatrix &A, const DenseMatrix &x, DenseMatrix &result,
517  bool diff_cache)
518 {
519  SYMENGINE_ASSERT(A.col_ == 1);
520  SYMENGINE_ASSERT(x.col_ == 1);
521  SYMENGINE_ASSERT(A.row_ == result.nrows() and x.row_ == result.ncols());
522  bool error = false;
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);
530  } else {
531  error = true;
532  break;
533  }
534  }
535  }
536  if (error) {
537  throw SymEngineException(
538  "'x' must contain Symbols only. "
539  "Use sjacobian for SymPy style differentiation");
540  }
541 }
542 
543 void sjacobian(const DenseMatrix &A, const DenseMatrix &x, DenseMatrix &result,
544  bool diff_cache)
545 {
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);
556  } else {
557  // TODO: Use a dummy symbol
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),
561  {{x_, x.m_[j]}});
562  }
563  }
564  }
565 }
566 
567 // ---------------------------- Diff -------------------------------------//
568 
569 void diff(const DenseMatrix &A, const RCP<const Symbol> &x, DenseMatrix &result,
570  bool diff_cache)
571 {
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);
578  }
579  }
580 }
581 
582 void sdiff(const DenseMatrix &A, const RCP<const Basic> &x, DenseMatrix &result,
583  bool diff_cache)
584 {
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);
593  } else {
594  // TODO: Use a dummy symbol
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),
599  {{x_, x}});
600  }
601  }
602  }
603 }
604 
605 // ----------------------------- Matrix Conjugate ----------------------------//
606 void conjugate_dense(const DenseMatrix &A, DenseMatrix &B)
607 {
608  SYMENGINE_ASSERT(B.col_ == A.col_ and B.row_ == A.row_);
609 
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]);
613 }
614 
615 // ----------------------------- Matrix Transpose ----------------------------//
616 void transpose_dense(const DenseMatrix &A, DenseMatrix &B)
617 {
618  SYMENGINE_ASSERT(B.row_ == A.col_ and B.col_ == A.row_);
619 
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];
623 }
624 
625 // ----------------------------- Matrix Conjugate Transpose -----------------//
626 void conjugate_transpose_dense(const DenseMatrix &A, DenseMatrix &B)
627 {
628  SYMENGINE_ASSERT(B.row_ == A.col_ and B.col_ == A.row_);
629 
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]);
633 }
634 
635 // ------------------------------- Submatrix ---------------------------------//
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)
639 {
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);
645 
646  unsigned row = B.row_, col = B.col_;
647 
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];
651 }
652 
653 // ------------------------------- Matrix Addition ---------------------------//
654 void add_dense_dense(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &C)
655 {
656  SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_ and A.row_ == C.row_
657  and A.col_ == C.col_);
658 
659  unsigned row = A.row_, col = A.col_;
660 
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]);
664  }
665  }
666 }
667 
668 void add_dense_scalar(const DenseMatrix &A, const RCP<const Basic> &k,
669  DenseMatrix &B)
670 {
671  SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_);
672 
673  unsigned row = A.row_, col = A.col_;
674 
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);
678  }
679  }
680 }
681 
682 // ------------------------------- Matrix Multiplication ---------------------//
683 void mul_dense_dense(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &C)
684 {
685  SYMENGINE_ASSERT(A.col_ == B.row_ and C.row_ == A.row_
686  and C.col_ == B.col_);
687 
688  unsigned row = A.row_, col = B.col_;
689 
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; // Integer Zero
694  for (unsigned k = 0; k < A.col_; k++)
695  C.m_[r * col + c]
696  = add(C.m_[r * col + c],
697  mul(A.m_[r * A.col_ + k], B.m_[k * col + c]));
698  }
699  }
700  } else {
701  DenseMatrix tmp = DenseMatrix(A.row_, B.col_);
702  mul_dense_dense(A, B, tmp);
703  C = tmp;
704  }
705 }
706 
707 void elementwise_mul_dense_dense(const DenseMatrix &A, const DenseMatrix &B,
708  DenseMatrix &C)
709 {
710  SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_ and A.row_ == C.row_
711  and A.col_ == C.col_);
712 
713  unsigned row = A.row_, col = A.col_;
714 
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]);
718  }
719  }
720 }
721 
722 void mul_dense_scalar(const DenseMatrix &A, const RCP<const Basic> &k,
723  DenseMatrix &B)
724 {
725  SYMENGINE_ASSERT(A.col_ == B.col_ and A.row_ == B.row_);
726 
727  unsigned row = A.row_, col = A.col_;
728 
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);
732  }
733  }
734 }
735 
736 // ---------------------------- Joining Operations ---------------------------//
737 void DenseMatrix::row_join(const DenseMatrix &B)
738 {
739  this->col_insert(B, col_);
740 }
741 
742 void DenseMatrix::col_join(const DenseMatrix &B)
743 {
744  this->row_insert(B, row_);
745 }
746 
747 void DenseMatrix::row_insert(const DenseMatrix &B, unsigned pos)
748 {
749  SYMENGINE_ASSERT(col_ == B.col_ and pos <= row_)
750 
751  unsigned row = row_, col = col_;
752  this->resize(row_ + B.row_, col_);
753 
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];
757  }
758  }
759 
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];
763  }
764  }
765 }
766 
767 void DenseMatrix::col_insert(const DenseMatrix &B, unsigned pos)
768 {
769  SYMENGINE_ASSERT(row_ == B.row_ and pos <= col_)
770 
771  unsigned row = row_, col = col_;
772  this->resize(row_, col_ + B.col_);
773 
774  for (unsigned i = row; i-- > 0;) {
775  for (unsigned j = col; j-- > 0;) {
776  if (j >= pos) {
777  this->m_[i * (col + B.col_) + j + B.col_]
778  = this->m_[i * col + j];
779  } else {
780  this->m_[i * (col + B.col_) + j] = this->m_[i * col + j];
781  }
782  }
783  }
784 
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];
788  }
789  }
790 }
791 
792 void DenseMatrix::row_del(unsigned k)
793 {
794  SYMENGINE_ASSERT(k < row_)
795 
796  if (row_ == 1)
797  this->resize(0, 0);
798  else {
799  for (unsigned i = k; i < row_ - 1; i++) {
800  row_exchange_dense(*this, i, i + 1);
801  }
802  this->resize(row_ - 1, col_);
803  }
804 }
805 
806 void DenseMatrix::col_del(unsigned k)
807 {
808  SYMENGINE_ASSERT(k < col_)
809 
810  if (col_ == 1)
811  this->resize(0, 0);
812  else {
813  unsigned row = row_, col = col_, m = 0;
814 
815  for (unsigned i = 0; i < row; i++) {
816  for (unsigned j = 0; j < col; j++) {
817  if (j != k) {
818  this->m_[m] = this->m_[i * col + j];
819  m++;
820  }
821  }
822  }
823  this->resize(row_, col_ - 1);
824  }
825 }
826 
827 // -------------------------------- Row Operations ---------------------------//
828 void row_exchange_dense(DenseMatrix &A, unsigned i, unsigned j)
829 {
830  SYMENGINE_ASSERT(i != j and i < A.row_ and j < A.row_);
831 
832  unsigned col = A.col_;
833 
834  for (unsigned k = 0; k < A.col_; k++)
835  std::swap(A.m_[i * col + k], A.m_[j * col + k]);
836 }
837 
838 void row_mul_scalar_dense(DenseMatrix &A, unsigned i, RCP<const Basic> &c)
839 {
840  SYMENGINE_ASSERT(i < A.row_);
841 
842  unsigned col = A.col_;
843 
844  for (unsigned j = 0; j < A.col_; j++)
845  A.m_[i * col + j] = mul(c, A.m_[i * col + j]);
846 }
847 
848 void row_add_row_dense(DenseMatrix &A, unsigned i, unsigned j,
849  RCP<const Basic> &c)
850 {
851  SYMENGINE_ASSERT(i != j and i < A.row_ and j < A.row_);
852 
853  unsigned col = A.col_;
854 
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]));
857 }
858 
859 void permuteFwd(DenseMatrix &A, permutelist &pl)
860 {
861  for (auto &p : pl) {
862  row_exchange_dense(A, p.first, p.second);
863  }
864 }
865 
866 // ----------------------------- Column Operations ---------------------------//
867 void column_exchange_dense(DenseMatrix &A, unsigned i, unsigned j)
868 {
869  SYMENGINE_ASSERT(i != j and i < A.col_ and j < A.col_);
870 
871  unsigned col = A.col_;
872 
873  for (unsigned k = 0; k < A.row_; k++)
874  std::swap(A.m_[k * col + i], A.m_[k * col + j]);
875 }
876 
877 // ------------------------------ Gaussian Elimination -----------------------//
878 void pivoted_gaussian_elimination(const DenseMatrix &A, DenseMatrix &B,
879  permutelist &pl)
880 {
881  SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_);
882 
883  unsigned row = A.row_, col = A.col_;
884  unsigned index = 0, i, j, k;
885  B.m_ = A.m_;
886 
887  RCP<const Basic> scale;
888 
889  for (i = 0; i < col - 1; i++) {
890  if (index == row)
891  break;
892 
893  k = pivot(B, index, i);
894  if (k == row)
895  continue;
896  if (k != index) {
897  row_exchange_dense(B, k, index);
898  pl.push_back({k, index});
899  }
900 
901  scale = div(one, B.m_[index * col + i]);
902  row_mul_scalar_dense(B, index, scale);
903 
904  for (j = i + 1; j < row; j++) {
905  for (k = i + 1; k < col; k++)
906  B.m_[j * 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;
910  }
911 
912  index++;
913  }
914 }
915 
916 // Algorithm 1, page 12, Nakos, G. C., Turner, P. R., Williams, R. M. (1997).
917 // Fraction-free algorithms for linear and polynomial equations.
918 // ACM SIGSAM Bulletin, 31(3), 11–19. doi:10.1145/271130.271133.
919 void fraction_free_gaussian_elimination(const DenseMatrix &A, DenseMatrix &B)
920 {
921  SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_);
922 
923  unsigned col = A.col_;
924  B.m_ = A.m_;
925 
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++) {
929  B.m_[j * 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]));
932  if (i > 0)
933  B.m_[j * col + k]
934  = div(B.m_[j * col + k], B.m_[i * col - col + i - 1]);
935  }
936  B.m_[j * col + i] = zero;
937  }
938 }
939 
940 // Pivoted version of `fraction_free_gaussian_elimination`
941 void pivoted_fraction_free_gaussian_elimination(const DenseMatrix &A,
942  DenseMatrix &B, permutelist &pl)
943 {
944  SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_);
945 
946  unsigned col = A.col_, row = A.row_;
947  unsigned index = 0, i, k, j;
948  B.m_ = A.m_;
949 
950  for (i = 0; i < col - 1; i++) {
951  if (index == row)
952  break;
953 
954  k = pivot(B, index, i);
955  if (k == row)
956  continue;
957  if (k != index) {
958  row_exchange_dense(B, k, index);
959  pl.push_back({k, index});
960  }
961 
962  for (j = i + 1; j < row; j++) {
963  for (k = i + 1; k < col; k++) {
964  B.m_[j * 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]));
967  if (i > 0)
968  B.m_[j * col + k]
969  = div(B.m_[j * col + k], B.m_[i * col - col + i - 1]);
970  }
971  B.m_[j * col + i] = zero;
972  }
973 
974  index++;
975  }
976 }
977 
978 void pivoted_gauss_jordan_elimination(const DenseMatrix &A, DenseMatrix &B,
979  permutelist &pl)
980 {
981  SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_);
982 
983  unsigned row = A.row_, col = A.col_;
984  unsigned index = 0, i, j, k;
985  RCP<const Basic> scale;
986  B.m_ = A.m_;
987 
988  for (i = 0; i < col; i++) {
989  if (index == row)
990  break;
991 
992  k = pivot(B, index, i);
993  if (k == row)
994  continue;
995  if (k != index) {
996  row_exchange_dense(B, k, index);
997  pl.push_back({k, index});
998  }
999 
1000  scale = div(one, B.m_[index * col + i]);
1001  row_mul_scalar_dense(B, index, scale);
1002 
1003  for (j = 0; j < row; j++) {
1004  if (j == index)
1005  continue;
1006 
1007  scale = mul(minus_one, B.m_[j * col + i]);
1008  row_add_row_dense(B, j, index, scale);
1009  }
1010 
1011  index++;
1012  }
1013 }
1014 
1015 // Algorithm 2, page 13, Nakos, G. C., Turner, P. R., Williams, R. M. (1997).
1016 // Fraction-free algorithms for linear and polynomial equations.
1017 // ACM SIGSAM Bulletin, 31(3), 11–19. doi:10.1145/271130.271133.
1018 void fraction_free_gauss_jordan_elimination(const DenseMatrix &A,
1019  DenseMatrix &B)
1020 {
1021  SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_);
1022 
1023  unsigned row = A.row_, col = A.col_;
1024  unsigned i, j, k;
1025  RCP<const Basic> d;
1026 
1027  B.m_ = A.m_;
1028 
1029  for (i = 0; i < col; i++) {
1030  if (i > 0)
1031  d = B.m_[i * col - col + i - 1];
1032  for (j = 0; j < row; j++)
1033  if (j != i)
1034  for (k = 0; k < col; k++) {
1035  if (k != i) {
1036  B.m_[j * 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]));
1039  if (i > 0)
1040  B.m_[j * col + k] = div(B.m_[j * col + k], d);
1041  }
1042  }
1043  for (j = 0; j < row; j++)
1044  if (j != i)
1045  B.m_[j * col + i] = zero;
1046  }
1047 }
1048 
1049 void pivoted_fraction_free_gauss_jordan_elimination(const DenseMatrix &A,
1050  DenseMatrix &B,
1051  permutelist &pl)
1052 {
1053  SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_);
1054 
1055  unsigned row = A.row_, col = A.col_;
1056  unsigned index = 0, i, k, j;
1057  RCP<const Basic> d;
1058 
1059  B.m_ = A.m_;
1060 
1061  for (i = 0; i < col; i++) {
1062  if (index == row)
1063  break;
1064 
1065  k = pivot(B, index, i);
1066  if (k == row)
1067  continue;
1068  if (k != index) {
1069  row_exchange_dense(B, k, index);
1070  pl.push_back({k, index});
1071  }
1072 
1073  for (j = 0; j < row; j++) {
1074  if (j == index)
1075  continue;
1076  for (k = 0; k < col; k++) {
1077  if (k == i)
1078  continue;
1079  B.m_[j * 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]));
1082  if (index == 0)
1083  continue;
1084  B.m_[j * col + k] = div(B.m_[j * col + k], d);
1085  }
1086  }
1087 
1088  d = B.m_[index * col + i];
1089 
1090  for (j = 0; j < row; j++) {
1091  if (j == index)
1092  continue;
1093  B.m_[j * col + i] = zero;
1094  }
1095  index++;
1096  }
1097 }
1098 
1099 unsigned pivot(DenseMatrix &B, unsigned r, unsigned c)
1100 {
1101  for (unsigned k = r; k < B.row_; k++) {
1102  if (!is_true(is_zero(*(B.m_[k * B.col_ + c])))) {
1103  return k;
1104  }
1105  }
1106  return B.row_;
1107 }
1108 
1109 void reduced_row_echelon_form(const DenseMatrix &A, DenseMatrix &b,
1110  vec_uint &pivot_cols, bool normalize_last)
1111 {
1112  permutelist pl;
1113  if (normalize_last) {
1114  pivoted_fraction_free_gauss_jordan_elimination(A, b, pl);
1115  } else {
1116  pivoted_gauss_jordan_elimination(A, b, pl);
1117  }
1118  unsigned row = 0;
1119  for (unsigned col = 0; col < b.col_ && row < b.row_; col++) {
1120  if (is_true(is_zero(*b.get(row, col))))
1121  continue;
1122  pivot_cols.push_back(col);
1123  if (row == 0 and normalize_last) {
1124  RCP<const Basic> m = div(one, b.get(row, col));
1125  b.mul_scalar(m, b);
1126  }
1127  row++;
1128  }
1129 }
1130 
1131 // --------------------------- Solve Ax = b ---------------------------------//
1132 // Assuming A is a diagonal square matrix
1133 void diagonal_solve(const DenseMatrix &A, const DenseMatrix &b, DenseMatrix &x)
1134 {
1135  SYMENGINE_ASSERT(A.row_ == A.col_);
1136  SYMENGINE_ASSERT(x.row_ == A.col_ and x.col_ == b.col_);
1137 
1138  const unsigned sys = b.col_;
1139 
1140  // No checks are done to see if the diagonal entries are zero
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]);
1144  }
1145  }
1146 }
1147 
1148 void back_substitution(const DenseMatrix &U, const DenseMatrix &b,
1149  DenseMatrix &x)
1150 {
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_);
1154 
1155  const unsigned col = U.col_;
1156  const unsigned sys = b.col_;
1157  x.m_ = b.m_;
1158 
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++)
1162  x.m_[i * sys + k]
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]);
1166  }
1167  }
1168 }
1169 
1170 void forward_substitution(const DenseMatrix &A, const DenseMatrix &b,
1171  DenseMatrix &x)
1172 {
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_);
1176 
1177  unsigned col = A.col_;
1178  const unsigned sys = b.col_;
1179  x.m_ = b.m_;
1180 
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++) {
1184  x.m_[j * sys + k]
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]));
1187  if (i > 0)
1188  x.m_[j * sys + k]
1189  = div(x.m_[j * sys + k], A.m_[(i - 1) * col + i - 1]);
1190  }
1191  }
1192  }
1193 }
1194 
1195 void fraction_free_gaussian_elimination_solve(const DenseMatrix &A,
1196  const DenseMatrix &b,
1197  DenseMatrix &x)
1198 {
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_);
1202 
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_);
1206 
1207  for (i = 0; i < col - 1; i++)
1208  for (j = i + 1; j < col; j++) {
1209  for (k = 0; k < bcol; k++) {
1210  b_.m_[j * 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]));
1213  if (i > 0)
1214  b_.m_[j * bcol + k] = div(b_.m_[j * bcol + k],
1215  A_.m_[i * col - col + i - 1]);
1216  }
1217 
1218  for (k = i + 1; k < col; k++) {
1219  A_.m_[j * 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]));
1222  if (i > 0)
1223  A_.m_[j * col + k]
1224  = div(A_.m_[j * col + k], A_.m_[i * col - col + i - 1]);
1225  }
1226  A_.m_[j * col + i] = zero;
1227  }
1228 
1229  for (i = 0; i < col * bcol; i++)
1230  x.m_[i] = zero; // Integer zero;
1231 
1232  for (k = 0; k < bcol; k++) {
1233  for (i = col - 1; i >= 0; i--) {
1234  for (j = i + 1; j < col; j++)
1235  b_.m_[i * bcol + k]
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]);
1239  }
1240  }
1241 }
1242 
1243 void fraction_free_gauss_jordan_solve(const DenseMatrix &A,
1244  const DenseMatrix &b, DenseMatrix &x,
1245  bool pivot)
1246 {
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_);
1250 
1251  unsigned i, j, k, p, col = A.col_, bcol = b.col_;
1252  RCP<const Basic> d;
1253  DenseMatrix A_ = DenseMatrix(A.row_, A.col_, A.m_);
1254  DenseMatrix b_ = DenseMatrix(b.row_, b.col_, b.m_);
1255 
1256  for (i = 0; i < col; i++) {
1257  if (i > 0)
1258  d = A_.m_[i * col - col + i - 1];
1259  if (pivot) {
1260  p = i;
1261  while (p < col and eq(*A_.m_[p * col + i], *zero)) {
1262  p++;
1263  }
1264  SYMENGINE_ASSERT(p != col);
1265  if (p != i) {
1266  // pivot A
1267  for (k = i; k < col; k++) {
1268  std::swap(A_.m_[p * col + k], A_.m_[i * col + k]);
1269  }
1270  for (k = 0; k < bcol; k++) {
1271  std::swap(b_.m_[p * bcol + k], b_.m_[i * bcol + k]);
1272  }
1273  }
1274  }
1275  for (j = 0; j < col; j++)
1276  if (j != i) {
1277  for (k = 0; k < bcol; k++) {
1278  b_.m_[j * 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]));
1281  if (i > 0)
1282  b_.m_[j * bcol + k] = div(b_.m_[j * bcol + k], d);
1283  }
1284 
1285  for (k = 0; k < col; k++) {
1286  if (k != i) {
1287  A_.m_[j * 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]));
1290  if (i > 0)
1291  A_.m_[j * col + k] = div(A_.m_[j * col + k], d);
1292  }
1293  }
1294  }
1295 
1296  for (j = 0; j < col; j++)
1297  if (j != i)
1298  A_.m_[j * col + i] = zero;
1299  }
1300 
1301  // No checks are done to see if the diagonal entries are 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]);
1305 }
1306 
1307 void fraction_free_LU_solve(const DenseMatrix &A, const DenseMatrix &b,
1308  DenseMatrix &x)
1309 {
1310  DenseMatrix LU = DenseMatrix(A.nrows(), A.ncols());
1311  DenseMatrix x_ = DenseMatrix(b.nrows(), b.ncols());
1312 
1313  fraction_free_LU(A, LU);
1314  forward_substitution(LU, b, x_);
1315  back_substitution(LU, x_, x);
1316 }
1317 
1318 void LU_solve(const DenseMatrix &A, const DenseMatrix &b, DenseMatrix &x)
1319 {
1320  DenseMatrix L = DenseMatrix(A.nrows(), A.ncols());
1321  DenseMatrix U = DenseMatrix(A.nrows(), A.ncols());
1322  DenseMatrix x_ = DenseMatrix(b.nrows(), b.ncols());
1323 
1324  LU(A, L, U);
1325  forward_substitution(L, b, x_);
1326  back_substitution(U, x_, x);
1327 }
1328 
1329 void pivoted_LU_solve(const DenseMatrix &A, const DenseMatrix &b,
1330  DenseMatrix &x)
1331 {
1332  DenseMatrix L = DenseMatrix(A.nrows(), A.ncols());
1333  DenseMatrix U = DenseMatrix(A.nrows(), A.ncols());
1334  DenseMatrix x_ = DenseMatrix(b);
1335  permutelist pl;
1336 
1337  pivoted_LU(A, L, U, pl);
1338  permuteFwd(x_, pl);
1339  forward_substitution(L, x_, x_);
1340  back_substitution(U, x_, x);
1341 }
1342 
1343 void LDL_solve(const DenseMatrix &A, const DenseMatrix &b, DenseMatrix &x)
1344 {
1345  DenseMatrix L = DenseMatrix(A.nrows(), A.ncols());
1346  DenseMatrix D = DenseMatrix(A.nrows(), A.ncols());
1347  DenseMatrix x_ = DenseMatrix(b.nrows(), b.ncols());
1348 
1349  if (not is_symmetric_dense(A))
1350  throw SymEngineException("Matrix must be symmetric");
1351 
1352  LDL(A, L, D);
1353  forward_substitution(L, b, x);
1354  diagonal_solve(D, x, x_);
1355  transpose_dense(L, D);
1356  back_substitution(D, x_, x);
1357 }
1358 
1359 // --------------------------- Matrix Decomposition --------------------------//
1360 
1361 // Algorithm 3, page 14, Nakos, G. C., Turner, P. R., Williams, R. M. (1997).
1362 // Fraction-free algorithms for linear and polynomial equations.
1363 // ACM SIGSAM Bulletin, 31(3), 11–19. doi:10.1145/271130.271133.
1364 // This algorithms is not a true factorization of the matrix A(i.e. A != LU))
1365 // but can be used to solve linear systems by applying forward and backward
1366 // substitutions respectively.
1367 void fraction_free_LU(const DenseMatrix &A, DenseMatrix &LU)
1368 {
1369  SYMENGINE_ASSERT(A.row_ == A.col_ and LU.row_ == LU.col_
1370  and A.row_ == LU.row_);
1371 
1372  unsigned n = A.row_;
1373  unsigned i, j, k;
1374 
1375  LU.m_ = A.m_;
1376 
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]));
1382  if (i)
1383  LU.m_[j * n + k]
1384  = div(LU.m_[j * n + k], LU.m_[i * n - n + i - 1]);
1385  }
1386 }
1387 
1388 // SymPy LUDecomposition algorithm, in
1389 // sympy.matrices.matrices.Matrix.LUdecomposition
1390 // with no pivoting
1391 void LU(const DenseMatrix &A, DenseMatrix &L, DenseMatrix &U)
1392 {
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_);
1396 
1397  unsigned n = A.row_;
1398  unsigned i, j, k;
1399  RCP<const Basic> scale;
1400 
1401  U.m_ = A.m_;
1402 
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]));
1408 
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]));
1413  }
1414 
1415  scale = div(one, U.m_[j * n + j]);
1416 
1417  for (i = j + 1; i < n; i++)
1418  U.m_[i * n + j] = mul(U.m_[i * n + j], scale);
1419  }
1420 
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; // Integer zero
1425  }
1426  L.m_[i * n + i] = one; // Integer one
1427  for (j = i + 1; j < n; j++)
1428  L.m_[i * n + j] = zero; // Integer zero
1429  }
1430 }
1431 
1432 // SymPy LUDecomposition_Simple algorithm, in
1433 // sympy.matrices.matrices.Matrix.LUdecomposition_Simple with pivoting.
1434 // P must be an initialized matrix and will be permuted.
1435 void pivoted_LU(const DenseMatrix &A, DenseMatrix &LU, permutelist &pl)
1436 {
1437  SYMENGINE_ASSERT(A.row_ == A.col_ and LU.row_ == LU.col_);
1438  SYMENGINE_ASSERT(A.row_ == LU.row_);
1439 
1440  unsigned n = A.row_;
1441  unsigned i, j, k;
1442  RCP<const Basic> scale;
1443  int pivot;
1444 
1445  LU.m_ = A.m_;
1446 
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]));
1452  pivot = -1;
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]));
1457  }
1458  if (pivot == -1 and !is_true(is_zero(*LU.m_[i * n + j])))
1459  pivot = i;
1460  }
1461  if (pivot == -1)
1462  throw SymEngineException("Matrix is rank deficient");
1463  if (pivot - j != 0) { // row must be swapped
1464  row_exchange_dense(LU, pivot, j);
1465  pl.push_back({pivot, j});
1466  }
1467  scale = div(one, LU.m_[j * n + j]);
1468 
1469  for (i = j + 1; i < n; i++)
1470  LU.m_[i * n + j] = mul(LU.m_[i * n + j], scale);
1471  }
1472 }
1473 
1474 // SymPy LUDecomposition algorithm, in
1475 // sympy.matrices.matrices.Matrix.LUdecomposition_Simple with pivoting.
1476 // P must be an initialized matrix and will be permuted.
1477 void pivoted_LU(const DenseMatrix &A, DenseMatrix &L, DenseMatrix &U,
1478  permutelist &pl)
1479 {
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_);
1483 
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; // Integer zero
1490  }
1491  L.m_[i * n + i] = one; // Integer one
1492  for (unsigned j = i + 1; j < n; j++)
1493  L.m_[i * n + j] = zero; // Integer zero
1494  }
1495 }
1496 
1497 // SymPy's fraction free LU decomposition, without pivoting
1498 // sympy.matrices.matrices.MatrixBase.LUdecompositionFF
1499 // W. Zhou & D.J. Jeffrey, "Fraction-free matrix factors: new forms for LU and
1500 // QR factors".
1501 // Frontiers in Computer Science in China, Vol 2, no. 1, pp. 67-80, 2008.
1502 void fraction_free_LDU(const DenseMatrix &A, DenseMatrix &L, DenseMatrix &D,
1503  DenseMatrix &U)
1504 {
1505  SYMENGINE_ASSERT(A.row_ == L.row_ and A.row_ == U.row_);
1506  SYMENGINE_ASSERT(A.col_ == L.col_ and A.col_ == U.col_);
1507 
1508  unsigned row = A.row_, col = A.col_;
1509  unsigned i, j, k;
1510  RCP<const Basic> old = integer(1);
1511 
1512  U.m_ = A.m_;
1513 
1514  // Initialize L
1515  for (i = 0; i < row; i++)
1516  for (j = 0; j < row; j++)
1517  if (i != j)
1518  L.m_[i * col + j] = zero;
1519  else
1520  L.m_[i * col + i] = one;
1521 
1522  // Initialize D
1523  for (i = 0; i < row * row; i++)
1524  D.m_[i] = zero; // Integer zero
1525 
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]);
1529 
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++)
1533  U.m_[i * 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])),
1536  old);
1537  U.m_[i * col + k] = zero; // Integer zero
1538  }
1539 
1540  old = U.m_[k * col + k];
1541  }
1542 
1543  D.m_[row * col - col + row - 1] = old;
1544 }
1545 
1546 // SymPy's QRecomposition in sympy.matrices.matrices.MatrixBase.QRdecomposition
1547 // Rank check is not performed
1548 void QR(const DenseMatrix &A, DenseMatrix &Q, DenseMatrix &R)
1549 {
1550  unsigned row = A.row_;
1551  unsigned col = A.col_;
1552 
1553  SYMENGINE_ASSERT(Q.row_ == row and Q.col_ == col and R.row_ == col
1554  and R.col_ == col);
1555 
1556  unsigned i, j, k;
1557  RCP<const Basic> t;
1558  std::vector<RCP<const Basic>> tmp(row);
1559 
1560  // Initialize Q
1561  for (i = 0; i < row * col; i++)
1562  Q.m_[i] = zero;
1563 
1564  // Initialize R
1565  for (i = 0; i < col * col; i++)
1566  R.m_[i] = zero;
1567 
1568  for (j = 0; j < col; j++) {
1569  // Use submatrix for this
1570  for (k = 0; k < row; k++)
1571  tmp[k] = A.m_[k * col + j];
1572 
1573  for (i = 0; i < j; i++) {
1574  t = zero;
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)));
1579  }
1580 
1581  // calculate norm
1582  t = zero;
1583  for (k = 0; k < row; k++)
1584  t = add(t, pow(tmp[k], integer(2)));
1585 
1586  t = pow(t, div(one, integer(2)));
1587 
1588  R.m_[j * col + j] = t;
1589  for (k = 0; k < row; k++)
1590  Q.m_[k * col + j] = div(tmp[k], t);
1591 
1592  for (i = 0; i < j; i++) {
1593  t = zero;
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;
1597  }
1598  }
1599 }
1600 
1601 // SymPy's LDL decomposition, Assuming A is a symmetric, square, positive
1602 // definite non singular matrix
1603 void LDL(const DenseMatrix &A, DenseMatrix &L, DenseMatrix &D)
1604 {
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_);
1608 
1609  unsigned col = A.col_;
1610  unsigned i, k, j;
1611  RCP<const Basic> sum;
1612  RCP<const Basic> i2 = integer(2);
1613 
1614  // Initialize D
1615  for (i = 0; i < col; i++)
1616  for (j = 0; j < col; j++)
1617  D.m_[i * col + j] = zero; // Integer zero
1618 
1619  // Initialize L
1620  for (i = 0; i < col; i++)
1621  for (j = 0; j < col; j++)
1622  L.m_[i * col + j] = (i != j) ? zero : one;
1623 
1624  for (i = 0; i < col; i++) {
1625  for (j = 0; j < i; j++) {
1626  sum = zero;
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]));
1630  L.m_[i * col + j]
1631  = mul(div(one, D.m_[j * col + j]), sub(A.m_[i * col + j], sum));
1632  }
1633  sum = zero;
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);
1637  }
1638 }
1639 
1640 // SymPy's cholesky decomposition
1641 void cholesky(const DenseMatrix &A, DenseMatrix &L)
1642 {
1643  SYMENGINE_ASSERT(A.row_ == A.col_);
1644  SYMENGINE_ASSERT(L.row_ == A.row_ and L.col_ == A.row_);
1645 
1646  unsigned col = A.col_;
1647  unsigned i, j, k;
1648  RCP<const Basic> sum;
1649  RCP<const Basic> i2 = integer(2);
1650  RCP<const Basic> half = div(one, i2);
1651 
1652  // Initialize L
1653  for (i = 0; i < col; i++)
1654  for (j = 0; j < col; j++)
1655  L.m_[i * col + j] = zero;
1656 
1657  for (i = 0; i < col; i++) {
1658  for (j = 0; j < i; j++) {
1659  sum = zero;
1660  for (k = 0; k < j; k++)
1661  sum = add(sum, mul(L.m_[i * col + k], L.m_[j * col + k]));
1662 
1663  L.m_[i * col + j]
1664  = mul(div(one, L.m_[j * col + j]), sub(A.m_[i * col + j], sum));
1665  }
1666  sum = zero;
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);
1670  }
1671 }
1672 
1673 // Matrix Queries
1674 bool is_symmetric_dense(const DenseMatrix &A)
1675 {
1676  if (A.col_ != A.row_)
1677  return false;
1678 
1679  unsigned col = A.col_;
1680  bool sym = true;
1681 
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]))) {
1685  sym = false;
1686  break;
1687  }
1688 
1689  return sym;
1690 }
1691 
1692 // ----------------------------- Determinant ---------------------------------//
1693 RCP<const Basic> det_bareis(const DenseMatrix &A)
1694 {
1695  SYMENGINE_ASSERT(A.row_ == A.col_);
1696 
1697  unsigned n = A.row_;
1698 
1699  if (n == 1) {
1700  return A.m_[0];
1701  } else if (n == 2) {
1702  // If A = [[a, b], [c, d]] then det(A) = ad - bc
1703  return sub(mul(A.m_[0], A.m_[3]), mul(A.m_[1], A.m_[2]));
1704  } else if (n == 3) {
1705  // if A = [[a, b, c], [d, e, f], [g, h, i]] then
1706  // det(A) = (aei + bfg + cdh) - (ceg + bdi + afh)
1707  return sub(add(add(mul(mul(A.m_[0], A.m_[4]), A.m_[8]),
1708  mul(mul(A.m_[1], A.m_[5]), A.m_[6])),
1709  mul(mul(A.m_[2], A.m_[3]), A.m_[7])),
1710  add(add(mul(mul(A.m_[2], A.m_[4]), A.m_[6]),
1711  mul(mul(A.m_[1], A.m_[3]), A.m_[8])),
1712  mul(mul(A.m_[0], A.m_[5]), A.m_[7])));
1713  } else {
1714 
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]);
1719  }
1720  return det;
1721  }
1722 
1723  DenseMatrix B = DenseMatrix(n, n, A.m_);
1724  unsigned i, sign = 1;
1725  RCP<const Basic> d;
1726 
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);
1732  sign *= -1;
1733  break;
1734  }
1735  if (i == n)
1736  return zero;
1737  }
1738 
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]));
1743  if (k > 0)
1744  d = div(d, B.m_[(k - 1) * n + k - 1]);
1745  B.m_[i * n + j] = d;
1746  }
1747  }
1748  }
1749 
1750  return (sign == 1) ? B.m_[n * n - 1] : mul(minus_one, B.m_[n * n - 1]);
1751  }
1752 }
1753 
1754 // Returns the coefficients of characterisitc polynomials of leading principal
1755 // minors of Matrix A as elements in `polys`. Principal leading minor of kth
1756 // order is the submatrix of A obtained by deleting last n-k rows and columns
1757 // from A. Here `n` is the dimension of the square matrix A.
1758 void berkowitz(const DenseMatrix &A, std::vector<DenseMatrix> &polys)
1759 {
1760  SYMENGINE_ASSERT(A.row_ == A.col_);
1761 
1762  unsigned col = A.col_;
1763  unsigned i, k, l, m;
1764 
1766  std::vector<DenseMatrix> transforms;
1768 
1769  for (unsigned n = col; n > 1; n--) {
1770  items.clear();
1771  k = n - 1;
1772  DenseMatrix T = DenseMatrix(n + 1, n);
1773  DenseMatrix C = DenseMatrix(k, 1);
1774 
1775  // Initialize T and C
1776  for (i = 0; i < n * (n + 1); i++)
1777  T.m_[i] = zero;
1778  for (i = 0; i < k; i++)
1779  C.m_[i] = A.m_[i * col + k];
1780  items.push_back(C);
1781 
1782  for (i = 0; i < n - 2; i++) {
1783  DenseMatrix B = DenseMatrix(k, 1);
1784  for (l = 0; l < k; l++) {
1785  B.m_[l] = zero;
1786  for (m = 0; m < k; m++)
1787  B.m_[l]
1788  = add(B.m_[l], mul(A.m_[l * col + m], items[i].m_[m]));
1789  }
1790  items.push_back(B);
1791  }
1792 
1793  items_.clear();
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]));
1798  items_.push_back(mul(minus_one, element));
1799  }
1800  items_.insert(items_.begin(), mul(minus_one, A.m_[k * col + k]));
1801  items_.insert(items_.begin(), one);
1802 
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];
1806  }
1807 
1808  transforms.push_back(T);
1809  }
1810 
1811  polys.push_back(DenseMatrix(2, 1, {one, mul(A.m_[0], minus_one)}));
1812 
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);
1817 
1818  for (l = 0; l < t_row; l++) {
1819  B.m_[l] = zero;
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],
1823  polys[i].m_[m]));
1824  B.m_[l] = expand(B.m_[l]);
1825  }
1826  }
1827  polys.push_back(B);
1828  }
1829 }
1830 
1831 RCP<const Basic> det_berkowitz(const DenseMatrix &A)
1832 {
1834 
1835  berkowitz(A, polys);
1836  DenseMatrix poly = polys[polys.size() - 1];
1837 
1838  if (polys.size() % 2 == 1)
1839  return mul(minus_one, poly.get(poly.nrows() - 1, 0));
1840 
1841  return poly.get(poly.nrows() - 1, 0);
1842 }
1843 
1844 void char_poly(const DenseMatrix &A, DenseMatrix &B)
1845 {
1846  SYMENGINE_ASSERT(B.ncols() == 1 and B.nrows() == A.nrows() + 1);
1847  SYMENGINE_ASSERT(A.nrows() == A.ncols());
1848 
1850 
1851  berkowitz(A, polys);
1852  B = polys[polys.size() - 1];
1853 }
1854 
1855 void inverse_fraction_free_LU(const DenseMatrix &A, DenseMatrix &B)
1856 {
1857  SYMENGINE_ASSERT(A.row_ == A.col_ and B.row_ == B.col_
1858  and B.row_ == A.row_);
1859 
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);
1865 
1866  // Initialize matrices
1867  for (i = 0; i < n * n; i++) {
1868  LU.m_[i] = zero;
1869  B.m_[i] = zero;
1870  }
1871  for (i = 0; i < n; i++) {
1872  e.m_[i] = zero;
1873  x.m_[i] = zero;
1874  x_.m_[i] = zero;
1875  }
1876 
1877  fraction_free_LU(A, LU);
1878 
1879  // We solve AX_{i} = e_{i} for i = 1, 2, .. n and combine the row vectors
1880  // X_{1}, X_{2}, ... X_{n} to form the inverse of A. Here, e_{i}'s are the
1881  // elements of the standard basis.
1882  for (unsigned j = 0; j < n; j++) {
1883  e.m_[j] = one;
1884 
1885  forward_substitution(LU, e, x_);
1886  back_substitution(LU, x_, x);
1887 
1888  for (i = 0; i < n; i++)
1889  B.m_[i * n + j] = x.m_[i];
1890 
1891  e.m_[j] = zero;
1892  }
1893 }
1894 
1895 void inverse_LU(const DenseMatrix &A, DenseMatrix &B)
1896 {
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_);
1900  eye(e);
1901  LU_solve(A, e, B);
1902 }
1903 
1904 void inverse_pivoted_LU(const DenseMatrix &A, DenseMatrix &B)
1905 {
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_);
1909  eye(e);
1910  pivoted_LU_solve(A, e, B);
1911 }
1912 
1913 void inverse_gauss_jordan(const DenseMatrix &A, DenseMatrix &B)
1914 {
1915  SYMENGINE_ASSERT(A.row_ == A.col_ and B.row_ == B.col_
1916  and B.row_ == A.row_);
1917 
1918  unsigned n = A.row_;
1919  DenseMatrix e = DenseMatrix(n, n);
1920 
1921  // Initialize matrices
1922  for (unsigned i = 0; i < n; i++)
1923  for (unsigned j = 0; j < n; j++) {
1924  if (i != j) {
1925  e.m_[i * n + j] = zero;
1926  } else {
1927  e.m_[i * n + i] = one;
1928  }
1929  B.m_[i * n + j] = zero;
1930  }
1931 
1932  fraction_free_gauss_jordan_solve(A, e, B);
1933 }
1934 
1935 // ----------------------- Vector-specific Methods --------------------------//
1936 
1937 void dot(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &C)
1938 {
1939  if (A.col_ == B.row_) {
1940  if (B.col_ != 1) {
1941  DenseMatrix tmp1 = DenseMatrix(A.col_, A.row_);
1942  A.transpose(tmp1);
1943  DenseMatrix tmp2 = DenseMatrix(B.col_, B.row_);
1944  B.transpose(tmp2);
1945  C.resize(tmp1.row_, tmp2.col_);
1946  mul_dense_dense(tmp1, tmp2, C);
1947  } else {
1948  C.resize(A.row_, B.col_);
1949  mul_dense_dense(A, B, C);
1950  }
1951  C.resize(1, C.row_ * C.col_);
1952  } else if (A.col_ == B.col_) {
1953  DenseMatrix tmp2 = DenseMatrix(B.col_, B.row_);
1954  B.transpose(tmp2);
1955  dot(A, tmp2, C);
1956  } else if (A.row_ == B.row_) {
1957  DenseMatrix tmp1 = DenseMatrix(A.col_, A.row_);
1958  A.transpose(tmp1);
1959  dot(tmp1, B, C);
1960  } else {
1961  throw SymEngineException("Dimensions incorrect for dot product");
1962  }
1963 }
1964 
1965 void cross(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &C)
1966 {
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_));
1969 
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]));
1973 }
1974 
1975 RCP<const Set> eigen_values(const DenseMatrix &A)
1976 {
1977  unsigned n = A.nrows();
1978  if (A.is_lower() or A.is_upper()) {
1979  RCP<const Set> eigenvals = emptyset();
1980  set_basic x;
1981  for (unsigned i = 0; i < n; ++i) {
1982  x.insert(A.get(i, i));
1983  }
1984  eigenvals = finiteset(x);
1985  return eigenvals;
1986  }
1987 
1988  DenseMatrix B = DenseMatrix(A.nrows() + 1, 1);
1989  char_poly(A, B);
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);
1996 }
1997 
1998 // ------------------------- NumPy-like functions ----------------------------//
1999 
2000 // Mimic `eye` function in NumPy
2001 void eye(DenseMatrix &A, int k)
2002 {
2003  if ((k >= 0 and (unsigned) k >= A.col_) or k + A.row_ <= 0) {
2004  zeros(A);
2005  }
2006 
2007  vec_basic v = vec_basic(k > 0 ? A.col_ - k : A.row_ + k, one);
2008 
2009  diag(A, v, k);
2010 }
2011 
2012 // Create diagonal matrices directly
2013 void diag(DenseMatrix &A, vec_basic &v, int k)
2014 {
2015  SYMENGINE_ASSERT(v.size() > 0);
2016 
2017  unsigned k_ = std::abs(k);
2018 
2019  if (k >= 0) {
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;
2024  } else {
2025  A.m_[i * A.col_ + j] = v[k - k_];
2026  }
2027  }
2028  k++;
2029  }
2030  } else {
2031  k = -k;
2032 
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;
2037  } else {
2038  A.m_[i * A.col_ + j] = v[k - k_];
2039  }
2040  }
2041  k++;
2042  }
2043  }
2044 }
2045 
2046 // Create a matrix filled with ones
2047 void ones(DenseMatrix &A)
2048 {
2049  for (unsigned i = 0; i < A.row_ * A.col_; i++) {
2050  A.m_[i] = one;
2051  }
2052 }
2053 
2054 // Create a matrix filled with zeros
2055 void zeros(DenseMatrix &A)
2056 {
2057  for (unsigned i = 0; i < A.row_ * A.col_; i++) {
2058  A.m_[i] = zero;
2059  }
2060 }
2061 
2062 std::string DenseMatrix::dumps() const
2063 {
2064 #if HAVE_SYMENGINE_RTTI
2065  std::ostringstream oss;
2066  unsigned short major = SYMENGINE_MAJOR_VERSION;
2067  unsigned short minor = SYMENGINE_MINOR_VERSION;
2069  major, minor, row_, col_, m_);
2070  return oss.str();
2071 #else
2072  throw NotImplementedError("Serialization not implemented in no-rtti mode");
2073 #endif
2074 }
2075 
2076 DenseMatrix DenseMatrix::loads(const std::string &serialized)
2077 {
2078 #if HAVE_SYMENGINE_RTTI
2079  unsigned short major, minor;
2080  unsigned row, col;
2081  vec_basic obj;
2082  std::istringstream iss(serialized);
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 << "."
2091  << minor << ".");
2092  }
2093  iarchive(row, col, obj);
2094  return DenseMatrix(row, col, std::move(obj));
2095 #else
2096  throw NotImplementedError("Serialization not implemented in no-rtti mode");
2097 #endif
2098 }
2099 
2100 } // namespace SymEngine
Classes and functions relating to the binary operation of addition.
T begin(T... args)
RCP< const Basic > sub(const RCP< const Basic > &a, const RCP< const Basic > &b)
Substracts b from a.
Definition: add.cpp:495
RCP< const Basic > add(const RCP< const Basic > &a, const RCP< const Basic > &b)
Adds two objects (safely).
Definition: add.cpp:425
T clear(T... args)
T div(T... args)
T get(T... args)
T insert(T... args)
T move(T... args)
Main namespace for SymEngine package.
Definition: add.cpp:19
std::enable_if< std::is_integral< T >::value, RCP< const Integer > >::type integer(T i)
Definition: integer.h:197
bool is_number_and_zero(const Basic &b)
Definition: number.h:139
RCP< const Symbol > symbol(const std::string &name)
inline version to return Symbol
Definition: symbol.h:82
bool eq(const Basic &a, const Basic &b)
Checks equality for a and b
Definition: basic-inl.h:21
RCP< const EmptySet > emptyset()
Definition: sets.h:590
RCP< const Basic > sign(const RCP< const Basic > &arg)
Canonicalize Sign.
Definition: functions.cpp:527
RCP< const Basic > abs(const RCP< const Basic > &arg)
Canonicalize Abs:
Definition: functions.cpp:3492
RCP< const Basic > mul(const RCP< const Basic > &a, const RCP< const Basic > &b)
Multiplication.
Definition: mul.cpp:352
void insert(T1 &m, const T2 &first, const T3 &second)
Definition: dict.h:83
tribool is_zero(const Basic &b, const Assumptions *assumptions=nullptr)
Check if a number is zero.
RCP< const Set > finiteset(const set_basic &container)
Definition: sets.h:602
RCP< const Basic > expand(const RCP< const Basic > &self, bool deep=true)
Expands self
Definition: expand.cpp:369
RCP< const Basic > ssubs(const RCP< const Basic > &x, const map_basic_basic &subs_dict, bool cache=true)
SymPy compatible subs.
Definition: subs.h:559
RCP< const Basic > conjugate(const RCP< const Basic > &arg)
Canonicalize Conjugate.
Definition: functions.cpp:149
T pow(T... args)
T push_back(T... args)
T size(T... args)
T str(T... args)
T swap(T... args)
Code