List of all members | Public Member Functions | Static Public Member Functions | Friends
casadi::Matrix< Scalar > Class Template Reference

Sparse matrix class. SX and DM are specializations. More...

#include <matrix_decl.hpp>

Detailed Description

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.

Inheritance diagram for casadi::Matrix< Scalar >:
Inheritance graph
[legend]

Public Member Functions

 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 Public Member Functions

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...
 
Construct symbolic primitives

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...
 

Friends

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 > &LT, 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 > &LT, 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...
 

Constructor & Destructor Documentation

◆ Matrix() [1/10]

template<typename Scalar >
casadi::Matrix< Scalar >::Matrix

Extra doc: https://github.com/casadi/casadi/wiki/L_18f empty 0-by-0 matrix constructor

Definition at line 576 of file matrix_impl.hpp.

◆ Matrix() [2/10]

template<typename Scalar >
casadi::Matrix< Scalar >::Matrix ( const Matrix< Scalar > &  m)

Definition at line 580 of file matrix_impl.hpp.

◆ Matrix() [3/10]

template<typename Scalar >
casadi::Matrix< Scalar >::Matrix ( casadi_int  nrow,
casadi_int  ncol 
)

◆ Matrix() [4/10]

template<typename Scalar >
casadi::Matrix< Scalar >::Matrix ( const Sparsity sp)
explicit

Same as Matrix::ones(sparsity)

Extra doc: https://github.com/casadi/casadi/wiki/L_18k

◆ Matrix() [5/10]

template<typename Scalar >
casadi::Matrix< Scalar >::Matrix ( const Sparsity sp,
const Matrix< Scalar > &  d 
)

◆ Matrix() [6/10]

template<typename Scalar >
casadi::Matrix< Scalar >::Matrix ( double  val)

◆ Matrix() [7/10]

template<typename Scalar >
casadi::Matrix< Scalar >::Matrix ( const std::vector< std::vector< double > > &  m)
explicit

◆ Matrix() [8/10]

template<typename Scalar >
template<typename A >
casadi::Matrix< Scalar >::Matrix ( const std::vector< A > &  x)
inline

Extra doc: https://github.com/casadi/casadi/wiki/L_18m

Definition at line 154 of file matrix_decl.hpp.

◆ Matrix() [9/10]

template<typename Scalar >
template<typename A >
casadi::Matrix< Scalar >::Matrix ( const Matrix< A > &  x)
inline

Assumes that the scalar conversion is valid.

Extra doc: https://github.com/casadi/casadi/wiki/L_18n

Definition at line 167 of file matrix_decl.hpp.

◆ Matrix() [10/10]

template<typename Scalar >
template<typename Scalar >
casadi::Matrix< Scalar >::Matrix ( const std::vector< Scalar > &  x)

Definition at line 584 of file matrix_impl.hpp.

Member Function Documentation

◆ __nonzero__()

bool casadi::SX::__nonzero__ ( ) const

Definition at line 70 of file matrix_impl.hpp.

◆ bilin()

static Matrix< Scalar > casadi::GenericMatrix< Matrix< Scalar > >::bilin ( const Matrix< Scalar > &  A,
const Matrix< Scalar > &  x,
const Matrix< Scalar > &  y 
)
staticinherited
Parameters
[in]ycan be omitted, in which case x^T A x is calculated

Extra doc: https://github.com/casadi/casadi/wiki/L_1bo

◆ clear()

template<typename Scalar >
void casadi::Matrix< Scalar >::clear ( )

◆ colind()

casadi_int casadi::GenericMatrix< Matrix< Scalar > >::colind ( casadi_int  col) const
inlineinherited

Extra doc: https://github.com/casadi/casadi/wiki/L_1b8

Definition at line 201 of file generic_matrix.hpp.

◆ columns()

casadi_int casadi::GenericMatrix< Matrix< Scalar > >::columns ( ) const
inlineinherited

Extra doc: https://github.com/casadi/casadi/wiki/L_1av

Definition at line 124 of file generic_matrix.hpp.

◆ dep()

SX casadi::SX::dep ( casadi_int  ch = 0) const

Only defined if symbolic scalar. Wraps SXElem SXElem::dep(casadi_int ch=0) const.

Extra doc: https://github.com/casadi/casadi/wiki/L_1a9

◆ deserialize() [1/3]

template<typename Scalar >
static Matrix<Scalar> casadi::Matrix< Scalar >::deserialize ( const std::string &  s)
static

◆ deserialize() [2/3]

template<typename Scalar >
static Matrix<Scalar> casadi::Matrix< Scalar >::deserialize ( DeserializingStream s)
static

◆ deserialize() [3/3]

template<typename Scalar >
static Matrix<Scalar> casadi::Matrix< Scalar >::deserialize ( std::istream &  stream)
static

◆ dim()

std::string casadi::GenericMatrix< Matrix< Scalar > >::dim ( bool  with_nz = false) const
inherited

The representation is e.g. "4x5" or "4x5,10nz"

Extra doc: https://github.com/casadi/casadi/wiki/L_1aw

◆ disp()

template<typename Scalar >
void casadi::Matrix< Scalar >::disp ( std::ostream &  stream,
bool  more = false 
) const

◆ element_hash()

casadi_int casadi::SX::element_hash ( ) const

Only defined if symbolic scalar.

Extra doc: https://github.com/casadi/casadi/wiki/L_19n

◆ enlarge()

template<typename Scalar >
void casadi::Matrix< Scalar >::enlarge ( casadi_int  nrow,
casadi_int  ncol,
const std::vector< casadi_int > &  rr,
const std::vector< casadi_int > &  cc,
bool  ind1 = false 
)

Make the matrix larger by inserting empty rows and columns, keeping the existing non-zeros

Extra doc: https://github.com/casadi/casadi/wiki/L_19i

◆ erase() [1/2]

template<typename Scalar >
void casadi::Matrix< Scalar >::erase ( const std::vector< casadi_int > &  rr,
bool  ind1 = false 
)

Erase elements of a matrix

Extra doc: https://github.com/casadi/casadi/wiki/L_19g

◆ erase() [2/2]

template<typename Scalar >
void casadi::Matrix< Scalar >::erase ( const std::vector< casadi_int > &  rr,
const std::vector< casadi_int > &  cc,
bool  ind1 = false 
)

Erase rows and/or columns of a matrix

Extra doc: https://github.com/casadi/casadi/wiki/L_19f

◆ export_code()

void casadi::DM::export_code ( const std::string &  lang,
std::ostream &  stream = casadi::uout(),
const Dict options = Dict() 
) const

lang: only 'matlab' supported for now

* options:
*   inline: Indicates if you want everything on a single line (default: False)
*   name: Name of exported variable (default: 'm')
*   indent_level: Level of indentation (default: 0)
*   spoof_zero: Replace numerical zero by a 1e-200 (default: false)
*               might be needed for matlab sparse construct,
*               which doesn't allow numerical zero
* 

Extra doc: https://github.com/casadi/casadi/wiki/L_1ac

◆ eye()

template<typename Scalar >
static Matrix<Scalar> casadi::Matrix< Scalar >::eye ( casadi_int  n)
static

◆ from_file()

DM casadi::DM::from_file ( const std::string &  filename,
const std::string &  format_hint = "" 
)
static

◆ get() [1/7]

template<typename Scalar >
void casadi::Matrix< Scalar >::get ( Matrix< Scalar > &  m,
bool  ind1,
const Matrix< casadi_int > &  rr 
) const

Get a submatrix, single argument

Definition at line 152 of file matrix_impl.hpp.

◆ get() [2/7]

template<typename Scalar >
void casadi::Matrix< Scalar >::get ( Matrix< Scalar > &  m,
bool  ind1,
const Matrix< casadi_int > &  rr,
const Matrix< casadi_int > &  cc 
) const

Get a submatrix, two arguments

Definition at line 111 of file matrix_impl.hpp.

◆ get() [3/7]

template<typename Scalar >
void casadi::Matrix< Scalar >::get ( Matrix< Scalar > &  m,
bool  ind1,
const Matrix< casadi_int > &  rr,
const Slice cc 
) const

Get a submatrix, two arguments

Definition at line 104 of file matrix_impl.hpp.

◆ get() [4/7]

template<typename Scalar >
void casadi::Matrix< Scalar >::get ( Matrix< Scalar > &  m,
bool  ind1,
const Slice rr 
) const

Get a submatrix, single argument

Definition at line 134 of file matrix_impl.hpp.

◆ get() [5/7]

template<typename Scalar >
void casadi::Matrix< Scalar >::get ( Matrix< Scalar > &  m,
bool  ind1,
const Slice rr,
const Matrix< casadi_int > &  cc 
) const

Get a submatrix, two arguments

Definition at line 97 of file matrix_impl.hpp.

◆ get() [6/7]

template<typename Scalar >
void casadi::Matrix< Scalar >::get ( Matrix< Scalar > &  m,
bool  ind1,
const Slice rr,
const Slice cc 
) const

Get a submatrix, two arguments

Definition at line 79 of file matrix_impl.hpp.

◆ get() [7/7]

template<typename Scalar >
void casadi::Matrix< Scalar >::get ( Matrix< Scalar > &  m,
bool  ind1,
const Sparsity sp 
) const

Get a submatrix, single argument

Definition at line 176 of file matrix_impl.hpp.

◆ get_colind()

std::vector<casadi_int> casadi::GenericMatrix< Matrix< Scalar > >::get_colind ( ) const
inlineinherited

Extra doc: https://github.com/casadi/casadi/wiki/L_1b8

Definition at line 195 of file generic_matrix.hpp.

◆ get_elements()

template<typename Scalar >
std::vector<Scalar> casadi::Matrix< Scalar >::get_elements ( ) const

◆ get_free()

std::vector< SX > casadi::SX::get_free ( const Function f)
static

◆ get_input()

std::vector< SX > casadi::SX::get_input ( const Function f)
static

◆ get_max_depth()

casadi_int casadi::SX::get_max_depth ( )
static

◆ get_nonzeros()

template<typename Scalar >
std::vector< A > casadi::Matrix< Scalar >::get_nonzeros ( ) const

Implementation of Matrix::get_nonzeros (in public API)

Extra doc: https://github.com/casadi/casadi/wiki/L_1a2

Definition at line 1325 of file matrix_decl.hpp.

◆ get_nz() [1/2]

template<typename Scalar >
void casadi::Matrix< Scalar >::get_nz ( Matrix< Scalar > &  m,
bool  ind1,
const Matrix< casadi_int > &  k 
) const

Get a set of nonzeros

Definition at line 406 of file matrix_impl.hpp.

◆ get_nz() [2/2]

template<typename Scalar >
void casadi::Matrix< Scalar >::get_nz ( Matrix< Scalar > &  m,
bool  ind1,
const Slice k 
) const

Get a set of nonzeros

Definition at line 394 of file matrix_impl.hpp.

◆ get_row()

std::vector<casadi_int> casadi::GenericMatrix< Matrix< Scalar > >::get_row ( ) const
inlineinherited

Extra doc: https://github.com/casadi/casadi/wiki/L_1b8

Definition at line 194 of file generic_matrix.hpp.

◆ get_sparsity()

template<typename Scalar >
Sparsity casadi::Matrix< Scalar >::get_sparsity ( ) const

◆ get_str()

template<typename Scalar >
std::string casadi::Matrix< Scalar >::get_str ( bool  more = false) const

◆ gradient()

static Matrix< Scalar > casadi::GenericMatrix< Matrix< Scalar > >::gradient ( const Matrix< Scalar > &  ex,
const Matrix< Scalar > &  arg,
const Dict opts = Dict() 
)
staticinherited

Functions called by friend functions defined here

◆ has_nz()

template<typename Scalar >
bool casadi::Matrix< Scalar >::has_nz ( casadi_int  rr,
casadi_int  cc 
) const

Definition at line 65 of file matrix_impl.hpp.

◆ has_zeros()

template<typename Scalar >
bool casadi::Matrix< Scalar >::has_zeros ( ) const

◆ inf() [1/3]

template<typename Scalar >
static Matrix<Scalar> casadi::Matrix< Scalar >::inf ( casadi_int  nrow = 1,
casadi_int  ncol = 1 
)
static

◆ inf() [2/3]

template<typename Scalar >
static Matrix<Scalar> casadi::Matrix< Scalar >::inf ( const Sparsity sp)
static

◆ inf() [3/3]

template<typename Scalar >
static Matrix<Scalar> casadi::Matrix< Scalar >::inf ( const std::pair< casadi_int, casadi_int > &  rc)
static

◆ info()

Dict CASADI_EXPORT casadi::SX::info ( ) const

Obtain information about sparsity

◆ is_column()

bool casadi::GenericMatrix< Matrix< Scalar > >::is_column ( ) const
inlineinherited

Extra doc: https://github.com/casadi/casadi/wiki/L_1b5

Definition at line 178 of file generic_matrix.hpp.

◆ is_commutative()

bool casadi::SX::is_commutative ( ) const

Only defined if symbolic scalar.

Extra doc: https://github.com/casadi/casadi/wiki/L_19q

◆ is_constant()

template<typename Scalar >
bool casadi::Matrix< Scalar >::is_constant ( ) const

◆ is_dense()

bool casadi::GenericMatrix< Matrix< Scalar > >::is_dense ( ) const
inlineinherited

Extra doc: https://github.com/casadi/casadi/wiki/L_1b0

Definition at line 153 of file generic_matrix.hpp.

◆ is_empty()

bool casadi::GenericMatrix< Matrix< Scalar > >::is_empty ( bool  both = false) const
inlineinherited

(or optionally both dimensions)

Extra doc: https://github.com/casadi/casadi/wiki/L_1az

Definition at line 148 of file generic_matrix.hpp.

◆ is_eye()

template<typename Scalar >
bool casadi::Matrix< Scalar >::is_eye ( ) const

◆ is_integer()

template<typename Scalar >
bool casadi::Matrix< Scalar >::is_integer ( ) const

(note that false negative answers are possible)

Extra doc: https://github.com/casadi/casadi/wiki/L_19w

◆ is_leaf()

bool casadi::SX::is_leaf ( ) const

Only defined if symbolic scalar.

Extra doc: https://github.com/casadi/casadi/wiki/L_19p

◆ is_minus_one()

template<typename Scalar >
bool casadi::Matrix< Scalar >::is_minus_one ( ) const

◆ is_one()

template<typename Scalar >
bool casadi::Matrix< Scalar >::is_one ( ) const

◆ is_op()

bool casadi::SX::is_op ( casadi_int  op) const

◆ is_regular()

bool casadi::SX::is_regular ( ) const

◆ is_row()

bool casadi::GenericMatrix< Matrix< Scalar > >::is_row ( ) const
inlineinherited

Extra doc: https://github.com/casadi/casadi/wiki/L_1b4

Definition at line 173 of file generic_matrix.hpp.

◆ is_scalar()

bool casadi::GenericMatrix< Matrix< Scalar > >::is_scalar ( bool  scalar_and_dense = false) const
inherited

◆ is_smooth()

bool casadi::SX::is_smooth ( ) const

◆ is_square()

bool casadi::GenericMatrix< Matrix< Scalar > >::is_square ( ) const
inlineinherited

Extra doc: https://github.com/casadi/casadi/wiki/L_1b2

Definition at line 163 of file generic_matrix.hpp.

◆ is_symbolic()

bool casadi::SX::is_symbolic ( ) const

Check if an expression is a pure symbol. Pure means that no operations should happen on the symbol (not even vec, transpose, index, concatenation, ...) By consequence, a slice of a vector-shaped MX symbol is not a symbol. However, the SX type is really a container format with scalar entries. Therefore, a slice of a vector-shaped SX symbol is still a symbol.

Sparse matrices invariable return false

\seealso is_valid_input

Extra doc: https://github.com/casadi/casadi/wiki/L_19r

◆ is_tril()

bool casadi::GenericMatrix< Matrix< Scalar > >::is_tril ( ) const
inlineinherited

Extra doc: https://github.com/casadi/casadi/wiki/L_1b7

Definition at line 188 of file generic_matrix.hpp.

◆ is_triu()

bool casadi::GenericMatrix< Matrix< Scalar > >::is_triu ( ) const
inlineinherited

Extra doc: https://github.com/casadi/casadi/wiki/L_1b6

Definition at line 183 of file generic_matrix.hpp.

◆ is_valid_input()

bool casadi::SX::is_valid_input ( ) const

is_valid_input is more forgiving than is_symbolic. Some compositions are allowed: vec, vertcat.

Sparse matrices can return true if all non-zero elements are symbolic

\seealso is_symbolic

Extra doc: https://github.com/casadi/casadi/wiki/L_19s

◆ is_vector()

bool casadi::GenericMatrix< Matrix< Scalar > >::is_vector ( ) const
inlineinherited

Extra doc: https://github.com/casadi/casadi/wiki/L_1b3

Definition at line 168 of file generic_matrix.hpp.

◆ is_zero()

template<typename Scalar >
bool casadi::Matrix< Scalar >::is_zero ( ) const

◆ jtimes()

static Matrix< Scalar > casadi::GenericMatrix< Matrix< Scalar > >::jtimes ( const Matrix< Scalar > &  ex,
const Matrix< Scalar > &  arg,
const Matrix< Scalar > &  v,
bool  tr = false,
const Dict opts = Dict() 
)
staticinherited

Functions called by friend functions defined here

◆ linearize()

static Matrix< Scalar > casadi::GenericMatrix< Matrix< Scalar > >::linearize ( const Matrix< Scalar > &  f,
const Matrix< Scalar > &  x,
const Matrix< Scalar > &  x0,
const Dict opts = Dict() 
)
staticinherited

Functions called by friend functions defined here

◆ logsumexp()

static Matrix< Scalar > casadi::GenericMatrix< Matrix< Scalar > >::logsumexp ( const Matrix< Scalar > &  x)
staticinherited

◆ mpower()

static Matrix< Scalar > casadi::GenericMatrix< Matrix< Scalar > >::mpower ( const Matrix< Scalar > &  x,
const Matrix< Scalar > &  y 
)
staticinherited

Functions called by friend functions defined here

◆ n_dep()

casadi_int casadi::SX::n_dep ( ) const

Only defined if symbolic scalar.

Extra doc: https://github.com/casadi/casadi/wiki/L_1aa

◆ name()

std::string casadi::SX::name ( ) const

◆ nan() [1/3]

template<typename Scalar >
static Matrix<Scalar> casadi::Matrix< Scalar >::nan ( casadi_int  nrow = 1,
casadi_int  ncol = 1 
)
static

◆ nan() [2/3]

template<typename Scalar >
static Matrix<Scalar> casadi::Matrix< Scalar >::nan ( const Sparsity sp)
static

◆ nan() [3/3]

template<typename Scalar >
static Matrix<Scalar> casadi::Matrix< Scalar >::nan ( const std::pair< casadi_int, casadi_int > &  rc)
static

◆ nnz()

casadi_int casadi::GenericMatrix< Matrix< Scalar > >::nnz ( ) const
inherited

◆ nnz_diag()

casadi_int casadi::GenericMatrix< Matrix< Scalar > >::nnz_diag ( ) const
inherited

◆ nnz_lower()

casadi_int casadi::GenericMatrix< Matrix< Scalar > >::nnz_lower ( ) const
inherited

◆ nnz_upper()

casadi_int casadi::GenericMatrix< Matrix< Scalar > >::nnz_upper ( ) const
inherited

◆ numel()

casadi_int casadi::GenericMatrix< Matrix< Scalar > >::numel ( ) const
inherited

◆ ones() [1/3]

static Matrix< Scalar > casadi::GenericMatrix< Matrix< Scalar > >::ones ( casadi_int  nrow = 1,
casadi_int  ncol = 1 
)
inlinestaticinherited

◆ ones() [2/3]

static Matrix< Scalar > casadi::GenericMatrix< Matrix< Scalar > >::ones ( const Sparsity sp)
inlinestaticinherited

◆ ones() [3/3]

static Matrix< Scalar > casadi::GenericMatrix< Matrix< Scalar > >::ones ( const std::pair< casadi_int, casadi_int > &  rc)
inlinestaticinherited

◆ op()

casadi_int casadi::SX::op ( ) const

◆ operator casadi_int()

template<typename Scalar >
casadi::Matrix< Scalar >::operator casadi_int ( ) const
explicit

◆ operator double()

template<typename Scalar >
casadi::Matrix< Scalar >::operator double ( ) const
explicit

◆ operator+()

template<typename Scalar >
Matrix<Scalar> casadi::Matrix< Scalar >::operator+ ( ) const
Examples
SX/SXaddition.cc.

◆ operator-()

template<typename Scalar >
Matrix<Scalar> casadi::Matrix< Scalar >::operator- ( ) const

◆ print_dense()

template<typename Scalar >
void casadi::Matrix< Scalar >::print_dense ( std::ostream &  stream,
bool  truncate = true 
) const

Definition at line 686 of file matrix_impl.hpp.

◆ print_scalar()

template<typename Scalar >
void casadi::Matrix< Scalar >::print_scalar ( std::ostream &  stream) const

Definition at line 609 of file matrix_impl.hpp.

◆ print_sparse()

template<typename Scalar >
void casadi::Matrix< Scalar >::print_sparse ( std::ostream &  stream,
bool  truncate = true 
) const

Definition at line 691 of file matrix_impl.hpp.

◆ print_split()

template<typename Scalar >
void casadi::Matrix< Scalar >::print_split ( std::vector< std::string > &  nz,
std::vector< std::string > &  inter 
) const

Definition at line 696 of file matrix_impl.hpp.

◆ print_vector()

template<typename Scalar >
void casadi::Matrix< Scalar >::print_vector ( std::ostream &  stream,
bool  truncate = true 
) const

Definition at line 636 of file matrix_impl.hpp.

◆ printme()

template<typename Scalar >
Matrix<Scalar> casadi::Matrix< Scalar >::printme ( const Matrix< Scalar > &  y) const

◆ rand() [1/3]

template<typename Scalar >
static Matrix<Scalar> casadi::Matrix< Scalar >::rand ( casadi_int  nrow = 1,
casadi_int  ncol = 1 
)
static

◆ rand() [2/3]

DM casadi::DM::rand ( const Sparsity sp)
static

◆ rand() [3/3]

template<typename Scalar >
static Matrix<Scalar> casadi::Matrix< Scalar >::rand ( const std::pair< casadi_int, casadi_int > &  rc)
static

◆ rank1()

static Matrix< Scalar > casadi::GenericMatrix< Matrix< Scalar > >::rank1 ( const Matrix< Scalar > &  A,
const Matrix< Scalar > &  alpha,
const Matrix< Scalar > &  x,
const Matrix< Scalar > &  y 
)
staticinherited

Calculates A + 1/2 * alpha * x*y'

Extra doc: https://github.com/casadi/casadi/wiki/L_1bp

◆ remove()

template<typename Scalar >
void casadi::Matrix< Scalar >::remove ( const std::vector< casadi_int > &  rr,
const std::vector< casadi_int > &  cc 
)

Remove/delete rows and/or columns of a matrix

Extra doc: https://github.com/casadi/casadi/wiki/L_19h

◆ reserve() [1/2]

template<typename Scalar >
void casadi::Matrix< Scalar >::reserve ( casadi_int  nnz)

◆ reserve() [2/2]

template<typename Scalar >
void casadi::Matrix< Scalar >::reserve ( casadi_int  nnz,
casadi_int  ncol 
)

◆ resize()

template<typename Scalar >
void casadi::Matrix< Scalar >::resize ( casadi_int  nrow,
casadi_int  ncol 
)

◆ rng()

template<typename Scalar >
void casadi::Matrix< Scalar >::rng ( casadi_int  seed)
static

Definition at line 60 of file matrix_impl.hpp.

◆ row()

casadi_int casadi::GenericMatrix< Matrix< Scalar > >::row ( casadi_int  el) const
inlineinherited

Extra doc: https://github.com/casadi/casadi/wiki/L_1b8

Definition at line 200 of file generic_matrix.hpp.

◆ rows()

casadi_int casadi::GenericMatrix< Matrix< Scalar > >::rows ( ) const
inlineinherited

Extra doc: https://github.com/casadi/casadi/wiki/L_1at

Definition at line 114 of file generic_matrix.hpp.

◆ serialize() [1/2]

template<typename Scalar >
std::string casadi::Matrix< Scalar >::serialize ( ) const

◆ serialize() [2/2]

template<typename Scalar >
void casadi::Matrix< Scalar >::serialize ( SerializingStream s) const

◆ set() [1/7]

template<typename Scalar >
void casadi::Matrix< Scalar >::set ( const Matrix< Scalar > &  m,
bool  ind1,
const Matrix< casadi_int > &  rr 
)

Set a submatrix, single argument

Definition at line 304 of file matrix_impl.hpp.

◆ set() [2/7]

template<typename Scalar >
void casadi::Matrix< Scalar >::set ( const Matrix< Scalar > &  m,
bool  ind1,
const Matrix< casadi_int > &  rr,
const Matrix< casadi_int > &  cc 
)

Set a submatrix, two arguments

Definition at line 218 of file matrix_impl.hpp.

◆ set() [3/7]

template<typename Scalar >
void casadi::Matrix< Scalar >::set ( const Matrix< Scalar > &  m,
bool  ind1,
const Matrix< casadi_int > &  rr,
const Slice cc 
)

Set a submatrix, two arguments

Definition at line 211 of file matrix_impl.hpp.

◆ set() [4/7]

template<typename Scalar >
void casadi::Matrix< Scalar >::set ( const Matrix< Scalar > &  m,
bool  ind1,
const Slice rr 
)

Set a submatrix, single argument

Definition at line 285 of file matrix_impl.hpp.

◆ set() [5/7]

template<typename Scalar >
void casadi::Matrix< Scalar >::set ( const Matrix< Scalar > &  m,
bool  ind1,
const Slice rr,
const Matrix< casadi_int > &  cc 
)

Set a submatrix, two arguments

Definition at line 204 of file matrix_impl.hpp.

◆ set() [6/7]

template<typename Scalar >
void casadi::Matrix< Scalar >::set ( const Matrix< Scalar > &  m,
bool  ind1,
const Slice rr,
const Slice cc 
)

Set a submatrix, two arguments

Definition at line 185 of file matrix_impl.hpp.

◆ set() [7/7]

template<typename Scalar >
void casadi::Matrix< Scalar >::set ( const Matrix< Scalar > &  m,
bool  ind1,
const Sparsity sp 
)

Set a submatrix, single argument

Definition at line 380 of file matrix_impl.hpp.

◆ set_max_depth()

void casadi::SX::set_max_depth ( casadi_int  eq_depth = 1)
static

◆ set_nz() [1/2]

template<typename Scalar >
void casadi::Matrix< Scalar >::set_nz ( const Matrix< Scalar > &  m,
bool  ind1,
const Matrix< casadi_int > &  k 
)

Set a set of nonzeros

Definition at line 447 of file matrix_impl.hpp.

◆ set_nz() [2/2]

template<typename Scalar >
void casadi::Matrix< Scalar >::set_nz ( const Matrix< Scalar > &  m,
bool  ind1,
const Slice k 
)

Set a set of nonzeros

Definition at line 435 of file matrix_impl.hpp.

◆ set_precision()

template<typename Scalar >
void casadi::Matrix< Scalar >::set_precision ( casadi_int  precision)
static

Definition at line 39 of file matrix_impl.hpp.

◆ set_scientific()

template<typename Scalar >
void casadi::Matrix< Scalar >::set_scientific ( bool  scientific)
static

Definition at line 45 of file matrix_impl.hpp.

◆ set_width()

template<typename Scalar >
void casadi::Matrix< Scalar >::set_width ( casadi_int  width)
static

Definition at line 42 of file matrix_impl.hpp.

◆ size() [1/2]

std::pair<casadi_int, casadi_int> casadi::GenericMatrix< Matrix< Scalar > >::size ( ) const
inherited

◆ size() [2/2]

casadi_int casadi::GenericMatrix< Matrix< Scalar > >::size ( casadi_int  axis) const
inherited

◆ size1()

casadi_int casadi::GenericMatrix< Matrix< Scalar > >::size1 ( ) const
inherited

◆ size2()

casadi_int casadi::GenericMatrix< Matrix< Scalar > >::size2 ( ) const
inherited

◆ soc()

static Matrix< Scalar > casadi::GenericMatrix< Matrix< Scalar > >::soc ( const Matrix< Scalar > &  x,
const Matrix< Scalar > &  y 
)
staticinherited

Functions called by friend functions defined here

◆ sparsity()

Sparsity casadi::GenericMatrix< Matrix< Scalar > >::sparsity ( ) const
inherited

◆ sym() [1/7]

static std::vector<Matrix< Scalar > > casadi::GenericMatrix< Matrix< Scalar > >::sym ( const std::string &  name,
casadi_int  nrow,
casadi_int  ncol,
casadi_int  p 
)
inlinestaticinherited

◆ sym() [2/7]

static std::vector<std::vector<Matrix< Scalar > > > casadi::GenericMatrix< Matrix< Scalar > >::sym ( const std::string &  name,
casadi_int  nrow,
casadi_int  ncol,
casadi_int  p,
casadi_int  r 
)
inlinestaticinherited

with nrow-by-ncol symbolic primitives

Extra doc: https://github.com/casadi/casadi/wiki/L_1dg

Definition at line 1074 of file generic_matrix.hpp.

◆ sym() [3/7]

static Matrix< Scalar > casadi::GenericMatrix< Matrix< Scalar > >::sym ( const std::string &  name,
casadi_int  nrow = 1,
casadi_int  ncol = 1 
)
inlinestaticinherited

◆ sym() [4/7]

static Matrix< Scalar > casadi::GenericMatrix< Matrix< Scalar > >::sym ( const std::string &  name,
const Sparsity sp 
)
inlinestaticinherited

◆ sym() [5/7]

std::vector< Matrix< Scalar > > casadi::GenericMatrix< Matrix< Scalar > >::sym ( const std::string &  name,
const Sparsity sp,
casadi_int  p 
)
staticinherited

with symbolic primitives of given sparsity

Extra doc: https://github.com/casadi/casadi/wiki/L_1dd

Definition at line 1050 of file generic_matrix.hpp.

◆ sym() [6/7]

std::vector< std::vector< Matrix< Scalar > > > casadi::GenericMatrix< Matrix< Scalar > >::sym ( const std::string &  name,
const Sparsity sp,
casadi_int  p,
casadi_int  r 
)
staticinherited

symbolic primitives with given sparsity

Extra doc: https://github.com/casadi/casadi/wiki/L_1df

Definition at line 1066 of file generic_matrix.hpp.

◆ sym() [7/7]

static Matrix< Scalar > casadi::GenericMatrix< Matrix< Scalar > >::sym ( const std::string &  name,
const std::pair< casadi_int, casadi_int > &  rc 
)
inlinestaticinherited

◆ T()

template<typename Scalar >
Matrix<Scalar> casadi::Matrix< Scalar >::T ( ) const

◆ tangent()

static Matrix< Scalar > casadi::GenericMatrix< Matrix< Scalar > >::tangent ( const Matrix< Scalar > &  ex,
const Matrix< Scalar > &  arg,
const Dict opts = Dict() 
)
staticinherited

Functions called by friend functions defined here

◆ to_file()

template<typename Scalar >
void casadi::Matrix< Scalar >::to_file ( const std::string &  filename,
const std::string &  format = "" 
) const

Export numerical matrix to file

Supported formats:

*   - .mtx   Matrix Market (sparse)
*   - .txt   Ascii full precision representation (sparse)
*            Whitespace separated, aligned.
*            Comments with # % or /
*            Uses C locale
*            Structural zeros represented by 00
*            Does not scale well for large sparse matrices
* 

◆ triplet() [1/3]

template<typename Scalar >
static Matrix<Scalar> casadi::Matrix< Scalar >::triplet ( const std::vector< casadi_int > &  row,
const std::vector< casadi_int > &  col,
const Matrix< Scalar > &  d 
)
static

Default matrix size is max(col) x max(row)

Extra doc: https://github.com/casadi/casadi/wiki/L_23t

◆ triplet() [2/3]

template<typename Scalar >
static Matrix<Scalar> casadi::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 
)
static

Default matrix size is max(col) x max(row)

Extra doc: https://github.com/casadi/casadi/wiki/L_23t

◆ triplet() [3/3]

template<typename Scalar >
static Matrix<Scalar> casadi::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 
)
static

Default matrix size is max(col) x max(row)

Extra doc: https://github.com/casadi/casadi/wiki/L_23t

◆ type_name()

template<typename Scalar >
std::string casadi::Matrix< Scalar >::type_name
static

Definition at line 606 of file matrix_impl.hpp.

◆ zeros() [1/3]

static Matrix< Scalar > casadi::GenericMatrix< Matrix< Scalar > >::zeros ( casadi_int  nrow = 1,
casadi_int  ncol = 1 
)
inlinestaticinherited

Extra doc: https://github.com/casadi/casadi/wiki/L_1dh

Examples
misc/all.cc.

Definition at line 1083 of file generic_matrix.hpp.

◆ zeros() [2/3]

static Matrix< Scalar > casadi::GenericMatrix< Matrix< Scalar > >::zeros ( const Sparsity sp)
inlinestaticinherited

◆ zeros() [3/3]

static Matrix< Scalar > casadi::GenericMatrix< Matrix< Scalar > >::zeros ( const std::pair< casadi_int, casadi_int > &  rc)
inlinestaticinherited

The documentation for this class was generated from the following files: