template<typename Scalar>
class casadi::Matrix< Scalar >
General sparse matrix class that is designed with the idea that "everything is a matrix", that is, also scalars and vectors.
This philosophy makes it easy to use and to interface in particularly with Python and Matlab/Octave.
Index starts with 0.
Index vec happens as follows: (rr, cc) -> k = rr+cc*size1()
Vectors are column vectors.
The storage format is Compressed Column Storage (CCS), similar to that used for sparse matrices in Matlab,
but unlike this format, we do allow for elements to be structurally non-zero but numerically zero.
Matrix<Scalar> is polymorphic with a std::vector<Scalar> that contain all non-identical-zero elements.
The sparsity can be accessed with Sparsity& sparsity()
- Author
- Joel Andersson
- Date
- 2010-2014
Extra doc: https://github.com/casadi/casadi/wiki/L_18e
Definition at line 88 of file matrix_decl.hpp.
|
| Matrix () |
| constructors More...
|
|
| Matrix (const Matrix< Scalar > &m) |
| Copy constructor. More...
|
|
| Matrix (casadi_int nrow, casadi_int ncol) |
| Create a sparse matrix with all structural zeros. More...
|
|
| Matrix (const Sparsity &sp) |
| Create a sparse matrix from a sparsity pattern. More...
|
|
| Matrix (const Sparsity &sp, const Matrix< Scalar > &d) |
| Construct matrix with a given sparsity and nonzeros. More...
|
|
| Matrix (double val) |
| This constructor enables implicit type conversion from a numeric type. More...
|
|
| Matrix (const std::vector< std::vector< double > > &m) |
| Dense matrix constructor with data given as vector of vectors. More...
|
|
template<typename A > |
| Matrix (const std::vector< A > &x) |
| Create an expression from a vector. More...
|
|
template<typename A > |
| Matrix (const Matrix< A > &x) |
| Create a matrix from another matrix with a different entry type. More...
|
|
bool | has_nz (casadi_int rr, casadi_int cc) const |
| Returns true if the matrix has a non-zero at location rr, cc. More...
|
|
bool | __nonzero__ () const |
| Returns the truth value of a Matrix. More...
|
|
Matrix< Scalar > | operator+ () const |
|
Matrix< Scalar > | operator- () const |
|
Matrix< Scalar > | printme (const Matrix< Scalar > &y) const |
|
Matrix< Scalar > | T () const |
| Transpose the matrix. More...
|
|
void | print_split (std::vector< std::string > &nz, std::vector< std::string > &inter) const |
| Get strings corresponding to the nonzeros and the interdependencies. More...
|
|
void | disp (std::ostream &stream, bool more=false) const |
| Print a representation of the object. More...
|
|
std::string | get_str (bool more=false) const |
| Get string representation. More...
|
|
void | print_scalar (std::ostream &stream) const |
| Print scalar. More...
|
|
void | print_vector (std::ostream &stream, bool truncate=true) const |
| Print vector-style. More...
|
|
void | print_dense (std::ostream &stream, bool truncate=true) const |
| Print dense matrix-stype. More...
|
|
void | print_sparse (std::ostream &stream, bool truncate=true) const |
| Print sparse matrix style. More...
|
|
void | clear () |
|
void | resize (casadi_int nrow, casadi_int ncol) |
|
void | reserve (casadi_int nnz) |
|
void | reserve (casadi_int nnz, casadi_int ncol) |
|
void | erase (const std::vector< casadi_int > &rr, const std::vector< casadi_int > &cc, bool ind1=false) |
| Erase a submatrix (leaving structural zeros in its place) More...
|
|
void | erase (const std::vector< casadi_int > &rr, bool ind1=false) |
| Erase a submatrix (leaving structural zeros in its place) More...
|
|
void | remove (const std::vector< casadi_int > &rr, const std::vector< casadi_int > &cc) |
| Remove columns and rows. More...
|
|
void | enlarge (casadi_int nrow, casadi_int ncol, const std::vector< casadi_int > &rr, const std::vector< casadi_int > &cc, bool ind1=false) |
| Enlarge matrix. More...
|
|
Sparsity | get_sparsity () const |
| Get an owning reference to the sparsity pattern. More...
|
|
casadi_int | element_hash () const |
| Returns a number that is unique for a given symbolic scalar. More...
|
|
bool | is_regular () const |
| Checks if expression does not contain NaN or Inf. More...
|
|
bool | is_smooth () const |
| Check if smooth. More...
|
|
bool | is_leaf () const |
| Check if SX is a leaf of the SX graph. More...
|
|
bool | is_commutative () const |
| Check whether a binary SX is commutative. More...
|
|
bool | is_symbolic () const |
| Check if symbolic (Dense) More...
|
|
bool | is_valid_input () const |
| Check if matrix can be used to define function inputs. More...
|
|
bool | is_constant () const |
| Check if the matrix is constant (note that false negative answers are possible) More...
|
|
bool | is_integer () const |
| Check if the matrix is integer-valued. More...
|
|
bool | is_zero () const |
| check if the matrix is 0 (note that false negative answers are possible) More...
|
|
bool | is_one () const |
| check if the matrix is 1 (note that false negative answers are possible) More...
|
|
bool | is_minus_one () const |
| check if the matrix is -1 (note that false negative answers are possible) More...
|
|
bool | is_eye () const |
| check if the matrix is an identity matrix (note that false negative answers More...
|
|
casadi_int | op () const |
| Get operation type. More...
|
|
bool | is_op (casadi_int op) const |
| Is it a certain operation. More...
|
|
bool | has_zeros () const |
| Check if the matrix has any zero entries which are not structural zeros. More...
|
|
std::vector< Scalar > | get_nonzeros () const |
| Get all nonzeros. More...
|
|
std::vector< Scalar > | get_elements () const |
| Get all elements. More...
|
|
| operator double () const |
| Type conversion to double. More...
|
|
| operator casadi_int () const |
| Type conversion to casadi_int. More...
|
|
std::string | name () const |
| Get name (only if symbolic scalar) More...
|
|
Matrix< Scalar > | dep (casadi_int ch=0) const |
| Get expressions of the children of the expression. More...
|
|
casadi_int | n_dep () const |
| Get the number of dependencies of a binary SXElem. More...
|
|
void | export_code (const std::string &lang, std::ostream &stream=casadi::uout(), const Dict &options=Dict()) const |
| Export matrix in specific language. More...
|
|
Dict | info () const |
|
std::string | serialize () const |
| Serialize. More...
|
|
void | serialize (SerializingStream &s) const |
| Serialize an object. More...
|
|
void | to_file (const std::string &filename, const std::string &format="") const |
|
template<typename Scalar > |
| Matrix (const std::vector< Scalar > &x) |
|
casadi_int | nnz () const |
| Get the number of (structural) non-zero elements. More...
|
|
casadi_int | nnz_lower () const |
| Get the number of non-zeros in the lower triangular half. More...
|
|
casadi_int | nnz_upper () const |
| Get the number of non-zeros in the upper triangular half. More...
|
|
casadi_int | nnz_diag () const |
| Get get the number of non-zeros on the diagonal. More...
|
|
casadi_int | numel () const |
| Get the number of elements. More...
|
|
casadi_int | size1 () const |
| Get the first dimension (i.e. number of rows) More...
|
|
casadi_int | rows () const |
| Get the number of rows, Octave-style syntax. More...
|
|
casadi_int | size2 () const |
| Get the second dimension (i.e. number of columns) More...
|
|
casadi_int | columns () const |
| Get the number of columns, Octave-style syntax. More...
|
|
std::string | dim (bool with_nz=false) const |
| Get string representation of dimensions. More...
|
|
std::pair< casadi_int, casadi_int > | size () const |
| Get the shape. More...
|
|
casadi_int | size (casadi_int axis) const |
| Get the size along a particular dimensions. More...
|
|
bool | is_empty (bool both=false) const |
| Check if the sparsity is empty, i.e. if one of the dimensions is zero. More...
|
|
bool | is_dense () const |
| Check if the matrix expression is dense. More...
|
|
bool | is_scalar (bool scalar_and_dense=false) const |
| Check if the matrix expression is scalar. More...
|
|
bool | is_square () const |
| Check if the matrix expression is square. More...
|
|
bool | is_vector () const |
| Check if the matrix is a row or column vector. More...
|
|
bool | is_row () const |
| Check if the matrix is a row vector (i.e. size1()==1) More...
|
|
bool | is_column () const |
| Check if the matrix is a column vector (i.e. size2()==1) More...
|
|
bool | is_triu () const |
| Check if the matrix is upper triangular. More...
|
|
bool | is_tril () const |
| Check if the matrix is lower triangular. More...
|
|
Sparsity | sparsity () const |
| Get the sparsity pattern. More...
|
|
|
void | get (Matrix< Scalar > &m, bool ind1, const Slice &rr) const |
|
void | get (Matrix< Scalar > &m, bool ind1, const Matrix< casadi_int > &rr) const |
|
void | get (Matrix< Scalar > &m, bool ind1, const Sparsity &sp) const |
|
|
void | get (Matrix< Scalar > &m, bool ind1, const Slice &rr, const Slice &cc) const |
|
void | get (Matrix< Scalar > &m, bool ind1, const Slice &rr, const Matrix< casadi_int > &cc) const |
|
void | get (Matrix< Scalar > &m, bool ind1, const Matrix< casadi_int > &rr, const Slice &cc) const |
|
void | get (Matrix< Scalar > &m, bool ind1, const Matrix< casadi_int > &rr, const Matrix< casadi_int > &cc) const |
|
|
void | set (const Matrix< Scalar > &m, bool ind1, const Slice &rr) |
|
void | set (const Matrix< Scalar > &m, bool ind1, const Matrix< casadi_int > &rr) |
|
void | set (const Matrix< Scalar > &m, bool ind1, const Sparsity &sp) |
|
|
void | set (const Matrix< Scalar > &m, bool ind1, const Slice &rr, const Slice &cc) |
|
void | set (const Matrix< Scalar > &m, bool ind1, const Slice &rr, const Matrix< casadi_int > &cc) |
|
void | set (const Matrix< Scalar > &m, bool ind1, const Matrix< casadi_int > &rr, const Slice &cc) |
|
void | set (const Matrix< Scalar > &m, bool ind1, const Matrix< casadi_int > &rr, const Matrix< casadi_int > &cc) |
|
|
void | get_nz (Matrix< Scalar > &m, bool ind1, const Slice &k) const |
|
void | get_nz (Matrix< Scalar > &m, bool ind1, const Matrix< casadi_int > &k) const |
|
|
void | set_nz (const Matrix< Scalar > &m, bool ind1, const Slice &k) |
|
void | set_nz (const Matrix< Scalar > &m, bool ind1, const Matrix< casadi_int > &k) |
|
|
std::vector< casadi_int > | get_row () const |
| Get the sparsity pattern. See the Sparsity class for details. More...
|
|
std::vector< casadi_int > | get_colind () const |
| Get the sparsity pattern. See the Sparsity class for details. More...
|
|
casadi_int | row (casadi_int el) const |
| Get the sparsity pattern. See the Sparsity class for details. More...
|
|
casadi_int | colind (casadi_int col) const |
| Get the sparsity pattern. See the Sparsity class for details. More...
|
|
|
static void | set_max_depth (casadi_int eq_depth=1) |
| Set or reset the depth to which equalities are being checked for simplifications. More...
|
|
static casadi_int | get_max_depth () |
| Get the depth to which equalities are being checked for simplifications. More...
|
|
static std::vector< Matrix< Scalar > > | get_input (const Function &f) |
| Get function input. More...
|
|
static std::vector< Matrix< Scalar > > | get_free (const Function &f) |
| Get free. More...
|
|
static std::string | type_name () |
| Get name of the class. More...
|
|
static Matrix< Scalar > | eye (casadi_int n) |
| create an n-by-n identity matrix More...
|
|
static void | set_precision (casadi_int precision) |
| Set the 'precision, width & scientific' used in printing and serializing to streams. More...
|
|
static void | set_width (casadi_int width) |
|
static void | set_scientific (bool scientific) |
|
static void | rng (casadi_int seed) |
| Seed the random number generator. More...
|
|
static Matrix< Scalar > | deserialize (std::istream &stream) |
| Build Sparsity from serialization. More...
|
|
static Matrix< Scalar > | deserialize (const std::string &s) |
| Build Sparsity from serialization. More...
|
|
static Matrix< Scalar > | deserialize (DeserializingStream &s) |
|
static Matrix< double > | from_file (const std::string &filename, const std::string &format_hint="") |
|
static Matrix< Scalar > | logsumexp (const Matrix< Scalar > &x) |
|
|
static Matrix< Scalar > | triplet (const std::vector< casadi_int > &row, const std::vector< casadi_int > &col, const Matrix< Scalar > &d) |
| Construct a sparse matrix from triplet form. More...
|
|
static Matrix< Scalar > | triplet (const std::vector< casadi_int > &row, const std::vector< casadi_int > &col, const Matrix< Scalar > &d, casadi_int nrow, casadi_int ncol) |
| Construct a sparse matrix from triplet form. More...
|
|
static Matrix< Scalar > | triplet (const std::vector< casadi_int > &row, const std::vector< casadi_int > &col, const Matrix< Scalar > &d, const std::pair< casadi_int, casadi_int > &rc) |
| Construct a sparse matrix from triplet form. More...
|
|
|
static Matrix< Scalar > | inf (const Sparsity &sp) |
| create a matrix with all inf More...
|
|
static Matrix< Scalar > | inf (casadi_int nrow=1, casadi_int ncol=1) |
| create a matrix with all inf More...
|
|
static Matrix< Scalar > | inf (const std::pair< casadi_int, casadi_int > &rc) |
| create a matrix with all inf More...
|
|
|
static Matrix< Scalar > | nan (const Sparsity &sp) |
| create a matrix with all nan More...
|
|
static Matrix< Scalar > | nan (casadi_int nrow=1, casadi_int ncol=1) |
| create a matrix with all nan More...
|
|
static Matrix< Scalar > | nan (const std::pair< casadi_int, casadi_int > &rc) |
| create a matrix with all nan More...
|
|
|
static Matrix< Scalar > | rand (casadi_int nrow=1, casadi_int ncol=1) |
| Create a matrix with uniformly distributed random numbers. More...
|
|
static Matrix< Scalar > | rand (const Sparsity &sp) |
| Create a matrix with uniformly distributed random numbers. More...
|
|
static Matrix< Scalar > | rand (const std::pair< casadi_int, casadi_int > &rc) |
| Create a matrix with uniformly distributed random numbers. More...
|
|
|
static Matrix< Scalar > | mpower (const Matrix< Scalar > &x, const Matrix< Scalar > &y) |
|
static Matrix< Scalar > | soc (const Matrix< Scalar > &x, const Matrix< Scalar > &y) |
|
static Matrix< Scalar > | linearize (const Matrix< Scalar > &f, const Matrix< Scalar > &x, const Matrix< Scalar > &x0, const Dict &opts=Dict()) |
|
static Matrix< Scalar > | gradient (const Matrix< Scalar > &ex, const Matrix< Scalar > &arg, const Dict &opts=Dict()) |
|
static Matrix< Scalar > | tangent (const Matrix< Scalar > &ex, const Matrix< Scalar > &arg, const Dict &opts=Dict()) |
|
static Matrix< Scalar > | jtimes (const Matrix< Scalar > &ex, const Matrix< Scalar > &arg, const Matrix< Scalar > &v, bool tr=false, const Dict &opts=Dict()) |
|
|
static Matrix< Scalar > | bilin (const Matrix< Scalar > &A, const Matrix< Scalar > &x, const Matrix< Scalar > &y) |
| Calculate bilinear/quadratic form x^T A y. More...
|
|
|
static Matrix< Scalar > | rank1 (const Matrix< Scalar > &A, const Matrix< Scalar > &alpha, const Matrix< Scalar > &x, const Matrix< Scalar > &y) |
| Make a rank-1 update to a matrix A. More...
|
|
|
The "sym" function is intended to work in a similar way as "sym" used in the Symbolic Toolbox for Matlab but instead creating a CasADi symbolic primitive.
|
static Matrix< Scalar > | sym (const std::string &name, casadi_int nrow=1, casadi_int ncol=1) |
| Create an nrow-by-ncol symbolic primitive. More...
|
|
static Matrix< Scalar > | sym (const std::string &name, const std::pair< casadi_int, casadi_int > &rc) |
| Construct a symbolic primitive with given dimensions. More...
|
|
static Matrix< Scalar > | sym (const std::string &name, const Sparsity &sp) |
| Create symbolic primitive with a given sparsity pattern. More...
|
|
static std::vector< Matrix< Scalar > > | sym (const std::string &name, const Sparsity &sp, casadi_int p) |
| Create a vector of length p with with matrices. More...
|
|
static std::vector< Matrix< Scalar > > | sym (const std::string &name, casadi_int nrow, casadi_int ncol, casadi_int p) |
| Create a vector of length p with nrow-by-ncol symbolic primitives. More...
|
|
static std::vector< std::vector< Matrix< Scalar > > > | sym (const std::string &name, const Sparsity &sp, casadi_int p, casadi_int r) |
| Create a vector of length r of vectors of length p with. More...
|
|
static std::vector< std::vector< Matrix< Scalar > > > | sym (const std::string &name, casadi_int nrow, casadi_int ncol, casadi_int p, casadi_int r) |
| Create a vector of length r of vectors of length p. More...
|
|
|
static Matrix< Scalar > | zeros (casadi_int nrow=1, casadi_int ncol=1) |
| Create a dense matrix or a matrix with specified sparsity with all entries zero. More...
|
|
static Matrix< Scalar > | zeros (const Sparsity &sp) |
| Create a dense matrix or a matrix with specified sparsity with all entries zero. More...
|
|
static Matrix< Scalar > | zeros (const std::pair< casadi_int, casadi_int > &rc) |
| Create a dense matrix or a matrix with specified sparsity with all entries zero. More...
|
|
|
static Matrix< Scalar > | ones (casadi_int nrow=1, casadi_int ncol=1) |
| Create a dense matrix or a matrix with specified sparsity with all entries one. More...
|
|
static Matrix< Scalar > | ones (const Sparsity &sp) |
| Create a dense matrix or a matrix with specified sparsity with all entries one. More...
|
|
static Matrix< Scalar > | ones (const std::pair< casadi_int, casadi_int > &rc) |
| Create a dense matrix or a matrix with specified sparsity with all entries one. More...
|
|
|
Matrix< Scalar > | adj (const Matrix< Scalar > &A) |
| Matrix adjoint. More...
|
|
Matrix< Scalar > | minor (const Matrix< Scalar > &x, casadi_int i, casadi_int j) |
| Get the (i,j) minor matrix. More...
|
|
Matrix< Scalar > | cofactor (const Matrix< Scalar > &x, casadi_int i, casadi_int j) |
| Get the (i,j) cofactor matrix. More...
|
|
void | qr (const Matrix< Scalar > &A, Matrix< Scalar > &Q, Matrix< Scalar > &R) |
| QR factorization using the modified Gram-Schmidt algorithm. More...
|
|
void | qr_sparse (const Matrix< Scalar > &A, Matrix< Scalar > &V, Matrix< Scalar > &R, Matrix< Scalar > &beta, std::vector< casadi_int > &prinv, std::vector< casadi_int > &pc, bool amd=true) |
| Sparse direct QR factorization. More...
|
|
Matrix< Scalar > | qr_solve (const Matrix< Scalar > &b, const Matrix< Scalar > &v, const Matrix< Scalar > &r, const Matrix< Scalar > &beta, const std::vector< casadi_int > &prinv, const std::vector< casadi_int > &pc, bool tr=false) |
| Solve using a sparse QR factorization. More...
|
|
Matrix< Scalar > | chol (const Matrix< Scalar > &A) |
| Obtain a Cholesky factorisation of a matrix. More...
|
|
void | ldl (const Matrix< Scalar > &A, Matrix< Scalar > &D, Matrix< Scalar > <, std::vector< casadi_int > &p, bool amd=true) |
| Sparse LDL^T factorization. More...
|
|
Matrix< Scalar > | ldl_solve (const Matrix< Scalar > &b, const Matrix< Scalar > &D, const Matrix< Scalar > <, const std::vector< casadi_int > &p) |
| Solve using a sparse LDL^T factorization. More...
|
|
Matrix< Scalar > | any (const Matrix< Scalar > &x) |
| Returns true only if any element in the matrix is true. More...
|
|
Matrix< Scalar > | all (const Matrix< Scalar > &x) |
| Returns true only if every element in the matrix is true. More...
|
|
Matrix< Scalar > | norm_inf_mul (const Matrix< Scalar > &x, const Matrix< Scalar > &y) |
| Inf-norm of a Matrix-Matrix product. More...
|
|
Matrix< Scalar > | sparsify (const Matrix< Scalar > &A, double tol=0) |
| Make a matrix sparse by removing numerical zeros. More...
|
|
void | expand (const Matrix< Scalar > &ex, Matrix< Scalar > &weights, Matrix< Scalar > &terms) |
| Expand the expression as a weighted sum (with constant weights) More...
|
|
Matrix< Scalar > | pw_const (const Matrix< Scalar > &t, const Matrix< Scalar > &tval, const Matrix< Scalar > &val) |
| Create a piecewise constant function. More...
|
|
Matrix< Scalar > | pw_lin (const Matrix< Scalar > &t, const Matrix< Scalar > &tval, const Matrix< Scalar > &val) |
| t a scalar variable (e.g. time) More...
|
|
Matrix< Scalar > | heaviside (const Matrix< Scalar > &x) |
| Heaviside function. More...
|
|
Matrix< Scalar > | rectangle (const Matrix< Scalar > &x) |
| rectangle function More...
|
|
Matrix< Scalar > | triangle (const Matrix< Scalar > &x) |
| triangle function More...
|
|
Matrix< Scalar > | ramp (const Matrix< Scalar > &x) |
| ramp function More...
|
|
Matrix< Scalar > | mtaylor (const Matrix< Scalar > &ex, const Matrix< Scalar > &x, const Matrix< Scalar > &a, casadi_int order=1) |
| multivariate Taylor series expansion More...
|
|
Matrix< Scalar > | mtaylor (const Matrix< Scalar > &ex, const Matrix< Scalar > &x, const Matrix< Scalar > &a, casadi_int order, const std::vector< casadi_int > &order_contributions) |
| multivariate Taylor series expansion More...
|
|
Matrix< Scalar > | poly_coeff (const Matrix< Scalar > &f, const Matrix< Scalar > &x) |
| extracts polynomial coefficients from an expression More...
|
|
Matrix< Scalar > | poly_roots (const Matrix< Scalar > &p) |
| Attempts to find the roots of a polynomial. More...
|
|
Matrix< Scalar > | eig_symbolic (const Matrix< Scalar > &m) |
| Attempts to find the eigenvalues of a symbolic matrix. More...
|
|
Matrix< double > | evalf (const Matrix< Scalar > &expr) |
| Evaluates the expression numerically. More...
|
|
|
Matrix< Scalar > | gauss_quadrature (const Matrix< Scalar > &f, const Matrix< Scalar > &x, const Matrix< Scalar > &a, const Matrix< Scalar > &b, casadi_int order=5) |
| Integrate f from a to b using Gaussian quadrature with n points. More...
|
|
Matrix< Scalar > | gauss_quadrature (const Matrix< Scalar > &f, const Matrix< Scalar > &x, const Matrix< Scalar > &a, const Matrix< Scalar > &b, casadi_int order, const Matrix< Scalar > &w) |
| Integrate f from a to b using Gaussian quadrature with n points. More...
|
|
|
Matrix< Scalar > | taylor (const Matrix< Scalar > &ex, const Matrix< Scalar > &x, const Matrix< Scalar > &a, casadi_int order=1) |
| univariate Taylor series expansion More...
|
|
Matrix< Scalar > | taylor (const Matrix< Scalar > &ex, const Matrix< Scalar > &x) |
| univariate Taylor series expansion More...
|
|