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