sparsity.hpp
1 /*
2  * This file is part of CasADi.
3  *
4  * CasADi -- A symbolic framework for dynamic optimization.
5  * Copyright (C) 2010-2023 Joel Andersson, Joris Gillis, Moritz Diehl,
6  * KU Leuven. All rights reserved.
7  * Copyright (C) 2011-2014 Greg Horn
8  * Copyright (C) 2005-2013 Timothy A. Davis
9  *
10  * CasADi is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public
12  * License as published by the Free Software Foundation; either
13  * version 3 of the License, or (at your option) any later version.
14  *
15  * CasADi is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with CasADi; if not, write to the Free Software
22  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
23  *
24  */
25 
26 
27 #ifndef CASADI_SPARSITY_HPP
28 #define CASADI_SPARSITY_HPP
29 
30 #include "shared_object.hpp"
31 #include "printable.hpp"
32 #include "casadi_common.hpp"
33 #include "sparsity_interface.hpp"
34 #include "generic_type.hpp"
35 #include <vector>
36 #include <list>
37 #include <limits>
38 #include <unordered_map>
39 
40 namespace casadi {
41  // Forward declaration
42  class SparsityInternal;
43  class SerializingStream;
44  class DeserializingStream;
45 
46  #ifndef SWIG
50  struct CASADI_EXPORT SparsityStruct {
51  casadi_int nrow;
52  casadi_int ncol;
53  const casadi_int* colind;
54  const casadi_int* row;
55  };
56  #endif // SWIG
57 
96  class CASADI_EXPORT Sparsity
97  : public SharedObject,
98  public SWIG_IF_ELSE(SparsityInterfaceCommon, SparsityInterface<Sparsity>),
99  public SWIG_IF_ELSE(PrintableCommon, Printable<Sparsity>) {
100  public:
101 
103  explicit Sparsity(casadi_int dummy=0);
104 
108  Sparsity(casadi_int nrow, casadi_int ncol);
109 
111  Sparsity(casadi_int nrow, casadi_int ncol,
112  const std::vector<casadi_int>& colind, const std::vector<casadi_int>& row,
113  bool order_rows=false);
114 
118  explicit Sparsity(const std::pair<casadi_int, casadi_int>& rc);
119 
120 #ifndef SWIG
122  Sparsity(casadi_int nrow, casadi_int ncol, const casadi_int* colind, const casadi_int* row,
123  bool order_rows=false);
124 
128  static Sparsity create(SparsityInternal *node);
129 
131  typedef SparsityInterface<Sparsity> B;
132 
134  using B::horzsplit;
135  using B::diagsplit;
136  using B::vertsplit;
137  using B::mtimes;
138 
139  SparsityInternal* get() const;
140 #endif
141 
146  static Sparsity scalar(bool dense_scalar=true)
147  { return dense_scalar ? dense(1, 1) : Sparsity(1, 1); }
149 
154  static Sparsity dense(casadi_int nrow, casadi_int ncol=1);
155  static Sparsity dense(const std::pair<casadi_int, casadi_int> &rc) {
156  return dense(rc.first, rc.second);
157  }
159 
166  static Sparsity unit(casadi_int n, casadi_int el);
168 
172  static Sparsity upper(casadi_int n);
173 
177  static Sparsity lower(casadi_int n);
178 
183  static Sparsity diag(casadi_int nrow) { return diag(nrow, nrow);}
184  static Sparsity diag(casadi_int nrow, casadi_int ncol);
185  static Sparsity diag(const std::pair<casadi_int, casadi_int> &rc) {
186  return diag(rc.first, rc.second);
187  }
189 
197  static Sparsity band(casadi_int n, casadi_int p);
198 
205  static Sparsity banded(casadi_int n, casadi_int p);
206 
210  static Sparsity rowcol(const std::vector<casadi_int>& row,
211  const std::vector<casadi_int>& col,
212  casadi_int nrow, casadi_int ncol);
213 
215 
218  static Sparsity triplet(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);
224 
230  static Sparsity nonzeros(casadi_int nrow, casadi_int ncol, const std::vector<casadi_int>& nz,
231  bool ind1=SWIG_IND1);
232 
241  static Sparsity compressed(const std::vector<casadi_int>& v, bool order_rows=false);
242 #ifndef SWIG
243  static Sparsity compressed(const casadi_int* v, bool order_rows=false);
244 #endif // SWIG
246 
257  static Sparsity permutation(const std::vector<casadi_int>& p, bool invert=false);
258 
262  const std::vector<casadi_int> permutation_vector(bool invert=false) const;
263 
269  Sparsity get_diag(std::vector<casadi_int>& SWIG_OUTPUT(mapping)) const;
270 
272  std::vector<casadi_int> compress() const;
273 
274 #ifndef SWIG
276  const SparsityInternal* operator->() const;
277 
279  const SparsityInternal& operator*() const;
280 #endif // SWIG
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;
286 #ifndef SWIG
287  bool is_equal(casadi_int nrow, casadi_int ncol,
288  const casadi_int* colind, const casadi_int* row) const;
289 #endif // SWIG
290 
291  bool operator==(const Sparsity& y) const { return is_equal(y);}
293 
295  bool operator!=(const Sparsity& y) const {return !is_equal(y);}
296 
298  bool is_stacked(const Sparsity& y, casadi_int n) const;
299 
300 #ifndef SWIG
307  operator const casadi_int*() const;
308 
312  operator const std::vector<casadi_int>&() const;
313 
317  operator SparsityStruct() const;
318 #endif // SWIG
319 
322 
324  casadi_int size1() const;
325 
327  casadi_int rows() const {return size1();}
328 
330  casadi_int size2() const;
331 
333  casadi_int columns() const {return size2();}
334 
341  casadi_int numel() const;
342 
348  double density() const;
349 
356  bool is_empty(bool both=false) const;
357 
363  casadi_int nnz() const;
364 
370  casadi_int nnz_upper(bool strictly=false) const;
371 
377  casadi_int nnz_lower(bool strictly=false) const;
378 
382  casadi_int nnz_diag() const;
383 
387  casadi_int bw_upper() const;
388 
392  casadi_int bw_lower() const;
393 
397  std::pair<casadi_int, casadi_int> size() const;
398 
402  casadi_int size(casadi_int axis) const;
404 
406  Dict info() const;
407 
413  void to_file(const std::string& filename, const std::string& format_hint="") const;
414 
415  static Sparsity from_file(const std::string& filename, const std::string& format_hint="");
416 
417 #ifndef SWIG
421  void serialize(std::ostream &stream) const;
422 #endif
423 
427  std::string serialize() const;
428 
432  static Sparsity deserialize(std::istream& stream);
433 
437  static Sparsity deserialize(const std::string& s);
438 
442  void serialize(SerializingStream& s) const;
443 
448 
449 #ifndef SWIG
455  const casadi_int* row() const;
456 
460  const casadi_int* colind() const;
461 #endif
462 
470  std::vector<casadi_int> get_row() const;
471 
478  std::vector<casadi_int> get_colind() const;
479 
483  casadi_int colind(casadi_int cc) const;
484 
488  casadi_int row(casadi_int el) const;
489 
496  std::vector<casadi_int> get_col() const;
497 
499  void resize(casadi_int nrow, casadi_int ncol);
500 
506  casadi_int add_nz(casadi_int rr, casadi_int cc);
507 
513  casadi_int get_nz(casadi_int rr, casadi_int cc) const;
514 
516  bool has_nz(casadi_int rr, casadi_int cc) const;
517 
523  std::vector<casadi_int> get_nz(const std::vector<casadi_int>& rr,
524  const std::vector<casadi_int>& cc) const;
525 
533  void get_nz(std::vector<casadi_int>& SWIG_INOUT(indices)) const;
534 
536  std::vector<casadi_int> get_lower() const;
537 
539  std::vector<casadi_int> get_upper() const;
540 
542  void get_ccs(std::vector<casadi_int>& SWIG_OUTPUT(colind),
543  std::vector<casadi_int>& SWIG_OUTPUT(row)) const;
544 
546  void get_crs(std::vector<casadi_int>& SWIG_OUTPUT(rowind),
547  std::vector<casadi_int>& SWIG_OUTPUT(col)) const;
548 
550  void get_triplet(std::vector<casadi_int>& SWIG_OUTPUT(row),
551  std::vector<casadi_int>& SWIG_OUTPUT(col)) const;
552 
559  Sparsity sub(const std::vector<casadi_int>& rr,
560  const std::vector<casadi_int>& cc,
561  std::vector<casadi_int>& SWIG_OUTPUT(mapping), bool ind1=false) const;
562 
569  Sparsity sub(const std::vector<casadi_int>& rr, const Sparsity& sp,
570  std::vector<casadi_int>& SWIG_OUTPUT(mapping), bool ind1=false) const;
571 
573  Sparsity T() const;
574 
581  Sparsity transpose(std::vector<casadi_int>& SWIG_OUTPUT(mapping),
582  bool invert_mapping=false) const;
583 
585  bool is_transpose(const Sparsity& y) const;
586 
588  bool is_reshape(const Sparsity& y) const;
589 
591 
601 #ifndef SWIG
602  Sparsity combine(const Sparsity& y, bool f0x_is_zero, bool function0_is_zero,
603  std::vector<unsigned char>& mapping) const;
604 #endif // SWIG
605  Sparsity combine(const Sparsity& y, bool f0x_is_zero, bool function0_is_zero) const;
607 
609 
612 #ifndef SWIG
613  Sparsity unite(const Sparsity& y, std::vector<unsigned char>& mapping) const;
614 #endif // SWIG
615  Sparsity unite(const Sparsity& y) const;
616  Sparsity operator+(const Sparsity& b) const;
618 
620 
628 #ifndef SWIG
629  Sparsity intersect(const Sparsity& y,
630  std::vector<unsigned char>& mapping) const;
631 #endif // SWIG
632  Sparsity intersect(const Sparsity& y) const;
633  Sparsity operator*(const Sparsity& b) const;
635 
637  bool is_subset(const Sparsity& rhs) const;
638 
652  Sparsity sparsity_cast_mod(const Sparsity& X, const Sparsity& Y) const;
653 
656 
657 #ifndef SWIG
664  static void mul_sparsityF(const bvec_t* x, const Sparsity& x_sp,
665  const bvec_t* y, const Sparsity& y_sp,
666  bvec_t* z, const Sparsity& z_sp,
667  bvec_t* w);
668 
675  static void mul_sparsityR(bvec_t* x, const Sparsity& x_sp,
676  bvec_t* y, const Sparsity& y_sp,
677  bvec_t* z, const Sparsity& z_sp,
678  bvec_t* w);
679 
682 
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>
694  diagsplit(const Sparsity& x,
695  const std::vector<casadi_int>& offset1,
696  const std::vector<casadi_int>& offset2);
697  static Sparsity mtimes(const Sparsity& x, const Sparsity& y);
698  static Sparsity mac(const Sparsity& x, const Sparsity& y, const Sparsity& z) { return z;}
699  static Sparsity reshape(const Sparsity& x, casadi_int nrow, casadi_int ncol);
700  static Sparsity reshape(const Sparsity& x, const Sparsity& sp);
701  static Sparsity sparsity_cast(const Sparsity& x, const Sparsity& sp);
702  static casadi_int sprank(const Sparsity& x);
703  static casadi_int norm_0_mul(const Sparsity& x, const Sparsity& A);
704  static Sparsity kron(const Sparsity& a, const Sparsity& b);
705  static Sparsity triu(const Sparsity& x, bool includeDiagonal=true);
706  static Sparsity tril(const Sparsity& x, bool includeDiagonal=true);
707 
708  static Sparsity sum2(const Sparsity &x);
709  static Sparsity sum1(const Sparsity &x);
710 #endif //SWIG
711 
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);
729 
733  void enlargeRows(casadi_int nrow, const std::vector<casadi_int>& rr, bool ind1=false);
734 
738  void enlargeColumns(casadi_int ncol, const std::vector<casadi_int>& cc, bool ind1=false);
739 
743  Sparsity makeDense(std::vector<casadi_int>& SWIG_OUTPUT(mapping)) const;
744 
748  std::vector<casadi_int> erase(const std::vector<casadi_int>& rr,
749  const std::vector<casadi_int>& cc, bool ind1=false);
750 
754  std::vector<casadi_int> erase(const std::vector<casadi_int>& rr, bool ind1=false);
755 
757  void append(const Sparsity& sp);
758 
760  void appendColumns(const Sparsity& sp);
761 
763  bool is_scalar(bool scalar_and_dense=false) const;
764 
766  bool is_dense() const;
767 
771  bool is_row() const;
772 
776  bool is_column() const;
777 
781  bool is_vector() const;
782 
786  bool is_diag() const;
787 
791  bool is_square() const;
792 
796  bool is_symmetric() const;
797 
801  bool is_triu(bool strictly = false) const;
802 
806  bool is_tril(bool strictly = false) const;
807 
811  bool is_singular() const;
812 
823  bool is_permutation() const;
824 
835  bool is_selection(bool allow_empty=false) const;
836 
842  bool is_orthonormal(bool allow_empty=false) const;
843 
849  bool is_orthonormal_rows(bool allow_empty=false) const;
850 
856  bool is_orthonormal_columns(bool allow_empty=false) const;
857 
863  bool rowsSequential(bool strictly=true) const;
864 
871  void removeDuplicates(std::vector<casadi_int>& SWIG_INOUT(mapping));
872 
873 #ifndef SWIG
874  typedef std::unordered_multimap<std::size_t, WeakRef> CachingMap;
875 
877  static CachingMap& getCache();
878 
880  static const Sparsity& getScalar();
881 
883  static const Sparsity& getScalarSparse();
884 
886  static const Sparsity& getEmpty();
887 
888 #endif //SWIG
889 
902  std::vector<casadi_int> etree(bool ata=false) const;
903 
913  Sparsity ldl(std::vector<casadi_int>& SWIG_OUTPUT(p), bool amd=true) const;
914 
924  void qr_sparse(Sparsity& SWIG_OUTPUT(V), Sparsity& SWIG_OUTPUT(R),
925  std::vector<casadi_int>& SWIG_OUTPUT(prinv),
926  std::vector<casadi_int>& SWIG_OUTPUT(pc), bool amd=true) const;
927 
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;
936 
960  casadi_int scc(std::vector<casadi_int>& SWIG_OUTPUT(index),
961  std::vector<casadi_int>& SWIG_OUTPUT(offset)) const;
962 
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;
990 
1003  std::vector<casadi_int> amd() const;
1004 
1005 #ifndef SWIG
1009  void spsolve(bvec_t* X, bvec_t* B, bool tr) const;
1010 #endif // SWIG
1011 
1023  std::vector<casadi_int> find(bool ind1=SWIG_IND1) const;
1024 
1025 #ifndef SWIG
1027  void find(std::vector<casadi_int>& loc, bool ind1=false) const;
1028 #endif // SWIG
1029 
1036  casadi_int cutoff = std::numeric_limits<casadi_int>::max()) const;
1037 
1049  Sparsity star_coloring(casadi_int ordering = 1,
1050  casadi_int cutoff = std::numeric_limits<casadi_int>::max()) const;
1051 
1063  Sparsity star_coloring2(casadi_int ordering = 1,
1064  casadi_int cutoff = std::numeric_limits<casadi_int>::max()) const;
1065 
1069  std::vector<casadi_int> largest_first() const;
1070 
1078  Sparsity pmult(const std::vector<casadi_int>& p,
1079  bool permute_rows=true, bool permute_columns=true,
1080  bool invert_permutation=false) const;
1081 
1083  std::string dim(bool with_nz=false) const;
1084 
1095  std::string postfix_dim() const;
1096 
1098  std::string repr_el(casadi_int k) const;
1099 
1103  void spy(std::ostream &stream=casadi::uout()) const;
1104 
1110  void spy_matlab(const std::string& mfile) const;
1111 
1124  void export_code(const std::string& lang, std::ostream &stream=casadi::uout(),
1125  const Dict& options=Dict()) const;
1126 
1128  static std::string type_name() {return "Sparsity";}
1129 
1130  // Hash the sparsity pattern
1131  std::size_t hash() const;
1132 
1134  static bool test_cast(const SharedObjectInternal* ptr);
1135 
1141  static Sparsity kkt(const Sparsity& H, const Sparsity& J,
1142  bool with_x_diag=true, bool with_lam_g_diag=true);
1143 
1144 #ifndef SWIG
1150  template<typename T>
1151  void set(T* data, const T* val_data, const Sparsity& val_sp) const;
1152 
1158  template<typename T>
1159  void add(T* data, const T* val_data, const Sparsity& val_sp) const;
1160 
1166  template<typename T>
1167  void bor(T* data, const T* val_data, const Sparsity& val_sp) const;
1168 
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;
1172  private:
1173 
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);
1177 
1179  void assign_cached(casadi_int nrow, casadi_int ncol,
1180  const casadi_int* colind, const casadi_int* row, bool order_rows=false);
1181 
1182 #endif //SWIG
1183  };
1184 
1188  template<typename T>
1189  inline size_t hash_value(T v) { return size_t(v);}
1190 
1194  template<typename T>
1195  inline void hash_combine(std::size_t& seed, T v) {
1196  seed ^= hash_value(v) + 0x9e3779b9 + (seed << 6) + (seed >> 2);
1197  }
1198 
1202  template<typename T>
1203  inline void hash_combine(std::size_t& seed, const T* v, std::size_t sz) {
1204  for (casadi_int i=0; i<sz; ++i) hash_combine(seed, v[i]);
1205  }
1206 
1210  template<typename T>
1211  inline void hash_combine(std::size_t& seed, const std::vector<T>& v) {
1212  hash_combine(seed, get_ptr(v), v.size());
1213  }
1214 
1215  template<>
1216  inline size_t hash_value(std::string v) {
1217  size_t seed = 0;
1218  hash_combine(seed, v.c_str(), v.size());
1219  return seed;
1220  }
1221 
1225  CASADI_EXPORT std::size_t hash_sparsity(casadi_int nrow, casadi_int ncol,
1226  const std::vector<casadi_int>& colind,
1227  const std::vector<casadi_int>& row);
1228 
1229  CASADI_EXPORT std::size_t hash_sparsity(casadi_int nrow, casadi_int ncol,
1230  const casadi_int* colind,
1231  const casadi_int* row);
1232 
1233 #ifndef SWIG
1234  // Template instantiations
1235  template<typename DataType>
1236  void Sparsity::set(DataType* data, const DataType* val_data, const Sparsity& val_sp) const {
1237  // Get dimensions of this
1238  const casadi_int sz = nnz();
1239  const casadi_int sz1 = size1();
1240  const casadi_int sz2 = size2();
1241 
1242  // Get dimensions of assigning matrix
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;
1247 
1248  // Check if sparsity matches
1249  if (val_sp==*this) {
1250  std::copy(val_data, val_data+sz, data);
1251  } else if (this->is_empty()) {
1252  // Quick return
1253  return;
1254  } else if (val_sp.is_empty()) {
1255  // Quick return
1256  return;
1257  } else if (val_nel==1) { // if scalar
1258  std::fill(data, data+sz, val_sz==0 ? DataType(0) : val_data[0]);
1259  } else if (sz2==val_sz2 && sz1==val_sz1) {
1260  // Matching dimensions
1261  // Sparsity
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();
1266 
1267  // For all columns
1268  for (casadi_int i=0; i<sz2; ++i) {
1269 
1270  // Nonzero of the assigning matrix
1271  casadi_int v_el = v_rind[i];
1272 
1273  // First nonzero of the following column
1274  casadi_int v_el_end = v_rind[i+1];
1275 
1276  // Next row of the assigning matrix
1277  casadi_int v_j = v_el<v_el_end ? v_c[v_el] : sz1;
1278 
1279  // Assign all nonzeros
1280  for (casadi_int el=rind[i]; el!=rind[i+1]; ++el) {
1281 
1282  // Get row
1283  casadi_int j=c[el];
1284 
1285  // Forward the assigning nonzero
1286  while (v_j<j) {
1287  v_el++;
1288  v_j = v_el<v_el_end ? v_c[v_el] : sz1;
1289  }
1290 
1291  // Assign nonzero
1292  if (v_j==j) {
1293  data[el] = val_data[v_el++];
1294  v_j = v_el<v_el_end ? v_c[v_el] : sz1;
1295  } else {
1296  data[el] = 0;
1297  }
1298  }
1299  }
1300  } else if (sz1==val_sz2 && sz2==val_sz1 && sz2 == 1) {
1301  // Assign transposed (this is column)
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]];
1307  }
1308  } else if (sz1==val_sz2 && sz2==val_sz1 && sz1 == 1) {
1309  // Assign transposed (this is row)
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];
1317  }
1318  }
1319  } else {
1320  // Make sure that dimension matches
1321  casadi_error("Sparsity::set<DataType>: shape mismatch. lhs is "
1322  + dim() + ", while rhs is " + val_sp.dim() + ".");
1323  }
1324  }
1325 
1326  template<typename DataType>
1327  void Sparsity::add(DataType* data, const DataType* val_data, const Sparsity& val_sp) const {
1328  // Get dimensions of this
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;
1333 
1334  // Get dimensions of assigning matrix
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;
1339 
1340  // Check if sparsity matches
1341  if (val_sp==*this) {
1342  for (casadi_int k=0; k<sz; ++k) {
1343  data[k] += val_data[k];
1344  }
1345  } else if (this->is_empty()) {
1346  // Quick return
1347  return;
1348  } else if (val_sp.is_empty()) {
1349  // Quick return
1350  return;
1351  } else if (val_nel==1) { // if scalar
1352  if (val_sz!=0) {
1353  for (casadi_int k=0; k<sz; ++k) {
1354  data[k] += val_data[0];
1355  }
1356  }
1357  } else {
1358  // Quick return if empty
1359  if (nel==0 && val_nel==0) return;
1360 
1361  // Make sure that dimension matches
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() + ".");
1365 
1366  // Sparsity
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();
1371 
1372  // For all columns
1373  for (casadi_int i=0; i<sz2; ++i) {
1374 
1375  // Nonzero of the assigning matrix
1376  casadi_int v_el = v_rind[i];
1377 
1378  // First nonzero of the following column
1379  casadi_int v_el_end = v_rind[i+1];
1380 
1381  // Next row of the assigning matrix
1382  casadi_int v_j = v_el<v_el_end ? v_c[v_el] : sz1;
1383 
1384  // Assign all nonzeros
1385  for (casadi_int el=rind[i]; el!=rind[i+1]; ++el) {
1386 
1387  // Get row
1388  casadi_int j=c[el];
1389 
1390  // Forward the assigning nonzero
1391  while (v_j<j) {
1392  v_el++;
1393  v_j = v_el<v_el_end ? v_c[v_el] : sz1;
1394  }
1395 
1396  // Assign nonzero
1397  if (v_j==j) {
1398  data[el] += val_data[v_el++];
1399  v_j = v_el<v_el_end ? v_c[v_el] : sz1;
1400  }
1401  }
1402  }
1403  }
1404  }
1405 
1406  template<typename DataType>
1407  void Sparsity::bor(DataType* data, const DataType* val_data, const Sparsity& val_sp) const {
1408  // Get dimensions of this
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;
1413 
1414  // Get dimensions of assigning matrix
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;
1419 
1420  // Check if sparsity matches
1421  if (val_sp==*this) {
1422  for (casadi_int k=0; k<sz; ++k) {
1423  data[k] |= val_data[k];
1424  }
1425  } else if (this->is_empty()) {
1426  // Quick return
1427  return;
1428  } else if (val_sp.is_empty()) {
1429  // Quick return
1430  return;
1431  } else if (val_nel==1) { // if scalar
1432  if (val_sz!=0) {
1433  for (casadi_int k=0; k<sz; ++k) {
1434  data[k] |= val_data[0];
1435  }
1436  }
1437  } else {
1438  // Quick return if empty
1439  if (nel==0 && val_nel==0) return;
1440 
1441  // Make sure that dimension matches
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() + ".");
1445 
1446  // Sparsity
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();
1451 
1452  // For all columns
1453  for (casadi_int i=0; i<sz2; ++i) {
1454 
1455  // Nonzero of the assigning matrix
1456  casadi_int v_el = v_rind[i];
1457 
1458  // First nonzero of the following column
1459  casadi_int v_el_end = v_rind[i+1];
1460 
1461  // Next row of the assigning matrix
1462  casadi_int v_j = v_el<v_el_end ? v_c[v_el] : sz1;
1463 
1464  // Assign all nonzeros
1465  for (casadi_int el=rind[i]; el!=rind[i+1]; ++el) {
1466 
1467  // Get row
1468  casadi_int j=c[el];
1469 
1470  // Forward the assigning nonzero
1471  while (v_j<j) {
1472  v_el++;
1473  v_j = v_el<v_el_end ? v_c[v_el] : sz1;
1474  }
1475 
1476  // Assign nonzero
1477  if (v_j==j) {
1478  data[el] |= val_data[v_el++];
1479  v_j = v_el<v_el_end ? v_c[v_el] : sz1;
1480  }
1481  }
1482  }
1483  }
1484  }
1485 
1486 #endif //SWIG
1487 
1490  typedef std::map<std::string, Sparsity> SpDict;
1492 
1493 } // namespace casadi
1494 
1495 #endif // CASADI_SPARSITY_HPP
Helper class for Serialization.
Helper class for Serialization.
SharedObject implements a reference counting framework similar for efficient and.
General sparsity class.
Definition: sparsity.hpp:99
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.
Definition: sparsity.hpp:333
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.
Definition: sparsity.hpp:1128
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.
Dict info() const
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 *.
Definition: sparsity.hpp:155
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 *.
Definition: sparsity.hpp:183
static Sparsity diag(const std::pair< casadi_int, casadi_int > &rc)
Create diagonal sparsity pattern *.
Definition: sparsity.hpp:185
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.
Definition: sparsity.hpp:327
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 *.
Definition: sparsity.hpp:146
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.
Definition: sparsity.hpp:295
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.
std::size_t hash() const
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
Definition: sparsity.hpp:291
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.
The casadi namespace.
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)
Definition: sparsity.hpp:1195
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.
Definition: sparsity.hpp:1189
std::map< std::string, Sparsity > SpDict
Definition: sparsity.hpp:1490