sparsity_internal.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_INTERNAL_HPP
28 #define CASADI_SPARSITY_INTERNAL_HPP
29 
30 #include "sparsity.hpp"
31 #include "shared_object.hpp"
32 #ifdef CASADI_WITH_THREAD
33 #ifdef CASADI_WITH_THREAD_MINGW
34 #include <mingw.mutex.h>
35 #else // CASADI_WITH_THREAD_MINGW
36 #include <mutex>
37 #endif // CASADI_WITH_THREAD_MINGW
38 #endif //CASADI_WITH_THREAD
39 
41 
42 namespace casadi {
43 
44  class CASADI_EXPORT SparsityInternal : public SharedObjectInternal {
45  private:
55  std::vector<casadi_int> sp_;
56 
60  struct Btf {
61  casadi_int nb;
62  std::vector<casadi_int> rowperm, colperm;
63  std::vector<casadi_int> rowblock, colblock;
64  std::vector<casadi_int> coarse_rowblock, coarse_colblock;
65  };
66 
72  mutable Btf* btf_;
73 
74 #ifdef CASADI_WITH_THREADSAFE_SYMBOLICS
76  mutable std::mutex btf_mtx_;
77 #endif // CASADI_WITH_THREADSAFE_SYMBOLICS
78 
79  public:
81  SparsityInternal(casadi_int nrow, casadi_int ncol,
82  const casadi_int* colind, const casadi_int* row);
83 
85  ~SparsityInternal() override;
86 
90  inline const std::vector<casadi_int>& sp() const { return sp_;}
91 
95  inline casadi_int size1() const { return sp_[0];}
96 
100  inline casadi_int size2() const { return sp_[1];}
101 
105  inline const casadi_int* colind() const { return &sp_.front()+2;}
106 
110  inline const casadi_int* row() const { return colind()+size2()+1;}
111 
113  inline casadi_int nnz() const { return colind()[size2()];}
114 
120  Sparsity get_diag(std::vector<casadi_int>& mapping) const;
121 
123  bool has_diag() const;
124 
126  Sparsity drop_diag() const;
127 
135  casadi_int dfs(casadi_int j, casadi_int top, std::vector<casadi_int>& xi,
136  std::vector<casadi_int>& pstack,
137  const std::vector<casadi_int>& pinv, std::vector<bool>& marked) const;
138 
146  casadi_int scc(std::vector<casadi_int>& p, std::vector<casadi_int>& r) const;
147 
155  std::vector<casadi_int> amd() const;
156 
167  static void etree(const casadi_int* sp, casadi_int* parent, casadi_int *w, casadi_int ata);
168 
177  static casadi_int postorder_dfs(casadi_int j, casadi_int k, casadi_int* head,
178  const casadi_int* next, casadi_int* post, casadi_int* stack);
179 
190  static void postorder(const casadi_int* parent, casadi_int n, casadi_int* post, casadi_int* w);
191 
200  static casadi_int leaf(casadi_int i, casadi_int j, const casadi_int* first,
201  casadi_int* maxfirst,
202  casadi_int* prevleaf, casadi_int* ancestor, casadi_int* jleaf);
203 
214  static casadi_int qr_counts(const casadi_int* tr_sp, const casadi_int* parent,
215  const casadi_int* post, casadi_int* counts, casadi_int* w);
216 
228  static casadi_int qr_nnz(const casadi_int* sp, casadi_int* pinv, casadi_int* leftmost,
229  const casadi_int* parent, casadi_int* nrow_ext, casadi_int* w);
230 
239  static void qr_init(const casadi_int* sp, const casadi_int* sp_tr,
240  casadi_int* leftmost, casadi_int* parent, casadi_int* pinv,
241  casadi_int* nrow_ext, casadi_int* v_nnz, casadi_int* r_nnz, casadi_int* w);
242 
259  static void qr_sparsities(
260  const casadi_int* sp_a, casadi_int nrow_ext, casadi_int* sp_v, casadi_int* sp_r,
261  const casadi_int* leftmost, const casadi_int* parent, const casadi_int* pinv,
262  casadi_int* iw);
263 
276  static void ldl_colind(const casadi_int* sp, casadi_int* parent,
277  casadi_int* l_colind, casadi_int* w);
278 
289  static void ldl_row(const casadi_int* sp, const casadi_int* parent,
290  casadi_int* l_colind, casadi_int* l_row, casadi_int *w);
291 
293  Sparsity T() const;
294 
300  Sparsity transpose(std::vector<casadi_int>& mapping, bool invert_mapping=false) const;
301 
303  bool is_transpose(const SparsityInternal& y) const;
304 
306  bool is_reshape(const SparsityInternal& y) const;
307 
315  void bfs(casadi_int n, std::vector<casadi_int>& wi, std::vector<casadi_int>& wj,
316  std::vector<casadi_int>& queue, const std::vector<casadi_int>& imatch,
317  const std::vector<casadi_int>& jmatch, casadi_int mark) const;
318 
326  static void matched(
327  casadi_int n, const std::vector<casadi_int>& wj, const std::vector<casadi_int>& imatch,
328  std::vector<casadi_int>& p, std::vector<casadi_int>& q, std::vector<casadi_int>& cc,
329  std::vector<casadi_int>& rr, casadi_int set, casadi_int mark);
330 
338  static void unmatched(casadi_int m, const std::vector<casadi_int>& wi,
339  std::vector<casadi_int>& p, std::vector<casadi_int>& rr, casadi_int set);
340 
348  static casadi_int rprune(casadi_int i, casadi_int j, double aij, void *other);
349 
357  static casadi_int drop(casadi_int (*fkeep)(casadi_int, casadi_int, double, void *), void *other,
358  casadi_int nrow, casadi_int ncol,
359  std::vector<casadi_int>& colind, std::vector<casadi_int>& row);
360 
362  casadi_int btf(std::vector<casadi_int>& rowperm, std::vector<casadi_int>& colperm,
363  std::vector<casadi_int>& rowblock, std::vector<casadi_int>& colblock,
364  std::vector<casadi_int>& coarse_rowblock,
365  std::vector<casadi_int>& coarse_colblock) const {
366  T()->dmperm(colperm, rowperm, colblock, rowblock,
367  coarse_colblock, coarse_rowblock);
368  return rowblock.size()-1;
369  }
370 
372  const Btf& btf() const;
373 
381  void dmperm(std::vector<casadi_int>& rowperm, std::vector<casadi_int>& colperm,
382  std::vector<casadi_int>& rowblock, std::vector<casadi_int>& colblock,
383  std::vector<casadi_int>& coarse_rowblock,
384  std::vector<casadi_int>& coarse_colblock) const;
385 
393  void maxtrans(std::vector<casadi_int>& imatch,
394  std::vector<casadi_int>& jmatch, Sparsity& trans, casadi_int seed) const;
395 
403  void augment(casadi_int k, std::vector<casadi_int>& jmatch,
404  casadi_int *cheap, std::vector<casadi_int>& w, casadi_int *js,
405  casadi_int *is, casadi_int *ps) const;
406 
416  static std::vector<casadi_int> randperm(casadi_int n, casadi_int seed);
417 
421  static std::vector<casadi_int> invertPermutation(const std::vector<casadi_int>& p);
422 
424  Sparsity permute(const std::vector<casadi_int>& pinv,
425  const std::vector<casadi_int>& q, casadi_int values) const;
426 
434  void permute(const std::vector<casadi_int>& pinv,
435  const std::vector<casadi_int>& q, casadi_int values,
436  std::vector<casadi_int>& colind_C,
437  std::vector<casadi_int>& row_C) const;
438 
446  static casadi_int wclear(casadi_int mark, casadi_int lemax, casadi_int *w, casadi_int n);
447 
455  static casadi_int diag(casadi_int i, casadi_int j, double aij, void *other);
456 
464  Sparsity multiply(const Sparsity& B) const;
465 
473  casadi_int scatter(casadi_int j, std::vector<casadi_int>& w, casadi_int mark,
474  casadi_int* Ci, casadi_int nz) const;
475 
477  std::vector<casadi_int> get_row() const;
478 
480  std::vector<casadi_int> get_colind() const;
481 
483  std::vector<casadi_int> get_col() const;
484 
486  Sparsity _resize(casadi_int nrow, casadi_int ncol) const;
487 
489  Sparsity _reshape(casadi_int nrow, casadi_int ncol) const;
490 
492  casadi_int numel() const;
493 
495  casadi_int nnz_lower(bool strictly=false) const;
496 
498  casadi_int nnz_upper(bool strictly=false) const;
499 
501  casadi_int nnz_diag() const;
502 
506  casadi_int bw_upper() const;
507 
511  casadi_int bw_lower() const;
512 
514  std::pair<casadi_int, casadi_int> size() const;
515 
517  bool is_scalar(bool scalar_and_dense) const;
518 
525  bool is_empty(bool both=false) const;
526 
528  bool is_dense() const;
529 
533  bool is_row() const;
534 
538  bool is_column() const;
539 
543  bool is_vector() const;
544 
546  bool is_diag() const;
547 
549  bool is_square() const;
550 
554  bool is_permutation() const;
555 
559  bool is_selection(bool allow_empty=false) const;
560 
564  bool is_orthonormal(bool allow_empty=false) const;
565 
569  bool is_orthonormal_rows(bool allow_empty=false) const;
570 
574  bool is_orthonormal_columns(bool allow_empty=false) const;
575 
577  bool is_symmetric() const;
578 
580  bool is_tril(bool strictly) const;
581 
583  bool is_triu(bool strictly) const;
584 
586  Sparsity _triu(bool includeDiagonal) const;
587 
589  Sparsity _tril(bool includeDiagonal) const;
590 
592  std::vector<casadi_int> get_lower() const;
593 
595  std::vector<casadi_int> get_upper() const;
596 
598  std::string dim(bool with_nz=false) const;
599 
601  std::string repr_el(casadi_int k) const;
602 
604  Sparsity _mtimes(const Sparsity& y) const;
605 
608  Sparsity combine(const Sparsity& y, bool f0x_is_zero, bool function0_is_zero,
609  std::vector<unsigned char>& mapping) const;
610  Sparsity combine(const Sparsity& y, bool f0x_is_zero, bool function0_is_zero) const;
611 
612  template<bool with_mapping>
613  Sparsity combineGen1(const Sparsity& y, bool f0x_is_zero, bool function0_is_zero,
614  std::vector<unsigned char>& mapping) const;
615 
616  template<bool with_mapping, bool f0x_is_zero, bool function0_is_zero>
617  Sparsity combineGen(const Sparsity& y, std::vector<unsigned char>& mapping) const;
619 
621  bool is_subset(const Sparsity& rhs) const;
622 
624  Sparsity pattern_inverse() const;
625 
627  bool is_equal(const Sparsity& y) const;
628 
630  bool is_equal(casadi_int y_nrow, casadi_int y_ncol, const std::vector<casadi_int>& y_colind,
631  const std::vector<casadi_int>& y_row) const;
632 
634  bool is_equal(casadi_int y_nrow, casadi_int y_ncol,
635  const casadi_int* y_colind, const casadi_int* y_row) const;
636 
638  bool is_stacked(const Sparsity& y, casadi_int n) const;
639 
641  Sparsity _enlargeRows(casadi_int nrow, const std::vector<casadi_int>& rr, bool ind1) const;
642 
644  Sparsity _enlargeColumns(casadi_int ncol, const std::vector<casadi_int>& cc, bool ind1) const;
645 
647  Sparsity makeDense(std::vector<casadi_int>& mapping) const;
648 
650  Sparsity _erase(const std::vector<casadi_int>& rr, const std::vector<casadi_int>& cc,
651  bool ind1, std::vector<casadi_int>& mapping) const;
652 
654  Sparsity _erase(const std::vector<casadi_int>& rr, bool ind1,
655  std::vector<casadi_int>& mapping) const;
656 
658  Sparsity _appendVector(const SparsityInternal& sp) const;
659 
661  Sparsity _appendColumns(const SparsityInternal& sp) const;
662 
669  Sparsity sub(const std::vector<casadi_int>& rr, const std::vector<casadi_int>& cc,
670  std::vector<casadi_int>& mapping, bool ind1) const;
671 
678  Sparsity sub(const std::vector<casadi_int>& rr, const SparsityInternal& sp,
679  std::vector<casadi_int>& mapping, bool ind1) const;
680 
682  casadi_int get_nz(casadi_int rr, casadi_int cc) const;
683 
685  std::vector<casadi_int> get_nz(const std::vector<casadi_int>& rr,
686  const std::vector<casadi_int>& cc) const;
687 
689  void get_nz(std::vector<casadi_int>& indices) const;
690 
692  bool rowsSequential(bool strictly) const;
693 
700  Sparsity _removeDuplicates(std::vector<casadi_int>& mapping) const;
701 
703  void find(std::vector<casadi_int>& loc, bool ind1) const;
704 
706  std::size_t hash() const;
707 
709  std::string class_name() const override {return "SparsityInternal";}
710 
712  void disp(std::ostream& stream, bool more) const override;
713 
720  Sparsity uni_coloring(const Sparsity& AT, casadi_int cutoff) const;
721 
727  Sparsity star_coloring(casadi_int ordering, casadi_int cutoff) const;
728 
734  Sparsity star_coloring2(casadi_int ordering, casadi_int cutoff) const;
735 
737  std::vector<casadi_int> largest_first() const;
738 
740  Sparsity pmult(const std::vector<casadi_int>& p, bool permute_rows=true, bool permute_cols=true,
741  bool invert_permutation=false) const;
742 
746  void spy(std::ostream &stream) const;
747 
749  void spy_matlab(const std::string& mfile) const;
750 
754  void export_code(const std::string& lang, std::ostream &stream,
755  const Dict& options) const;
756 
758  void spsolve(bvec_t* X, bvec_t* B, bool tr) const;
759 };
760 
761 } // namespace casadi
763 
764 #endif // CASADI_SPARSITY_INTERNAL_HPP
The casadi namespace.
Definition: archiver.hpp:32
GenericType::Dict Dict
C++ equivalent of Python's dict or MATLAB's struct.