Class DenseMatrix

Inheritance Relationships

Base Type

Class Documentation

class SymEngine::DenseMatrix : public SymEngine::MatrixBase

Public Functions

DenseMatrix()
DenseMatrix(const DenseMatrix&) = default
DenseMatrix(unsigned row, unsigned col)
DenseMatrix(unsigned row, unsigned col, const vec_basic &l)
DenseMatrix(const vec_basic &column_elements)
DenseMatrix &operator=(const DenseMatrix &other) = default
void resize(unsigned i, unsigned j)
RCP<const Basic> get(unsigned i, unsigned j) const
void set(unsigned i, unsigned j, const RCP<const Basic> &e)
vec_basic as_vec_basic() const
unsigned nrows() const
unsigned ncols() const
bool is_lower() const
bool is_upper() const
unsigned rank() const
RCP<const Basic> det() const
void inv(MatrixBase &result) const
void add_matrix(const MatrixBase &other, MatrixBase &result) const
void mul_matrix(const MatrixBase &other, MatrixBase &result) const
void elementwise_mul_matrix(const MatrixBase &other, MatrixBase &result) const
void add_scalar(const RCP<const Basic> &k, MatrixBase &result) const
void mul_scalar(const RCP<const Basic> &k, MatrixBase &result) const
void conjugate(MatrixBase &result) const
void transpose(MatrixBase &result) const
void conjugate_transpose(MatrixBase &result) const
void submatrix(MatrixBase &result, unsigned row_start, unsigned col_start, unsigned row_end, unsigned col_end, unsigned row_step = 1, unsigned col_step = 1) const
void LU(MatrixBase &L, MatrixBase &U) const
void LDL(MatrixBase &L, MatrixBase &D) const
void LU_solve(const MatrixBase &b, MatrixBase &x) const
void FFLU(MatrixBase &LU) const
void FFLDU(MatrixBase &L, MatrixBase &D, MatrixBase &U) const
void QR(MatrixBase &Q, MatrixBase &R) const
void cholesky(MatrixBase &L) const
void row_join(const DenseMatrix &B)
void col_join(const DenseMatrix &B)
void row_insert(const DenseMatrix &B, unsigned pos)
void col_insert(const DenseMatrix &B, unsigned pos)
void row_del(unsigned k)
void col_del(unsigned k)

Public Members

friend CSRMatrix

Friends

friend void jacobian(const DenseMatrix &A, const DenseMatrix &x, DenseMatrix &result, bool diff_cache)
friend void sjacobian(const DenseMatrix &A, const DenseMatrix &x, DenseMatrix &result, bool diff_cache)
friend void diff(const DenseMatrix &A, const RCP<const Symbol> &x, DenseMatrix &result, bool diff_cache)
friend void sdiff(const DenseMatrix &A, const RCP<const Basic> &x, DenseMatrix &result, bool diff_cache)
friend void add_dense_dense(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &C)
friend void add_dense_scalar(const DenseMatrix &A, const RCP<const Basic> &k, DenseMatrix &B)
friend void mul_dense_dense(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &C)
friend void elementwise_mul_dense_dense(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &C)
friend void mul_dense_scalar(const DenseMatrix &A, const RCP<const Basic> &k, DenseMatrix &C)
friend void conjugate_dense(const DenseMatrix &A, DenseMatrix &B)
friend void transpose_dense(const DenseMatrix &A, DenseMatrix &B)
friend void conjugate_transpose_dense(const DenseMatrix &A, DenseMatrix &B)
friend void submatrix_dense(const DenseMatrix &A, DenseMatrix &B, unsigned row_start, unsigned col_start, unsigned row_end, unsigned col_end, unsigned row_step, unsigned col_step)
friend void row_exchange_dense(DenseMatrix &A, unsigned i, unsigned j)
friend void row_mul_scalar_dense(DenseMatrix &A, unsigned i, RCP<const Basic> &c)
friend void row_add_row_dense(DenseMatrix &A, unsigned i, unsigned j, RCP<const Basic> &c)
friend void permuteFwd(DenseMatrix &A, permutelist &pl)
friend void column_exchange_dense(DenseMatrix &A, unsigned i, unsigned j)
friend void pivoted_gaussian_elimination(const DenseMatrix &A, DenseMatrix &B, permutelist &pivotlist)
friend void fraction_free_gaussian_elimination(const DenseMatrix &A, DenseMatrix &B)
friend void pivoted_fraction_free_gaussian_elimination(const DenseMatrix &A, DenseMatrix &B, permutelist &pivotlist)
friend void pivoted_gauss_jordan_elimination(const DenseMatrix &A, DenseMatrix &B, permutelist &pivotlist)
friend void fraction_free_gauss_jordan_elimination(const DenseMatrix &A, DenseMatrix &B)
friend void pivoted_fraction_free_gauss_jordan_elimination(const DenseMatrix &A, DenseMatrix &B, permutelist &pivotlist)
friend unsigned pivot(DenseMatrix &B, unsigned r, unsigned c)
friend void reduced_row_echelon_form(const DenseMatrix &A, DenseMatrix &B, vec_uint &pivot_cols, bool normalize_last)
friend void diagonal_solve(const DenseMatrix &A, const DenseMatrix &b, DenseMatrix &x)
friend void back_substitution(const DenseMatrix &U, const DenseMatrix &b, DenseMatrix &x)
friend void forward_substitution(const DenseMatrix &A, const DenseMatrix &b, DenseMatrix &x)
friend void fraction_free_gaussian_elimination_solve(const DenseMatrix &A, const DenseMatrix &b, DenseMatrix &x)
friend void fraction_free_gauss_jordan_solve(const DenseMatrix &A, const DenseMatrix &b, DenseMatrix &x)
friend void fraction_free_LU(const DenseMatrix &A, DenseMatrix &LU)
friend void LU(const DenseMatrix &A, DenseMatrix &L, DenseMatrix &U)
friend void pivoted_LU(const DenseMatrix &A, DenseMatrix &LU, permutelist &pl)
friend void pivoted_LU(const DenseMatrix &A, DenseMatrix &L, DenseMatrix &U, permutelist &pl)
friend void fraction_free_LDU(const DenseMatrix &A, DenseMatrix &L, DenseMatrix &D, DenseMatrix &U)
friend void QR(const DenseMatrix &A, DenseMatrix &Q, DenseMatrix &R)
friend void LDL(const DenseMatrix &A, DenseMatrix &L, DenseMatrix &D)
friend void cholesky(const DenseMatrix &A, DenseMatrix &L)
friend bool is_symmetric_dense(const DenseMatrix &A)
friend RCP<const Basic> det_bareis(const DenseMatrix &A)
friend void berkowitz(const DenseMatrix &A, std::vector<DenseMatrix> &polys)
friend void inverse_fraction_free_LU(const DenseMatrix &A, DenseMatrix &B)
friend void inverse_LU(const DenseMatrix &A, DenseMatrix &B)
friend void inverse_pivoted_LU(const DenseMatrix &A, DenseMatrix &B)
friend void inverse_gauss_jordan(const DenseMatrix &A, DenseMatrix &B)
friend void dot(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &C)
friend void cross(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &C)
friend void eye(DenseMatrix &A, int k)
friend void diag(DenseMatrix &A, vec_basic &v, int k)
friend void ones(DenseMatrix &A)
friend void zeros(DenseMatrix &A)