suanPan
SparseMat.hpp
Go to the documentation of this file.
1 /*******************************************************************************
2  * Copyright (C) 2017-2023 Theodore Chang
3  *
4  * This program is free software: you can redistribute it and/or modify
5  * it under the terms of the GNU General Public License as published by
6  * the Free Software Foundation, either version 3 of the License, or
7  * (at your option) any later version.
8  *
9  * This program is distributed in the hope that it will be useful,
10  * but WITHOUT ANY WARRANTY; without even the implied warranty of
11  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12  * GNU General Public License for more details.
13  *
14  * You should have received a copy of the GNU General Public License
15  * along with this program. If not, see <http://www.gnu.org/licenses/>.
16  ******************************************************************************/
29 #ifndef SPARSEMAT_HPP
30 #define SPARSEMAT_HPP
31 
32 #include "MetaMat.hpp"
33 
34 template<sp_d T> class SparseMat : public MetaMat<T> {
35 protected:
37 
38  int direct_solve(Mat<T>& X, Mat<T>&& B) override { return this->direct_solve(X, B); }
39 
40 public:
41  SparseMat(const uword in_row, const uword in_col, const uword in_elem = 0)
42  : MetaMat<T>(in_row, in_col, 0) { this->triplet_mat.init(in_elem); }
43 
44  [[nodiscard]] bool is_empty() const override { return this->triplet_mat.is_empty(); }
45 
46  void zeros() override {
47  this->factored = false;
48  this->triplet_mat.zeros();
49  }
50 
51  void nullify(const uword idx) override {
52  this->factored = false;
53  suanpan_for(static_cast<uword>(0), this->triplet_mat.n_elem, [&](const uword I) { if(this->triplet_mat.row(I) == idx || this->triplet_mat.col(I) == idx) this->triplet_mat.val_mem()[I] = T(0); });
54  }
55 
56  [[nodiscard]] T max() const override { return this->triplet_mat.max(); }
57 
58  [[nodiscard]] Col<T> diag() const override { return this->triplet_mat.diag(); }
59 
60  T operator()(const uword in_row, const uword in_col) const override { return this->triplet_mat(in_row, in_col); }
61 
62  T& at(const uword in_row, const uword in_col) override {
63  this->factored = false;
64  return this->triplet_mat.at(in_row, in_col);
65  }
66 
67  [[nodiscard]] const T* memptr() const override { throw invalid_argument("not supported"); }
68 
69  T* memptr() override { throw invalid_argument("not supported"); }
70 
71  void operator+=(const shared_ptr<MetaMat<T>>& in_mat) override {
72  if(nullptr == in_mat) return;
73  if(!in_mat->triplet_mat.is_empty()) return this->operator+=(in_mat->triplet_mat);
74  this->factored = false;
75  for(uword I = 0llu; I < in_mat->n_rows; ++I) for(uword J = 0llu; J < in_mat->n_cols; ++J) if(const auto t_val = in_mat->operator()(I, J); !suanpan::approx_equal(T(0), t_val)) at(I, J) = t_val;
76  }
77 
78  void operator-=(const shared_ptr<MetaMat<T>>& in_mat) override {
79  if(nullptr == in_mat) return;
80  if(!in_mat->triplet_mat.is_empty()) return this->operator-=(in_mat->triplet_mat);
81  this->factored = false;
82  for(uword I = 0llu; I < in_mat->n_rows; ++I) for(uword J = 0llu; J < in_mat->n_cols; ++J) if(const auto t_val = in_mat->operator()(I, J); !suanpan::approx_equal(T(0), t_val)) at(I, J) = -t_val;
83  }
84 
85  void operator+=(const triplet_form<T, uword>& in_mat) override {
86  this->factored = false;
87  this->triplet_mat += in_mat;
88  }
89 
90  void operator-=(const triplet_form<T, uword>& in_mat) override {
91  this->factored = false;
92  this->triplet_mat -= in_mat;
93  }
94 
95  Mat<T> operator*(const Mat<T>& in_mat) const override { return this->triplet_mat * in_mat; }
96 
97  void operator*=(const T scalar) override {
98  this->factored = false;
99  this->triplet_mat *= scalar;
100  }
101 
102  [[nodiscard]] int sign_det() const override { throw invalid_argument("not supported"); }
103 
104  void csc_condense() override { this->triplet_mat.csc_condense(); }
105 
106  void csr_condense() override { this->triplet_mat.csr_condense(); }
107 };
108 
109 #endif
110 
A MetaMat class that holds matrices.
Definition: MetaMat.hpp:39
triplet_form< T, uword > triplet_mat
Definition: MetaMat.hpp:83
bool factored
Definition: MetaMat.hpp:41
A SparseMat class that holds matrices.
Definition: SparseMat.hpp:34
void zeros() override
Definition: SparseMat.hpp:46
int sign_det() const override
Definition: SparseMat.hpp:102
Mat< T > operator*(const Mat< T > &in_mat) const override
Definition: SparseMat.hpp:95
T operator()(const uword in_row, const uword in_col) const override
Access element (read-only), returns zero if out-of-bound.
Definition: SparseMat.hpp:60
void operator+=(const shared_ptr< MetaMat< T >> &in_mat) override
Definition: SparseMat.hpp:71
Col< T > diag() const override
Definition: SparseMat.hpp:58
void operator-=(const shared_ptr< MetaMat< T >> &in_mat) override
Definition: SparseMat.hpp:78
SparseMat(const uword in_row, const uword in_col, const uword in_elem=0)
Definition: SparseMat.hpp:41
T & at(const uword in_row, const uword in_col) override
Access element with bound check.
Definition: SparseMat.hpp:62
T max() const override
Definition: SparseMat.hpp:56
const T * memptr() const override
Definition: SparseMat.hpp:67
void csc_condense() override
Definition: SparseMat.hpp:104
void nullify(const uword idx) override
Definition: SparseMat.hpp:51
void operator+=(const triplet_form< T, uword > &in_mat) override
Definition: SparseMat.hpp:85
void csr_condense() override
Definition: SparseMat.hpp:106
int direct_solve(Mat< T > &X, Mat< T > &&B) override
Definition: SparseMat.hpp:38
void operator-=(const triplet_form< T, uword > &in_mat) override
Definition: SparseMat.hpp:90
bool is_empty() const override
Definition: SparseMat.hpp:44
void operator*=(const T scalar) override
Definition: SparseMat.hpp:97
T * memptr() override
Definition: SparseMat.hpp:69
void zeros()
Definition: triplet_form.hpp:176
data_t max() const
Definition: triplet_form.hpp:171
void csc_condense()
Definition: triplet_form.hpp:209
bool is_empty() const
Definition: triplet_form.hpp:169
data_t & at(index_t, index_t)
Definition: triplet_form.hpp:384
void csr_condense()
Definition: triplet_form.hpp:204
const index_t n_elem
Definition: triplet_form.hpp:130
void init(const index_t in_elem)
Definition: triplet_form.hpp:181
Col< data_t > diag() const
Definition: triplet_form.hpp:571
std::enable_if_t<!std::numeric_limits< T >::is_integer, bool > approx_equal(T x, T y, int ulp=2)
Definition: utility.h:58
void suanpan_for(const IT start, const IT end, F &&FN)
Definition: utility.h:27