27 #ifndef CASADI_SPARSITY_HPP
28 #define CASADI_SPARSITY_HPP
30 #include "shared_object.hpp"
31 #include "printable.hpp"
32 #include "casadi_common.hpp"
33 #include "sparsity_interface.hpp"
34 #include "generic_type.hpp"
38 #include <unordered_map>
42 class SparsityInternal;
43 class SerializingStream;
44 class DeserializingStream;
50 struct CASADI_EXPORT SparsityStruct {
53 const casadi_int* colind;
54 const casadi_int* row;
98 public SWIG_IF_ELSE(SparsityInterfaceCommon, SparsityInterface<Sparsity>),
99 public SWIG_IF_ELSE(PrintableCommon, Printable<Sparsity>) {
112 const std::vector<casadi_int>& colind,
const std::vector<casadi_int>& row,
113 bool order_rows=
false);
118 explicit Sparsity(
const std::pair<casadi_int, casadi_int>& rc);
122 Sparsity(casadi_int nrow, casadi_int ncol,
const casadi_int* colind,
const casadi_int* row,
123 bool order_rows=
false);
128 static Sparsity create(SparsityInternal *node);
131 typedef SparsityInterface<Sparsity> B;
139 SparsityInternal* get()
const;
147 {
return dense_scalar ? dense(1, 1) :
Sparsity(1, 1); }
156 return dense(rc.first, rc.second);
186 return diag(rc.first, rc.second);
211 const std::vector<casadi_int>& col,
212 casadi_int nrow, casadi_int ncol);
219 const std::vector<casadi_int>& row,
const std::vector<casadi_int>& col,
220 std::vector<casadi_int>& SWIG_OUTPUT(mapping),
bool invert_mapping);
221 static Sparsity triplet(casadi_int nrow, casadi_int ncol,
const std::vector<casadi_int>& row,
222 const std::vector<casadi_int>& col);
230 static Sparsity nonzeros(casadi_int nrow, casadi_int ncol,
const std::vector<casadi_int>& nz,
231 bool ind1=SWIG_IND1);
243 static Sparsity compressed(
const casadi_int* v,
bool order_rows=
false);
276 const SparsityInternal* operator->()
const;
279 const SparsityInternal& operator*()
const;
283 bool is_equal(const Sparsity& y) const;
284 bool is_equal(casadi_int nrow, casadi_int ncol,
const std::vector<casadi_int>& colind,
285 const std::vector<casadi_int>& row)
const;
287 bool is_equal(casadi_int nrow, casadi_int ncol,
288 const casadi_int* colind,
const casadi_int* row)
const;
307 operator const casadi_int*()
const;
312 operator const std::vector<casadi_int>&()
const;
317 operator SparsityStruct()
const;
327 casadi_int
rows()
const {
return size1();}
397 std::pair<casadi_int, casadi_int>
size()
const;
402 casadi_int
size(casadi_int axis)
const;
413 void to_file(
const std::string& filename,
const std::string& format_hint=
"")
const;
421 void serialize(std::ostream &stream)
const;
455 const casadi_int* row()
const;
460 const casadi_int* colind()
const;
488 casadi_int
row(casadi_int el)
const;
499 void resize(casadi_int nrow, casadi_int ncol);
506 casadi_int
add_nz(casadi_int rr, casadi_int cc);
513 casadi_int
get_nz(casadi_int rr, casadi_int cc)
const;
516 bool has_nz(casadi_int rr, casadi_int cc)
const;
523 std::vector<casadi_int>
get_nz(
const std::vector<casadi_int>& rr,
524 const std::vector<casadi_int>& cc)
const;
533 void get_nz(std::vector<casadi_int>& SWIG_INOUT(indices))
const;
542 void get_ccs(std::vector<casadi_int>& SWIG_OUTPUT(colind),
543 std::vector<casadi_int>& SWIG_OUTPUT(row))
const;
546 void get_crs(std::vector<casadi_int>& SWIG_OUTPUT(rowind),
547 std::vector<casadi_int>& SWIG_OUTPUT(col))
const;
551 std::vector<casadi_int>& SWIG_OUTPUT(col))
const;
560 const std::vector<casadi_int>& cc,
561 std::vector<casadi_int>& SWIG_OUTPUT(mapping),
bool ind1=
false)
const;
570 std::vector<casadi_int>& SWIG_OUTPUT(mapping),
bool ind1=
false)
const;
582 bool invert_mapping=
false)
const;
602 Sparsity combine(
const Sparsity& y,
bool f0x_is_zero,
bool function0_is_zero,
603 std::vector<unsigned char>& mapping)
const;
613 Sparsity unite(
const Sparsity& y, std::vector<unsigned char>& mapping)
const;
630 std::vector<unsigned char>& mapping)
const;
664 static void mul_sparsityF(
const bvec_t* x,
const Sparsity& x_sp,
665 const bvec_t* y,
const Sparsity& y_sp,
675 static void mul_sparsityR(bvec_t* x,
const Sparsity& x_sp,
685 static Sparsity horzcat(
const std::vector<Sparsity> & sp);
686 static Sparsity vertcat(
const std::vector<Sparsity> & sp);
687 static Sparsity blockcat(
const std::vector< std::vector< Sparsity > > &v);
688 static Sparsity diagcat(
const std::vector< Sparsity > &v);
689 static std::vector<Sparsity>
690 horzsplit(
const Sparsity& x,
const std::vector<casadi_int>& offset);
691 static std::vector<Sparsity>
692 vertsplit(
const Sparsity& x,
const std::vector<casadi_int>& offset);
693 static std::vector<Sparsity>
695 const std::vector<casadi_int>& offset1,
696 const std::vector<casadi_int>& offset2);
699 static Sparsity reshape(
const Sparsity& x, casadi_int nrow, casadi_int ncol);
702 static casadi_int sprank(
const Sparsity& x);
726 void enlarge(casadi_int nrow, casadi_int ncol,
727 const std::vector<casadi_int>& rr,
728 const std::vector<casadi_int>& cc,
bool ind1=
false);
733 void enlargeRows(casadi_int nrow,
const std::vector<casadi_int>& rr,
bool ind1=
false);
738 void enlargeColumns(casadi_int ncol,
const std::vector<casadi_int>& cc,
bool ind1=
false);
748 std::vector<casadi_int>
erase(
const std::vector<casadi_int>& rr,
749 const std::vector<casadi_int>& cc,
bool ind1=
false);
754 std::vector<casadi_int>
erase(
const std::vector<casadi_int>& rr,
bool ind1=
false);
874 typedef std::unordered_multimap<std::size_t, WeakRef> CachingMap;
877 static CachingMap& getCache();
883 static const Sparsity& getScalarSparse();
902 std::vector<casadi_int>
etree(
bool ata=
false)
const;
913 Sparsity ldl(std::vector<casadi_int>& SWIG_OUTPUT(p),
bool amd=
true)
const;
925 std::vector<casadi_int>& SWIG_OUTPUT(prinv),
926 std::vector<casadi_int>& SWIG_OUTPUT(pc),
bool amd=
true)
const;
933 casadi_int
dfs(casadi_int j, casadi_int top, std::vector<casadi_int>& SWIG_INOUT(xi),
934 std::vector<casadi_int>& SWIG_INOUT(pstack),
935 const std::vector<casadi_int>& pinv, std::vector<bool>& SWIG_INOUT(marked))
const;
960 casadi_int
scc(std::vector<casadi_int>& SWIG_OUTPUT(index),
961 std::vector<casadi_int>& SWIG_OUTPUT(offset))
const;
984 casadi_int
btf(std::vector<casadi_int>& SWIG_OUTPUT(rowperm),
985 std::vector<casadi_int>& SWIG_OUTPUT(colperm),
986 std::vector<casadi_int>& SWIG_OUTPUT(rowblock),
987 std::vector<casadi_int>& SWIG_OUTPUT(colblock),
988 std::vector<casadi_int>& SWIG_OUTPUT(coarse_rowblock),
989 std::vector<casadi_int>& SWIG_OUTPUT(coarse_colblock))
const;
1003 std::vector<casadi_int>
amd()
const;
1009 void spsolve(bvec_t* X, bvec_t* B,
bool tr)
const;
1023 std::vector<casadi_int>
find(
bool ind1=SWIG_IND1)
const;
1027 void find(std::vector<casadi_int>& loc,
bool ind1=
false)
const;
1036 casadi_int cutoff = std::numeric_limits<casadi_int>::max())
const;
1050 casadi_int cutoff = std::numeric_limits<casadi_int>::max())
const;
1064 casadi_int cutoff = std::numeric_limits<casadi_int>::max())
const;
1079 bool permute_rows=
true,
bool permute_columns=
true,
1080 bool invert_permutation=
false)
const;
1083 std::string
dim(
bool with_nz=
false)
const;
1142 bool with_x_diag=
true,
bool with_lam_g_diag=
true);
1150 template<
typename T>
1151 void set(T* data,
const T* val_data,
const Sparsity& val_sp)
const;
1158 template<
typename T>
1159 void add(T* data,
const T* val_data,
const Sparsity& val_sp)
const;
1166 template<
typename T>
1167 void bor(T* data,
const T* val_data,
const Sparsity& val_sp)
const;
1169 static std::string file_format(
const std::string& filename,
1170 const std::string& format_hint,
const std::set<std::string>& file_formats);
1171 static std::set<std::string> file_formats;
1175 void assign_cached(casadi_int nrow, casadi_int ncol,
const std::vector<casadi_int>& colind,
1176 const std::vector<casadi_int>& row,
bool order_rows=
false);
1179 void assign_cached(casadi_int nrow, casadi_int ncol,
1180 const casadi_int* colind,
const casadi_int* row,
bool order_rows=
false);
1188 template<
typename T>
1194 template<
typename T>
1196 seed ^=
hash_value(v) + 0x9e3779b9 + (seed << 6) + (seed >> 2);
1202 template<
typename T>
1204 for (casadi_int i=0; i<sz; ++i)
hash_combine(seed, v[i]);
1210 template<
typename T>
1226 const std::vector<casadi_int>& colind,
1227 const std::vector<casadi_int>& row);
1230 const casadi_int* colind,
1231 const casadi_int* row);
1235 template<
typename DataType>
1236 void Sparsity::set(DataType* data,
const DataType* val_data,
const Sparsity& val_sp)
const {
1238 const casadi_int sz =
nnz();
1239 const casadi_int sz1 =
size1();
1240 const casadi_int sz2 =
size2();
1243 const casadi_int val_sz = val_sp.
nnz();
1244 const casadi_int val_sz1 = val_sp.
size1();
1245 const casadi_int val_sz2 = val_sp.
size2();
1246 const casadi_int val_nel = val_sz1*val_sz2;
1249 if (val_sp==*
this) {
1250 std::copy(val_data, val_data+sz, data);
1257 }
else if (val_nel==1) {
1258 std::fill(data, data+sz, val_sz==0 ? DataType(0) : val_data[0]);
1259 }
else if (sz2==val_sz2 && sz1==val_sz1) {
1262 const casadi_int* c =
row();
1263 const casadi_int* rind =
colind();
1264 const casadi_int* v_c = val_sp.
row();
1265 const casadi_int* v_rind = val_sp.
colind();
1268 for (casadi_int i=0; i<sz2; ++i) {
1271 casadi_int v_el = v_rind[i];
1274 casadi_int v_el_end = v_rind[i+1];
1277 casadi_int v_j = v_el<v_el_end ? v_c[v_el] : sz1;
1280 for (casadi_int el=rind[i]; el!=rind[i+1]; ++el) {
1288 v_j = v_el<v_el_end ? v_c[v_el] : sz1;
1293 data[el] = val_data[v_el++];
1294 v_j = v_el<v_el_end ? v_c[v_el] : sz1;
1300 }
else if (sz1==val_sz2 && sz2==val_sz1 && sz2 == 1) {
1302 const casadi_int* v_cind = val_sp.
colind();
1303 const casadi_int* r =
row();
1304 for (casadi_int el=0; el<sz; ++el) {
1305 casadi_int rr=r[el];
1306 data[el] = v_cind[rr]==v_cind[rr+1] ? 0 : val_data[v_cind[rr]];
1308 }
else if (sz1==val_sz2 && sz2==val_sz1 && sz1 == 1) {
1310 for (casadi_int el=0; el<sz; ++el) data[el] = 0;
1311 const casadi_int* cind =
colind();
1312 const casadi_int* v_r = val_sp.
row();
1313 for (casadi_int el=0; el<val_sz; ++el) {
1314 casadi_int rr=v_r[el];
1315 if (cind[rr]!=cind[rr+1]) {
1316 data[cind[rr]] = val_data[el];
1321 casadi_error(
"Sparsity::set<DataType>: shape mismatch. lhs is "
1322 +
dim() +
", while rhs is " + val_sp.
dim() +
".");
1326 template<
typename DataType>
1327 void Sparsity::add(DataType* data,
const DataType* val_data,
const Sparsity& val_sp)
const {
1329 const casadi_int sz =
nnz();
1330 const casadi_int sz1 =
size1();
1331 const casadi_int sz2 =
size2();
1332 const casadi_int nel = sz1*sz2;
1335 const casadi_int val_sz = val_sp.nnz();
1336 const casadi_int val_sz1 = val_sp.size1();
1337 const casadi_int val_sz2 = val_sp.size2();
1338 const casadi_int val_nel = val_sz1*val_sz2;
1341 if (val_sp==*
this) {
1342 for (casadi_int k=0; k<sz; ++k) {
1343 data[k] += val_data[k];
1348 }
else if (val_sp.is_empty()) {
1351 }
else if (val_nel==1) {
1353 for (casadi_int k=0; k<sz; ++k) {
1354 data[k] += val_data[0];
1359 if (nel==0 && val_nel==0)
return;
1362 casadi_assert(sz2==val_sz2 && sz1==val_sz1,
1363 "Sparsity::add<DataType>: shape mismatch. lhs is "
1364 +
dim() +
", while rhs is " + val_sp.dim() +
".");
1367 const casadi_int* c =
row();
1368 const casadi_int* rind =
colind();
1369 const casadi_int* v_c = val_sp.row();
1370 const casadi_int* v_rind = val_sp.colind();
1373 for (casadi_int i=0; i<sz2; ++i) {
1376 casadi_int v_el = v_rind[i];
1379 casadi_int v_el_end = v_rind[i+1];
1382 casadi_int v_j = v_el<v_el_end ? v_c[v_el] : sz1;
1385 for (casadi_int el=rind[i]; el!=rind[i+1]; ++el) {
1393 v_j = v_el<v_el_end ? v_c[v_el] : sz1;
1398 data[el] += val_data[v_el++];
1399 v_j = v_el<v_el_end ? v_c[v_el] : sz1;
1406 template<
typename DataType>
1407 void Sparsity::bor(DataType* data,
const DataType* val_data,
const Sparsity& val_sp)
const {
1409 const casadi_int sz =
nnz();
1410 const casadi_int sz1 =
size1();
1411 const casadi_int sz2 =
size2();
1412 const casadi_int nel = sz1*sz2;
1415 const casadi_int val_sz = val_sp.nnz();
1416 const casadi_int val_sz1 = val_sp.size1();
1417 const casadi_int val_sz2 = val_sp.size2();
1418 const casadi_int val_nel = val_sz1*val_sz2;
1421 if (val_sp==*
this) {
1422 for (casadi_int k=0; k<sz; ++k) {
1423 data[k] |= val_data[k];
1428 }
else if (val_sp.is_empty()) {
1431 }
else if (val_nel==1) {
1433 for (casadi_int k=0; k<sz; ++k) {
1434 data[k] |= val_data[0];
1439 if (nel==0 && val_nel==0)
return;
1442 casadi_assert(sz2==val_sz2 && sz1==val_sz1,
1443 "Sparsity::add<DataType>: shape mismatch. lhs is "
1444 +
dim() +
", while rhs is " + val_sp.dim() +
".");
1447 const casadi_int* c =
row();
1448 const casadi_int* rind =
colind();
1449 const casadi_int* v_c = val_sp.row();
1450 const casadi_int* v_rind = val_sp.colind();
1453 for (casadi_int i=0; i<sz2; ++i) {
1456 casadi_int v_el = v_rind[i];
1459 casadi_int v_el_end = v_rind[i+1];
1462 casadi_int v_j = v_el<v_el_end ? v_c[v_el] : sz1;
1465 for (casadi_int el=rind[i]; el!=rind[i+1]; ++el) {
1473 v_j = v_el<v_el_end ? v_c[v_el] : sz1;
1478 data[el] |= val_data[v_el++];
1479 v_j = v_el<v_el_end ? v_c[v_el] : sz1;
1490 typedef std::map<std::string, Sparsity>
SpDict;
Helper class for Serialization.
Helper class for Serialization.
SharedObject implements a reference counting framework similar for efficient and.
static Sparsity upper(casadi_int n)
Create a upper triangular square sparsity pattern *.
static Sparsity kkt(const Sparsity &H, const Sparsity &J, bool with_x_diag=true, bool with_lam_g_diag=true)
Get KKT system sparsity.
void get_nz(std::vector< casadi_int > &indices) const
Get the nonzero index for a set of elements.
casadi_int get_nz(casadi_int rr, casadi_int cc) const
Get the index of an existing non-zero element.
std::vector< casadi_int > find(bool ind1=SWIG_IND1) const
Get the location of all non-zero elements as they would appear in a Dense matrix.
casadi_int columns() const
Get the number of columns, Octave-style syntax.
casadi_int colind(casadi_int cc) const
Get a reference to the colindex of column cc (see class description)
Sparsity combine(const Sparsity &y, bool f0x_is_zero, bool function0_is_zero) const
Combine two sparsity patterns.
Sparsity pmult(const std::vector< casadi_int > &p, bool permute_rows=true, bool permute_columns=true, bool invert_permutation=false) const
Permute rows and/or columns.
Sparsity(casadi_int dummy=0)
Default constructor.
std::vector< casadi_int > get_upper() const
Get nonzeros in upper triangular part.
std::vector< casadi_int > erase(const std::vector< casadi_int > &rr, bool ind1=false)
Erase elements of a matrix.
bool is_subset(const Sparsity &rhs) const
Is subset?
static Sparsity band(casadi_int n, casadi_int p)
Create a single band in a square sparsity pattern.
std::vector< casadi_int > compress() const
Compress a sparsity pattern.
static Sparsity permutation(const std::vector< casadi_int > &p, bool invert=false)
Construct a permutation matrix P from a permutation vector p.
std::vector< casadi_int > etree(bool ata=false) const
Calculate the elimination tree.
bool is_stacked(const Sparsity &y, casadi_int n) const
Check if pattern is horizontal repeat of another.
Sparsity makeDense(std::vector< casadi_int > &mapping) const
Make a patten dense.
static bool test_cast(const SharedObjectInternal *ptr)
Check if a particular cast is allowed.
static std::string type_name()
Readable name of the public class.
bool is_vector() const
Check if the pattern is a row or column vector.
static Sparsity unit(casadi_int n, casadi_int el)
Create the sparsity pattern for a unit vector of length n and a nonzero on.
std::vector< casadi_int > get_lower() const
Get nonzeros in lower triangular part.
Sparsity sub(const std::vector< casadi_int > &rr, const std::vector< casadi_int > &cc, std::vector< casadi_int > &mapping, bool ind1=false) const
Get a submatrix.
casadi_int numel() const
The total number of elements, including structural zeros, i.e. size2()*size1()
casadi_int size1() const
Get the number of rows.
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.
casadi_int dfs(casadi_int j, casadi_int top, std::vector< casadi_int > &xi, std::vector< casadi_int > &pstack, const std::vector< casadi_int > &pinv, std::vector< bool > &marked) const
Depth-first search on the adjacency graph of the sparsity.
static Sparsity dense(casadi_int nrow, casadi_int ncol=1)
Create a dense rectangular sparsity pattern *.
Sparsity(casadi_int nrow, casadi_int ncol)
Pattern with all structural zeros.
Sparsity operator*(const Sparsity &b) const
Intersection of two sparsity patterns.
std::vector< casadi_int > largest_first() const
Order the columns by decreasing degree.
Sparsity star_coloring(casadi_int ordering=1, casadi_int cutoff=std::numeric_limits< casadi_int >::max()) const
Perform a star coloring of a symmetric matrix:
std::vector< casadi_int > get_colind() const
Get the column index for each column.
static Sparsity dense(const std::pair< casadi_int, casadi_int > &rc)
Create a dense rectangular sparsity pattern *.
Sparsity pattern_inverse() const
Take the inverse of a sparsity pattern; flip zeros and non-zeros.
static Sparsity diag(casadi_int nrow)
Create diagonal sparsity pattern *.
static Sparsity diag(const std::pair< casadi_int, casadi_int > &rc)
Create diagonal sparsity pattern *.
static Sparsity from_file(const std::string &filename, const std::string &format_hint="")
bool is_orthonormal(bool allow_empty=false) const
Are both rows and columns orthonormal ?
static Sparsity deserialize(DeserializingStream &s)
Deserialize.
casadi_int nnz_lower(bool strictly=false) const
Number of non-zeros in the lower triangular half,.
casadi_int scc(std::vector< casadi_int > &index, std::vector< casadi_int > &offset) const
Find the strongly connected components of the bigraph defined by the sparsity pattern.
Sparsity(casadi_int nrow, casadi_int ncol, const std::vector< casadi_int > &colind, const std::vector< casadi_int > &row, bool order_rows=false)
Construct from sparsity pattern vectors given in compressed column storage format.
std::string dim(bool with_nz=false) const
Get the dimension as a string.
std::vector< casadi_int > amd() const
Approximate minimal degree preordering.
Sparsity transpose(std::vector< casadi_int > &mapping, bool invert_mapping=false) const
Transpose the matrix and get the reordering of the non-zero entries.
Sparsity T() const
Transpose the matrix.
bool is_column() const
Check if the pattern is a column vector (i.e. size2()==1)
void spy(std::ostream &stream=casadi::uout()) const
Print a textual representation of sparsity.
bool is_reshape(const Sparsity &y) const
Check if the sparsity is a reshape of another.
void get_crs(std::vector< casadi_int > &rowind, std::vector< casadi_int > &col) const
Get the sparsity in compressed row storage (CRS) format.
void removeDuplicates(std::vector< casadi_int > &mapping)
Remove duplicate entries.
bool is_scalar(bool scalar_and_dense=false) const
Is scalar?
static Sparsity nonzeros(casadi_int nrow, casadi_int ncol, const std::vector< casadi_int > &nz, bool ind1=SWIG_IND1)
Create a sparsity from nonzeros.
void serialize(SerializingStream &s) const
Serialize an object.
Sparsity(const std::pair< casadi_int, casadi_int > &rc)
Create a sparse matrix with all structural zeros.
Sparsity unite(const Sparsity &y) const
Union of two sparsity patterns.
static Sparsity triplet(casadi_int nrow, casadi_int ncol, const std::vector< casadi_int > &row, const std::vector< casadi_int > &col, std::vector< casadi_int > &mapping, bool invert_mapping)
Create a sparsity pattern given the nonzeros in sparse triplet form *.
bool has_nz(casadi_int rr, casadi_int cc) const
Returns true if the pattern has a non-zero at location rr, cc.
bool is_orthonormal_rows(bool allow_empty=false) const
Are the rows of the pattern orthonormal ?
static Sparsity rowcol(const std::vector< casadi_int > &row, const std::vector< casadi_int > &col, casadi_int nrow, casadi_int ncol)
Construct a block sparsity pattern from (row, col) vectors.
void enlargeColumns(casadi_int ncol, const std::vector< casadi_int > &cc, bool ind1=false)
Enlarge the matrix along the second dimension (i.e. insert columns)
void append(const Sparsity &sp)
Append another sparsity patten vertically (NOTE: only efficient if vector)
casadi_int add_nz(casadi_int rr, casadi_int cc)
Get the index of a non-zero element.
bool is_row() const
Check if the pattern is a row vector (i.e. size1()==1)
const std::vector< casadi_int > permutation_vector(bool invert=false) const
Construct permutation vector from permutation matrix.
casadi_int size(casadi_int axis) const
Get the size along a particular dimensions.
Sparsity get_diag(std::vector< casadi_int > &mapping) const
std::vector< casadi_int > get_nz(const std::vector< casadi_int > &rr, const std::vector< casadi_int > &cc) const
Get a set of non-zero element.
Sparsity sub(const std::vector< casadi_int > &rr, const Sparsity &sp, std::vector< casadi_int > &mapping, bool ind1=false) const
Get a set of elements.
bool is_diag() const
Is diagonal?
Sparsity intersect(const Sparsity &y) const
Intersection of two sparsity patterns.
static Sparsity banded(casadi_int n, casadi_int p)
Create banded square sparsity pattern.
void enlargeRows(casadi_int nrow, const std::vector< casadi_int > &rr, bool ind1=false)
Enlarge the matrix along the first dimension (i.e. insert rows)
bool is_tril(bool strictly=false) const
Is lower triangular?
casadi_int rows() const
Get the number of rows, Octave-style syntax.
void spy_matlab(const std::string &mfile) const
Generate a script for Matlab or Octave which visualizes.
casadi_int nnz_upper(bool strictly=false) const
Number of non-zeros in the upper triangular half,.
casadi_int nnz() const
Get the number of (structural) non-zeros.
std::vector< casadi_int > erase(const std::vector< casadi_int > &rr, const std::vector< casadi_int > &cc, bool ind1=false)
Erase rows and/or columns of a matrix.
casadi_int size2() const
Get the number of columns.
casadi_int btf(std::vector< casadi_int > &rowperm, std::vector< casadi_int > &colperm, std::vector< casadi_int > &rowblock, std::vector< casadi_int > &colblock, std::vector< casadi_int > &coarse_rowblock, std::vector< casadi_int > &coarse_colblock) const
Calculate the block triangular form (BTF)
std::string serialize() const
Serialize.
bool is_equal(casadi_int nrow, casadi_int ncol, const std::vector< casadi_int > &colind, const std::vector< casadi_int > &row) const
bool is_selection(bool allow_empty=false) const
Is this a selection matrix?
Sparsity sparsity_cast_mod(const Sparsity &X, const Sparsity &Y) const
Propagates subset according to sparsity cast.
bool is_permutation() const
Is this a permutation matrix?
casadi_int row(casadi_int el) const
Get the row of a non-zero element.
Sparsity ldl(std::vector< casadi_int > &p, bool amd=true) const
Symbolic LDL factorization.
static Sparsity scalar(bool dense_scalar=true)
Create a scalar sparsity pattern *.
bool is_empty(bool both=false) const
Check if the sparsity is empty.
casadi_int bw_upper() const
Upper half-bandwidth.
bool is_singular() const
Check whether the sparsity-pattern indicates structural singularity.
void export_code(const std::string &lang, std::ostream &stream=casadi::uout(), const Dict &options=Dict()) const
Export matrix in specific language.
void to_file(const std::string &filename, const std::string &format_hint="") const
bool operator!=(const Sparsity &y) const
Check if two sparsity patterns are difference.
static Sparsity triplet(casadi_int nrow, casadi_int ncol, const std::vector< casadi_int > &row, const std::vector< casadi_int > &col)
Create a sparsity pattern given the nonzeros in sparse triplet form *.
static Sparsity lower(casadi_int n)
Create a lower triangular square sparsity pattern *.
void get_triplet(std::vector< casadi_int > &row, std::vector< casadi_int > &col) const
Get the sparsity in sparse triplet format.
double density() const
The percentage of nonzero.
Sparsity uni_coloring(const Sparsity &AT=Sparsity(), casadi_int cutoff=std::numeric_limits< casadi_int >::max()) const
Perform a unidirectional coloring: A greedy distance-2 coloring algorithm.
Sparsity operator+(const Sparsity &b) const
Union of two sparsity patterns.
casadi_int nnz_diag() const
Number of non-zeros on the diagonal, i.e. the number of elements (i, j) with j==i.
bool is_transpose(const Sparsity &y) const
Check if the sparsity is the transpose of another.
static Sparsity diag(casadi_int nrow, casadi_int ncol)
Create diagonal sparsity pattern *.
void appendColumns(const Sparsity &sp)
Append another sparsity patten horizontally.
void get_ccs(std::vector< casadi_int > &colind, std::vector< casadi_int > &row) const
Get the sparsity in compressed column storage (CCS) format.
bool is_dense() const
Is dense?
bool is_orthonormal_columns(bool allow_empty=false) const
Are the columns of the pattern orthonormal ?
static Sparsity compressed(const std::vector< casadi_int > &v, bool order_rows=false)
void resize(casadi_int nrow, casadi_int ncol)
Resize.
std::string postfix_dim() const
Dimension string as a postfix to a name.
bool operator==(const Sparsity &y) const
static Sparsity deserialize(std::istream &stream)
Build Sparsity from serialization.
bool is_square() const
Is square?
void qr_sparse(Sparsity &V, Sparsity &R, std::vector< casadi_int > &prinv, std::vector< casadi_int > &pc, bool amd=true) const
Symbolic QR factorization.
bool rowsSequential(bool strictly=true) const
Do the rows appear sequentially on each column.
bool is_triu(bool strictly=false) const
Is upper triangular?
std::pair< casadi_int, casadi_int > size() const
Get the shape.
std::vector< casadi_int > get_row() const
Get the row for each non-zero entry.
std::string repr_el(casadi_int k) const
Describe the nonzero location k as a string.
Sparsity star_coloring2(casadi_int ordering=1, casadi_int cutoff=std::numeric_limits< casadi_int >::max()) const
Perform a star coloring of a symmetric matrix:
bool is_symmetric() const
Is symmetric?
casadi_int bw_lower() const
Lower half-bandwidth.
std::vector< casadi_int > get_col() const
Get the column for each non-zero entry.
static Sparsity deserialize(const std::string &s)
Build Sparsity from serialization.
CASADI_EXPORT std::ostream & uout()
GenericType::Dict Dict
C++ equivalent of Python's dict or MATLAB's struct.
void hash_combine(std::size_t &seed, T v)
Generate a hash value incrementally (function taken from boost)
CASADI_EXPORT std::size_t hash_sparsity(casadi_int nrow, casadi_int ncol, const std::vector< casadi_int > &colind, const std::vector< casadi_int > &row)
Hash a sparsity pattern.
size_t hash_value(T v)
Hash value of an integer.
std::map< std::string, Sparsity > SpDict