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