matrix.h
1 #ifndef SYMENGINE_MATRIX_H
2 #define SYMENGINE_MATRIX_H
3 
4 #include <symengine/basic.h>
5 #include <symengine/sets.h>
6 
7 namespace SymEngine
8 {
9 
10 enum MatrixTypeID {
11  SYMENGINE_DENSE_MATRIX,
12  SYMENGINE_CSR_MATRIX,
13 };
14 
15 // Base class for matrices
17 {
18 public:
19  virtual ~MatrixBase(){};
20 
21  virtual MatrixTypeID get_type_code() const = 0;
22 
23  bool is_square() const
24  {
25  return ncols() == nrows();
26  }
27 
28  // Below methods should be implemented by the derived classes. If not
29  // applicable, raise an exception
30 
31  // Get the # of rows and # of columns
32  virtual unsigned nrows() const = 0;
33  virtual unsigned ncols() const = 0;
34  virtual bool eq(const MatrixBase &other) const;
35 
36  virtual tribool is_real(const Assumptions *assumptions = nullptr) const = 0;
37 
38  // Get and set elements
39  virtual RCP<const Basic> get(unsigned i, unsigned j) const = 0;
40  virtual void set(unsigned i, unsigned j, const RCP<const Basic> &e) = 0;
41 
42  // Print Matrix, very mundane version, should be overriden derived
43  // class if better printing is available
44  virtual std::string __str__() const;
45 
46  virtual unsigned rank() const = 0;
47  virtual RCP<const Basic> det() const = 0;
48  virtual void inv(MatrixBase &result) const = 0;
49 
50  // Matrix addition
51  virtual void add_matrix(const MatrixBase &other,
52  MatrixBase &result) const = 0;
53 
54  // Matrix Multiplication
55  virtual void mul_matrix(const MatrixBase &other,
56  MatrixBase &result) const = 0;
57 
58  // Matrix elementwise Multiplication
59  virtual void elementwise_mul_matrix(const MatrixBase &other,
60  MatrixBase &result) const = 0;
61 
62  // Add a scalar
63  virtual void add_scalar(const RCP<const Basic> &k,
64  MatrixBase &result) const = 0;
65 
66  // Multiply by a scalar
67  virtual void mul_scalar(const RCP<const Basic> &k,
68  MatrixBase &result) const = 0;
69 
70  // Matrix conjugate
71  virtual void conjugate(MatrixBase &result) const = 0;
72 
73  // Matrix transpose
74  virtual void transpose(MatrixBase &result) const = 0;
75 
76  // Matrix conjugate transpose
77  virtual void conjugate_transpose(MatrixBase &result) const = 0;
78 
79  // Extract out a submatrix
80  virtual void submatrix(MatrixBase &result, unsigned row_start,
81  unsigned col_start, unsigned row_end,
82  unsigned col_end, unsigned row_step = 1,
83  unsigned col_step = 1) const = 0;
84  // LU factorization
85  virtual void LU(MatrixBase &L, MatrixBase &U) const = 0;
86 
87  // LDL factorization
88  virtual void LDL(MatrixBase &L, MatrixBase &D) const = 0;
89 
90  // Fraction free LU factorization
91  virtual void FFLU(MatrixBase &LU) const = 0;
92 
93  // Fraction free LDU factorization
94  virtual void FFLDU(MatrixBase &L, MatrixBase &D, MatrixBase &U) const = 0;
95 
96  // QR factorization
97  virtual void QR(MatrixBase &Q, MatrixBase &R) const = 0;
98 
99  // Cholesky decomposition
100  virtual void cholesky(MatrixBase &L) const = 0;
101 
102  // Solve Ax = b using LU factorization
103  virtual void LU_solve(const MatrixBase &b, MatrixBase &x) const = 0;
104 };
105 
107 
108 class CSRMatrix;
109 
110 // ----------------------------- Dense Matrix --------------------------------//
111 class DenseMatrix : public MatrixBase
112 {
113 public:
114  // Constructors
115  DenseMatrix();
116  DenseMatrix(const DenseMatrix &) = default;
117  DenseMatrix(unsigned row, unsigned col);
118  DenseMatrix(unsigned row, unsigned col, const vec_basic &l);
119  DenseMatrix(const vec_basic &column_elements);
120  DenseMatrix &operator=(const DenseMatrix &other) = default;
121  // type_code
122  const static MatrixTypeID type_code_id = SYMENGINE_DENSE_MATRIX;
123  virtual MatrixTypeID get_type_code() const
124  {
125  return SYMENGINE_DENSE_MATRIX;
126  }
127  // Resize
128  void resize(unsigned i, unsigned j);
129 
130  // Should implement all the virtual methods from MatrixBase
131  // and throw an exception if a method is not applicable.
132 
133  // Get and set elements
134  RCP<const Basic> get(unsigned i, unsigned j) const override;
135  void set(unsigned i, unsigned j, const RCP<const Basic> &e) override;
136  virtual vec_basic as_vec_basic() const;
137 
138  unsigned nrows() const override
139  {
140  return row_;
141  }
142  unsigned ncols() const override
143  {
144  return col_;
145  }
146 
148  std::string dumps() const;
149 
151  static DenseMatrix loads(const std::string &);
152 
153  virtual bool is_lower() const;
154  virtual bool is_upper() const;
155  virtual tribool is_zero() const;
156  virtual tribool is_diagonal() const;
157  tribool is_real(const Assumptions *assumptions = nullptr) const override;
158  virtual tribool is_symmetric() const;
159  virtual tribool is_hermitian() const;
160  virtual tribool is_weakly_diagonally_dominant() const;
161  virtual tribool is_strictly_diagonally_dominant() const;
162  virtual tribool is_positive_definite() const;
163  virtual tribool is_negative_definite() const;
164 
165  RCP<const Basic> trace() const;
166  unsigned rank() const override;
167  RCP<const Basic> det() const override;
168  void inv(MatrixBase &result) const override;
169 
170  // Matrix addition
171  void add_matrix(const MatrixBase &other, MatrixBase &result) const override;
172 
173  // Matrix multiplication
174  void mul_matrix(const MatrixBase &other, MatrixBase &result) const override;
175 
176  // Matrix elementwise Multiplication
177  void elementwise_mul_matrix(const MatrixBase &other,
178  MatrixBase &result) const override;
179 
180  // Add a scalar
181  void add_scalar(const RCP<const Basic> &k,
182  MatrixBase &result) const override;
183 
184  // Multiply by a scalar
185  void mul_scalar(const RCP<const Basic> &k,
186  MatrixBase &result) const override;
187 
188  // Matrix conjugate
189  void conjugate(MatrixBase &result) const override;
190 
191  // Matrix transpose
192  void transpose(MatrixBase &result) const override;
193 
194  // Matrix conjugate transpose
195  void conjugate_transpose(MatrixBase &result) const override;
196 
197  // Extract out a submatrix
198  void submatrix(MatrixBase &result, unsigned row_start, unsigned col_start,
199  unsigned row_end, unsigned col_end, unsigned row_step = 1,
200  unsigned col_step = 1) const override;
201 
202  // LU factorization
203  void LU(MatrixBase &L, MatrixBase &U) const override;
204 
205  // LDL factorization
206  void LDL(MatrixBase &L, MatrixBase &D) const override;
207 
208  // Solve Ax = b using LU factorization
209  void LU_solve(const MatrixBase &b, MatrixBase &x) const override;
210 
211  // Fraction free LU factorization
212  void FFLU(MatrixBase &LU) const override;
213 
214  // Fraction free LDU factorization
215  void FFLDU(MatrixBase &L, MatrixBase &D, MatrixBase &U) const override;
216 
217  // QR factorization
218  void QR(MatrixBase &Q, MatrixBase &R) const override;
219 
220  // Cholesky decomposition
221  void cholesky(MatrixBase &L) const override;
222 
223  // Return the Jacobian of the matrix
224  friend void jacobian(const DenseMatrix &A, const DenseMatrix &x,
225  DenseMatrix &result, bool diff_cache);
226  // Return the Jacobian of the matrix using sdiff
227  friend void sjacobian(const DenseMatrix &A, const DenseMatrix &x,
228  DenseMatrix &result, bool diff_cache);
229 
230  // Differentiate the matrix element-wise
231  friend void diff(const DenseMatrix &A, const RCP<const Symbol> &x,
232  DenseMatrix &result, bool diff_cache);
233  // Differentiate the matrix element-wise using SymPy compatible diff
234  friend void sdiff(const DenseMatrix &A, const RCP<const Basic> &x,
235  DenseMatrix &result, bool diff_cache);
236 
237  // Friend functions related to Matrix Operations
238  friend void add_dense_dense(const DenseMatrix &A, const DenseMatrix &B,
239  DenseMatrix &C);
240  friend void add_dense_scalar(const DenseMatrix &A,
241  const RCP<const Basic> &k, DenseMatrix &B);
242  friend void mul_dense_dense(const DenseMatrix &A, const DenseMatrix &B,
243  DenseMatrix &C);
244  friend void elementwise_mul_dense_dense(const DenseMatrix &A,
245  const DenseMatrix &B,
246  DenseMatrix &C);
247  friend void mul_dense_scalar(const DenseMatrix &A,
248  const RCP<const Basic> &k, DenseMatrix &C);
249  friend void conjugate_dense(const DenseMatrix &A, DenseMatrix &B);
250  friend void transpose_dense(const DenseMatrix &A, DenseMatrix &B);
251  friend void conjugate_transpose_dense(const DenseMatrix &A, DenseMatrix &B);
252  friend void submatrix_dense(const DenseMatrix &A, DenseMatrix &B,
253  unsigned row_start, unsigned col_start,
254  unsigned row_end, unsigned col_end,
255  unsigned row_step, unsigned col_step);
256  void row_join(const DenseMatrix &B);
257  void col_join(const DenseMatrix &B);
258  void row_insert(const DenseMatrix &B, unsigned pos);
259  void col_insert(const DenseMatrix &B, unsigned pos);
260  void row_del(unsigned k);
261  void col_del(unsigned k);
262 
263  // Row operations
264  friend void row_exchange_dense(DenseMatrix &A, unsigned i, unsigned j);
265  friend void row_mul_scalar_dense(DenseMatrix &A, unsigned i,
266  RCP<const Basic> &c);
267  friend void row_add_row_dense(DenseMatrix &A, unsigned i, unsigned j,
268  RCP<const Basic> &c);
269  friend void permuteFwd(DenseMatrix &A, permutelist &pl);
270 
271  // Column operations
272  friend void column_exchange_dense(DenseMatrix &A, unsigned i, unsigned j);
273 
274  // Gaussian elimination
275  friend void pivoted_gaussian_elimination(const DenseMatrix &A,
276  DenseMatrix &B,
277  permutelist &pivotlist);
278  friend void fraction_free_gaussian_elimination(const DenseMatrix &A,
279  DenseMatrix &B);
280  friend void pivoted_fraction_free_gaussian_elimination(
281  const DenseMatrix &A, DenseMatrix &B, permutelist &pivotlist);
282  friend void pivoted_gauss_jordan_elimination(const DenseMatrix &A,
283  DenseMatrix &B,
284  permutelist &pivotlist);
285  friend void fraction_free_gauss_jordan_elimination(const DenseMatrix &A,
286  DenseMatrix &B);
287  friend void pivoted_fraction_free_gauss_jordan_elimination(
288  const DenseMatrix &A, DenseMatrix &B, permutelist &pivotlist);
289  friend unsigned pivot(DenseMatrix &B, unsigned r, unsigned c);
290 
291  friend void reduced_row_echelon_form(const DenseMatrix &A, DenseMatrix &B,
292  vec_uint &pivot_cols,
293  bool normalize_last);
294 
295  // Ax = b
296  friend void diagonal_solve(const DenseMatrix &A, const DenseMatrix &b,
297  DenseMatrix &x);
298  friend void back_substitution(const DenseMatrix &U, const DenseMatrix &b,
299  DenseMatrix &x);
300  friend void forward_substitution(const DenseMatrix &A, const DenseMatrix &b,
301  DenseMatrix &x);
302  friend void fraction_free_gaussian_elimination_solve(const DenseMatrix &A,
303  const DenseMatrix &b,
304  DenseMatrix &x);
305  friend void fraction_free_gauss_jordan_solve(const DenseMatrix &A,
306  const DenseMatrix &b,
307  DenseMatrix &x, bool pivot);
308 
309  // Matrix Decomposition
310  friend void fraction_free_LU(const DenseMatrix &A, DenseMatrix &LU);
311  friend void LU(const DenseMatrix &A, DenseMatrix &L, DenseMatrix &U);
312  friend void pivoted_LU(const DenseMatrix &A, DenseMatrix &LU,
313  permutelist &pl);
314  friend void pivoted_LU(const DenseMatrix &A, DenseMatrix &L, DenseMatrix &U,
315  permutelist &pl);
316  friend void fraction_free_LDU(const DenseMatrix &A, DenseMatrix &L,
317  DenseMatrix &D, DenseMatrix &U);
318  friend void QR(const DenseMatrix &A, DenseMatrix &Q, DenseMatrix &R);
319  friend void LDL(const DenseMatrix &A, DenseMatrix &L, DenseMatrix &D);
320  friend void cholesky(const DenseMatrix &A, DenseMatrix &L);
321 
322  // Matrix queries
323  friend bool is_symmetric_dense(const DenseMatrix &A);
324 
325  // Determinant
326  friend RCP<const Basic> det_bareis(const DenseMatrix &A);
327  friend void berkowitz(const DenseMatrix &A,
328  std::vector<DenseMatrix> &polys);
329 
330  // Inverse
331  friend void inverse_fraction_free_LU(const DenseMatrix &A, DenseMatrix &B);
332  friend void inverse_LU(const DenseMatrix &A, DenseMatrix &B);
333  friend void inverse_pivoted_LU(const DenseMatrix &A, DenseMatrix &B);
334  friend void inverse_gauss_jordan(const DenseMatrix &A, DenseMatrix &B);
335 
336  // Vector-specific methods
337  friend void dot(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &C);
338  friend void cross(const DenseMatrix &A, const DenseMatrix &B,
339  DenseMatrix &C);
340 
341  // NumPy-like functions
342  friend void eye(DenseMatrix &A, int k);
343  friend void diag(DenseMatrix &A, vec_basic &v, int k);
344  friend void ones(DenseMatrix &A);
345  friend void zeros(DenseMatrix &A);
346 
347  friend CSRMatrix;
348 
349 private:
350  // Matrix elements are stored in row-major order
351  vec_basic m_;
352  // Stores the dimension of the Matrix
353  unsigned row_;
354  unsigned col_;
355 
356  tribool shortcut_to_posdef() const;
357  tribool is_positive_definite_GE();
358 };
359 
360 // ----------------------------- Sparse Matrices -----------------------------//
361 class CSRMatrix : public MatrixBase
362 {
363 public:
364  CSRMatrix();
365  CSRMatrix(unsigned row, unsigned col);
366  CSRMatrix(unsigned row, unsigned col, const std::vector<unsigned> &p,
367  const std::vector<unsigned> &j, const vec_basic &x);
368  CSRMatrix(unsigned row, unsigned col, std::vector<unsigned> &&p,
370  CSRMatrix &operator=(CSRMatrix &&other);
371  CSRMatrix(const CSRMatrix &) = default;
372  const static MatrixTypeID type_code_id = SYMENGINE_CSR_MATRIX;
373  virtual MatrixTypeID get_type_code() const
374  {
375  return SYMENGINE_CSR_MATRIX;
376  }
378  as_vectors() const;
379 
380  bool is_canonical() const;
381 
382  bool eq(const MatrixBase &other) const override;
383 
384  // Get and set elements
385  RCP<const Basic> get(unsigned i, unsigned j) const override;
386  void set(unsigned i, unsigned j, const RCP<const Basic> &e) override;
387 
388  unsigned nrows() const override
389  {
390  return row_;
391  }
392  unsigned ncols() const override
393  {
394  return col_;
395  }
396 
397  tribool is_real(const Assumptions *assumptions = nullptr) const override;
398  unsigned rank() const override;
399  RCP<const Basic> det() const override;
400  void inv(MatrixBase &result) const override;
401 
402  // Matrix addition
403  void add_matrix(const MatrixBase &other, MatrixBase &result) const override;
404 
405  // Matrix Multiplication
406  void mul_matrix(const MatrixBase &other, MatrixBase &result) const override;
407 
408  // Matrix elementwise Multiplication
409  void elementwise_mul_matrix(const MatrixBase &other,
410  MatrixBase &result) const override;
411 
412  // Add a scalar
413  void add_scalar(const RCP<const Basic> &k,
414  MatrixBase &result) const override;
415 
416  // Multiply by a scalar
417  void mul_scalar(const RCP<const Basic> &k,
418  MatrixBase &result) const override;
419 
420  // Matrix conjugate
421  void conjugate(MatrixBase &result) const override;
422 
423  // Matrix transpose
424  void transpose(MatrixBase &result) const override;
425  CSRMatrix transpose(bool conjugate = false) const;
426 
427  // Matrix conjugate transpose
428  void conjugate_transpose(MatrixBase &result) const override;
429 
430  // Extract out a submatrix
431  void submatrix(MatrixBase &result, unsigned row_start, unsigned col_start,
432  unsigned row_end, unsigned col_end, unsigned row_step = 1,
433  unsigned col_step = 1) const override;
434 
435  // LU factorization
436  void LU(MatrixBase &L, MatrixBase &U) const override;
437 
438  // LDL factorization
439  void LDL(MatrixBase &L, MatrixBase &D) const override;
440 
441  // Solve Ax = b using LU factorization
442  void LU_solve(const MatrixBase &b, MatrixBase &x) const override;
443 
444  // Fraction free LU factorization
445  void FFLU(MatrixBase &LU) const override;
446 
447  // Fraction free LDU factorization
448  void FFLDU(MatrixBase &L, MatrixBase &D, MatrixBase &U) const override;
449 
450  // QR factorization
451  void QR(MatrixBase &Q, MatrixBase &R) const override;
452 
453  // Cholesky decomposition
454  void cholesky(MatrixBase &L) const override;
455 
456  static void csr_sum_duplicates(std::vector<unsigned> &p_,
458  unsigned row_);
459 
460  static void csr_sort_indices(std::vector<unsigned> &p_,
462  unsigned row_);
463 
464  static bool csr_has_sorted_indices(const std::vector<unsigned> &p_,
465  const std::vector<unsigned> &j_,
466  unsigned row_);
467 
468  static bool csr_has_duplicates(const std::vector<unsigned> &p_,
469  const std::vector<unsigned> &j_,
470  unsigned row_);
471 
472  static bool csr_has_canonical_format(const std::vector<unsigned> &p_,
473  const std::vector<unsigned> &j_,
474  unsigned row_);
475 
476  static CSRMatrix from_coo(unsigned row, unsigned col,
477  const std::vector<unsigned> &i,
478  const std::vector<unsigned> &j,
479  const vec_basic &x);
480  static CSRMatrix jacobian(const vec_basic &exprs, const vec_sym &x,
481  bool diff_cache = true);
482  static CSRMatrix jacobian(const DenseMatrix &A, const DenseMatrix &x,
483  bool diff_cache = true);
484 
485  friend void csr_matmat_pass1(const CSRMatrix &A, const CSRMatrix &B,
486  CSRMatrix &C);
487  friend void csr_matmat_pass2(const CSRMatrix &A, const CSRMatrix &B,
488  CSRMatrix &C);
489  friend void csr_diagonal(const CSRMatrix &A, DenseMatrix &D);
490  friend void csr_scale_rows(CSRMatrix &A, const DenseMatrix &X);
491  friend void csr_scale_columns(CSRMatrix &A, const DenseMatrix &X);
492 
493  friend void csr_binop_csr_canonical(
494  const CSRMatrix &A, const CSRMatrix &B, CSRMatrix &C,
495  RCP<const Basic> (&bin_op)(const RCP<const Basic> &,
496  const RCP<const Basic> &));
497 
498 private:
501  vec_basic x_;
502  // Stores the dimension of the Matrix
503  unsigned row_;
504  unsigned col_;
505 };
506 
507 // Return the Jacobian of the matrix
508 void jacobian(const DenseMatrix &A, const DenseMatrix &x, DenseMatrix &result,
509  bool diff_cache = true);
510 // Return the Jacobian of the matrix using sdiff
511 void sjacobian(const DenseMatrix &A, const DenseMatrix &x, DenseMatrix &result,
512  bool diff_cache = true);
513 
514 // Differentiate all the elements
515 void diff(const DenseMatrix &A, const RCP<const Symbol> &x, DenseMatrix &result,
516  bool diff_cache = true);
517 // Differentiate all the elements using SymPy compatible diff
518 void sdiff(const DenseMatrix &A, const RCP<const Basic> &x, DenseMatrix &result,
519  bool diff_cache = true);
520 
521 // Get submatrix from a DenseMatrix
522 void submatrix_dense(const DenseMatrix &A, DenseMatrix &B, unsigned row_start,
523  unsigned col_start, unsigned row_end, unsigned col_end,
524  unsigned row_step = 1, unsigned col_step = 1);
525 
526 // Row operations
527 void row_exchange_dense(DenseMatrix &A, unsigned i, unsigned j);
528 void row_mul_scalar_dense(DenseMatrix &A, unsigned i, RCP<const Basic> &c);
529 void row_add_row_dense(DenseMatrix &A, unsigned i, unsigned j,
530  RCP<const Basic> &c);
531 
532 // Column operations
533 void column_exchange_dense(DenseMatrix &A, unsigned i, unsigned j);
534 
535 // Vector-specific methods
536 void dot(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &C);
537 void cross(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &C);
538 
539 // Matrix Factorization
540 void LU(const DenseMatrix &A, DenseMatrix &L, DenseMatrix &U);
541 void LDL(const DenseMatrix &A, DenseMatrix &L, DenseMatrix &D);
542 void QR(const DenseMatrix &A, DenseMatrix &Q, DenseMatrix &R);
543 void cholesky(const DenseMatrix &A, DenseMatrix &L);
544 
545 // Inverse
546 void inverse_fraction_free_LU(const DenseMatrix &A, DenseMatrix &B);
547 
548 void inverse_gauss_jordan(const DenseMatrix &A, DenseMatrix &B);
549 
550 // Solving Ax = b
551 void fraction_free_LU_solve(const DenseMatrix &A, const DenseMatrix &b,
552  DenseMatrix &x);
553 
554 void fraction_free_gauss_jordan_solve(const DenseMatrix &A,
555  const DenseMatrix &b, DenseMatrix &x,
556  bool pivot = true);
557 
558 void LU_solve(const DenseMatrix &A, const DenseMatrix &b, DenseMatrix &x);
559 void pivoted_LU_solve(const DenseMatrix &A, const DenseMatrix &b,
560  DenseMatrix &x);
561 
562 void LDL_solve(const DenseMatrix &A, const DenseMatrix &b, DenseMatrix &x);
563 
564 // Determinant
565 RCP<const Basic> det_berkowitz(const DenseMatrix &A);
566 
567 // Characteristic polynomial: Only the coefficients of monomials in decreasing
568 // order of monomial powers is returned, i.e. if `B = transpose([1, -2, 3])`
569 // then the corresponding polynomial is `x**2 - 2x + 3`.
570 void char_poly(const DenseMatrix &A, DenseMatrix &B);
571 
572 // returns a finiteset of eigenvalues of a matrix
573 RCP<const Set> eigen_values(const DenseMatrix &A);
574 
575 // Mimic `eye` function in NumPy
576 void eye(DenseMatrix &A, int k = 0);
577 
578 // Create diagonal matrices directly
579 void diag(DenseMatrix &A, vec_basic &v, int k = 0);
580 
581 // Create a matrix filled with ones
582 void ones(DenseMatrix &A);
583 
584 // Create a matrix filled with zeros
585 void zeros(DenseMatrix &A);
586 
587 // Reduced row echelon form and returns the cols with pivots
588 void reduced_row_echelon_form(const DenseMatrix &A, DenseMatrix &B,
589  vec_uint &pivot_cols,
590  bool normalize_last = false);
591 
592 // Returns true if `b` is exactly the type T.
593 // Here T can be a DenseMatrix, CSRMatrix, etc.
594 template <class T>
595 inline bool is_a(const MatrixBase &b)
596 {
597  return T::type_code_id == b.get_type_code();
598 }
599 
600 // Test two matrices for equality
601 inline bool operator==(const SymEngine::MatrixBase &lhs,
602  const SymEngine::MatrixBase &rhs)
603 {
604  return lhs.eq(rhs);
605 }
606 
607 // Test two matrices for equality
608 inline bool operator!=(const SymEngine::MatrixBase &lhs,
609  const SymEngine::MatrixBase &rhs)
610 {
611  return not lhs.eq(rhs);
612 }
613 
614 } // namespace SymEngine
615 
616 // Print Matrix
618  const SymEngine::MatrixBase &A)
619 {
620  return out << A.__str__();
621 }
622 
623 #endif
The base class for SymEngine.
static DenseMatrix loads(const std::string &)
Creates an instance of a serialized string.
std::string dumps() const
Returns a string of the instance serialized.
Main namespace for SymEngine package.
Definition: add.cpp:19
bool is_a(const Basic &b)
Templatised version to check is_a type.
Definition: basic-inl.h:36
std::ostream & operator<<(std::ostream &out, const SymEngine::Basic &p)
<< Operator
Definition: basic-inl.h:55
RCP< const Basic > conjugate(const RCP< const Basic > &arg)
Canonicalize Conjugate.
Definition: functions.cpp:149
T operator!=(T... args)
Code