suanPan
MetaMat.hpp
Go to the documentation of this file.
1 /*******************************************************************************
2  * Copyright (C) 2017-2022 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 METAMAT_HPP
30 #define METAMAT_HPP
31 
32 #include "triplet_form.hpp"
33 #include "IterativeSolver.hpp"
34 #include "ILU.hpp"
35 #include "Jacobi.hpp"
36 
37 template<typename T, typename U> concept ArmaContainer = std::is_floating_point_v<U> && (std::is_convertible_v<T, Mat<U>> || std::is_convertible_v<T, SpMat<U>>) ;
38 
39 template<sp_d T> class MetaMat {
40 protected:
41  bool factored = false;
42 
44 
45 public:
47 
48  const uword n_rows;
49  const uword n_cols;
50  const uword n_elem;
51 
52  MetaMat(uword, uword, uword);
53  MetaMat(const MetaMat&) = default;
54  MetaMat(MetaMat&&) noexcept = delete;
55  MetaMat& operator=(const MetaMat&) = delete;
56  MetaMat& operator=(MetaMat&&) noexcept = delete;
57  virtual ~MetaMat() = default;
58 
59  void set_solver_setting(const SolverSetting<T>&);
60  [[nodiscard]] SolverSetting<T>& get_solver_setting();
61 
62  void set_factored(bool);
63 
64  [[nodiscard]] virtual bool is_empty() const = 0;
65  virtual void zeros() = 0;
66 
67  virtual unique_ptr<MetaMat> make_copy() = 0;
68 
69  virtual void unify(uword) = 0;
70  virtual void nullify(uword) = 0;
71 
72  [[nodiscard]] virtual T max() const = 0;
73  [[nodiscard]] virtual Col<T> diag() const = 0;
74 
75  virtual const T& operator()(uword, uword) const = 0;
76  virtual T& at(uword, uword) = 0;
77 
78  [[nodiscard]] virtual const T* memptr() const = 0;
79  virtual T* memptr() = 0;
80 
81  virtual void operator+=(const shared_ptr<MetaMat>&) = 0;
82  virtual void operator-=(const shared_ptr<MetaMat>&) = 0;
83 
84  virtual void operator+=(const triplet_form<T, uword>&) = 0;
85  virtual void operator-=(const triplet_form<T, uword>&) = 0;
86 
87  virtual Mat<T> operator*(const Mat<T>&) const = 0;
88 
89  virtual void operator*=(T) = 0;
90 
91  template<ArmaContainer<T> C> Mat<T> solve(const C& B) {
92  if(IterativeSolver::NONE == this->setting.iterative_solver) return this->direct_solve(B);
93  return this->iterative_solve(B);
94  }
95 
96  template<ArmaContainer<T> C> Mat<T> solve(C&& B) {
97  if(IterativeSolver::NONE == this->setting.iterative_solver) return this->direct_solve(std::forward<C>(B));
98  return this->iterative_solve(std::forward<C>(B));
99  }
100 
101  template<ArmaContainer<T> C> int solve(Mat<T>& X, const C& B) {
102  if(IterativeSolver::NONE == this->setting.iterative_solver) return this->direct_solve(X, B);
103  return this->iterative_solve(X, B);
104  }
105 
106  template<ArmaContainer<T> C> int solve(Mat<T>& X, C&& B) {
107  if(IterativeSolver::NONE == this->setting.iterative_solver) return this->direct_solve(X, std::forward<C>(B));
108  return this->iterative_solve(X, std::forward<C>(B));
109  }
110 
111  template<ArmaContainer<T> C> Mat<T> direct_solve(const C& B) {
112  Mat<T> X;
113  if(0 != this->direct_solve(X, B)) X.reset();
114  return X;
115  }
116 
117  template<ArmaContainer<T> C> Mat<T> direct_solve(C&& B) {
118  Mat<T> X;
119  if(0 != this->direct_solve(X, std::forward<C>(B))) X.reset();
120  return X;
121  }
122 
123  virtual int direct_solve(Mat<T>&, const Mat<T>&) = 0;
124  virtual int direct_solve(Mat<T>&, const SpMat<T>&);
125  virtual int direct_solve(Mat<T>&, Mat<T>&&);
126  virtual int direct_solve(Mat<T>&, SpMat<T>&&);
127 
128  [[nodiscard]] virtual int sign_det() const = 0;
129 
130  void save(const char*);
131 
132  virtual void csc_condense();
133  virtual void csr_condense();
134 
135  Mat<T> iterative_solve(const Mat<T>&);
136  Mat<T> iterative_solve(const SpMat<T>&);
137 
138  virtual int iterative_solve(Mat<T>&, const Mat<T>&);
139  int iterative_solve(Mat<T>&, const SpMat<T>&);
140 
141  [[nodiscard]] Col<T> evaluate(const Col<T>&) const;
142 };
143 
144 template<sp_d T> MetaMat<T>::MetaMat(const uword in_rows, const uword in_cols, const uword in_elem)
145  : triplet_mat(in_rows, in_cols)
146  , n_rows(in_rows)
147  , n_cols(in_cols)
148  , n_elem(in_elem) {}
149 
150 template<sp_d T> void MetaMat<T>::set_solver_setting(const SolverSetting<T>& SS) { setting = SS; }
151 
152 template<sp_d T> SolverSetting<T>& MetaMat<T>::get_solver_setting() { return setting; }
153 
154 template<sp_d T> void MetaMat<T>::set_factored(const bool F) { factored = F; }
155 
156 template<sp_d T> int MetaMat<T>::direct_solve(Mat<T>& X, const SpMat<T>& B) { return this->direct_solve(X, Mat<T>(B)); }
157 
158 template<sp_d T> int MetaMat<T>::direct_solve(Mat<T>& X, Mat<T>&& B) { return this->direct_solve(X, B); }
159 
160 template<sp_d T> int MetaMat<T>::direct_solve(Mat<T>& X, SpMat<T>&& B) { return this->direct_solve(X, B); }
161 
162 template<sp_d T> void MetaMat<T>::save(const char* name) { if(!to_mat(*this).save(name)) suanpan_error("cannot save matrix to file.\n"); }
163 
164 template<sp_d T> void MetaMat<T>::csc_condense() {}
165 
166 template<sp_d T> void MetaMat<T>::csr_condense() {}
167 
168 template<sp_d T> Mat<T> MetaMat<T>::iterative_solve(const Mat<T>& B) {
169  Mat<T> X;
170  if(SUANPAN_SUCCESS != this->iterative_solve(X, B)) X.reset();
171  return X;
172 }
173 
174 template<sp_d T> Mat<T> MetaMat<T>::iterative_solve(const SpMat<T>& B) { return this->iterative_solve(mat(B)); }
175 
176 template<sp_d T> int MetaMat<T>::iterative_solve(Mat<T>& X, const Mat<T>& B) {
177  X.zeros(arma::size(B));
178 
179  unique_ptr<Preconditioner<T>> preconditioner;
180  if(PreconditionerType::JACOBI == this->setting.preconditioner_type) preconditioner = std::make_unique<Jacobi<T>>(this->diag());
181 #ifndef SUANPAN_SUPERLUMT
182  else if(PreconditionerType::ILU == this->setting.preconditioner_type) {
183  if(this->triplet_mat.is_empty()) preconditioner = std::make_unique<ILU<T>>(to_triplet_form<T, int>(this));
184  else preconditioner = std::make_unique<ILU<T>>(this->triplet_mat);
185  }
186 #endif
187  else if(PreconditionerType::NONE == this->setting.preconditioner_type) preconditioner = std::make_unique<UnityPreconditioner<T>>();
188 
189  if(SUANPAN_SUCCESS != preconditioner->init()) return SUANPAN_FAIL;
190 
191  this->setting.preconditioner = preconditioner.get();
192 
193  std::atomic_int code = 0;
194 
195  if(IterativeSolver::GMRES == setting.iterative_solver)
196  suanpan_for(0llu, B.n_cols, [&](const uword I) {
197  Col<T> sub_x(X.colptr(I), X.n_rows, false, true);
198  const Col<T> sub_b(B.colptr(I), B.n_rows);
199  auto col_setting = setting;
200  code += GMRES(this, sub_x, sub_b, col_setting);
201  });
202  else if(IterativeSolver::BICGSTAB == setting.iterative_solver)
203  suanpan_for(0llu, B.n_cols, [&](const uword I) {
204  Col<T> sub_x(X.colptr(I), X.n_rows, false, true);
205  const Col<T> sub_b(B.colptr(I), B.n_rows);
206  auto col_setting = setting;
207  code += BiCGSTAB(this, sub_x, sub_b, col_setting);
208  });
209  else throw invalid_argument("no proper iterative solver assigned but somehow iterative solving is called");
210 
211  return 0 == code ? SUANPAN_SUCCESS : SUANPAN_FAIL;
212 }
213 
214 template<sp_d T> int MetaMat<T>::iterative_solve(Mat<T>& X, const SpMat<T>& B) { return this->iterative_solve(X, mat(B)); }
215 
216 template<sp_d T> Col<T> MetaMat<T>::evaluate(const Col<T>& X) const { return this->operator*(X); }
217 
218 template<sp_d T> Mat<T> to_mat(const MetaMat<T>& in_mat) {
219  Mat<T> out_mat(in_mat.n_rows, in_mat.n_cols);
220  for(uword J = 0; J < in_mat.n_cols; ++J) for(uword I = 0; I < in_mat.n_rows; ++I) out_mat(I, J) = in_mat(I, J);
221  return out_mat;
222 }
223 
224 template<sp_d T> Mat<T> to_mat(const shared_ptr<MetaMat<T>>& in_mat) { return to_mat(*in_mat); }
225 
226 template<sp_d data_t, sp_i index_t> Mat<data_t> to_mat(const triplet_form<data_t, index_t>& in_mat) {
227  Mat<data_t> out_mat(in_mat.n_rows, in_mat.n_cols, fill::zeros);
228  for(index_t I = 0; I < in_mat.n_elem; ++I) out_mat(in_mat.row(I), in_mat.col(I)) += in_mat.val(I);
229  return out_mat;
230 }
231 
232 template<sp_d data_t, sp_i index_t> Mat<data_t> to_mat(const csr_form<data_t, index_t>& in_mat) {
233  Mat<data_t> out_mat(in_mat.n_rows, in_mat.n_cols, fill::zeros);
234 
235  index_t c_idx = 1;
236  for(index_t I = 0; I < in_mat.n_elem; ++I) {
237  if(I >= in_mat.row_mem()[c_idx]) ++c_idx;
238  out_mat(c_idx - 1, in_mat.col_mem()[I]) += in_mat.val_mem()[I];
239  }
240 
241  return out_mat;
242 }
243 
244 template<sp_d data_t, sp_i index_t> Mat<data_t> to_mat(const csc_form<data_t, index_t>& in_mat) {
245  Mat<data_t> out_mat(in_mat.n_rows, in_mat.n_cols, fill::zeros);
246 
247  index_t c_idx = 1;
248  for(index_t I = 0; I < in_mat.n_elem; ++I) {
249  if(I >= in_mat.col_mem()[c_idx]) ++c_idx;
250  out_mat(in_mat.row_mem()[I], c_idx - 1) += in_mat.val_mem()[I];
251  }
252 
253  return out_mat;
254 }
255 
256 template<sp_d data_t, sp_i index_t> triplet_form<data_t, index_t> to_triplet_form(MetaMat<data_t>* in_mat) {
257  if(!in_mat->triplet_mat.is_empty()) return triplet_form<data_t, index_t>(in_mat->triplet_mat);
258 
259  const sp_i auto n_rows = index_t(in_mat->n_rows);
260  const sp_i auto n_cols = index_t(in_mat->n_cols);
261  const sp_i auto n_elem = index_t(in_mat->n_elem);
262 
263  triplet_form<data_t, index_t> out_mat(n_rows, n_cols, n_elem);
264  for(index_t J = 0; J < n_cols; ++J) for(index_t I = 0; I < n_rows; ++I) out_mat.at(I, J) = in_mat->operator()(I, J);
265 
266  return out_mat;
267 }
268 
269 template<sp_d data_t, sp_i index_t> triplet_form<data_t, index_t> to_triplet_form(const shared_ptr<MetaMat<data_t>>& in_mat) { return to_triplet_form<data_t, index_t>(in_mat.get()); }
270 
271 #endif
272 
A ILU class.
Definition: ILU.hpp:40
A MetaMat class that holds matrices.
Definition: MetaMat.hpp:39
triplet_form< T, uword > triplet_mat
Definition: MetaMat.hpp:46
MetaMat(const MetaMat &)=default
virtual void unify(uword)=0
virtual unique_ptr< MetaMat > make_copy()=0
virtual T max() const =0
virtual int sign_det() const =0
const uword n_cols
Definition: MetaMat.hpp:49
int solve(Mat< T > &X, C &&B)
Definition: MetaMat.hpp:106
Mat< T > solve(C &&B)
Definition: MetaMat.hpp:96
virtual T & at(uword, uword)=0
virtual Col< T > diag() const =0
virtual bool is_empty() const =0
const uword n_rows
Definition: MetaMat.hpp:48
virtual const T * memptr() const =0
virtual void nullify(uword)=0
MetaMat(MetaMat &&) noexcept=delete
Mat< T > direct_solve(C &&B)
Definition: MetaMat.hpp:117
bool factored
Definition: MetaMat.hpp:41
virtual int direct_solve(Mat< T > &, const Mat< T > &)=0
Mat< T > solve(const C &B)
Definition: MetaMat.hpp:91
int solve(Mat< T > &X, const C &B)
Definition: MetaMat.hpp:101
const uword n_elem
Definition: MetaMat.hpp:50
SolverSetting< T > setting
Definition: MetaMat.hpp:43
Mat< T > direct_solve(const C &B)
Definition: MetaMat.hpp:111
virtual void zeros()=0
Definition: csc_form.hpp:25
const index_t n_rows
Definition: csc_form.hpp:50
const index_t n_cols
Definition: csc_form.hpp:51
const data_t * val_mem() const
Definition: csc_form.hpp:65
const index_t n_elem
Definition: csc_form.hpp:52
const index_t * col_mem() const
Definition: csc_form.hpp:63
const index_t * row_mem() const
Definition: csc_form.hpp:61
Definition: csr_form.hpp:25
const index_t * col_mem() const
Definition: csr_form.hpp:63
const index_t n_rows
Definition: csr_form.hpp:50
const index_t n_cols
Definition: csr_form.hpp:51
const index_t * row_mem() const
Definition: csr_form.hpp:61
const data_t * val_mem() const
Definition: csr_form.hpp:65
const index_t n_elem
Definition: csr_form.hpp:52
const index_t n_rows
Definition: triplet_form.hpp:128
bool is_empty() const
Definition: triplet_form.hpp:169
data_t & at(index_t, index_t)
Definition: triplet_form.hpp:382
const index_t n_cols
Definition: triplet_form.hpp:129
const index_t n_elem
Definition: triplet_form.hpp:130
index_t col(const index_t I) const
Definition: triplet_form.hpp:161
data_t val(const index_t I) const
Definition: triplet_form.hpp:163
index_t row(const index_t I) const
Definition: triplet_form.hpp:159
void set_solver_setting(const SolverSetting< T > &)
Definition: MetaMat.hpp:150
void set_factored(bool)
Definition: MetaMat.hpp:154
Col< T > evaluate(const Col< T > &) const
Definition: MetaMat.hpp:216
void save(const char *)
Definition: MetaMat.hpp:162
Mat< T > iterative_solve(const Mat< T > &)
Definition: MetaMat.hpp:168
SolverSetting< T > & get_solver_setting()
Definition: MetaMat.hpp:152
virtual void csr_condense()
Definition: MetaMat.hpp:166
virtual void csc_condense()
Definition: MetaMat.hpp:164
MetaMat(uword, uword, uword)
Definition: MetaMat.hpp:144
triplet_form< data_t, index_t > to_triplet_form(MetaMat< data_t > *in_mat)
Definition: MetaMat.hpp:256
concept ArmaContainer
Definition: MetaMat.hpp:37
Mat< T > to_mat(const MetaMat< T > &in_mat)
Definition: MetaMat.hpp:218
unique_ptr< MetaMat< T > > operator*(const T value, const unique_ptr< MetaMat< T >> &M)
Definition: operator_times.hpp:24
void suanpan_error(const char *M,...)
Definition: print.cpp:116
IterativeSolver iterative_solver
Definition: SolverSetting.hpp:46
constexpr auto SUANPAN_SUCCESS
Definition: suanPan.h:161
constexpr auto SUANPAN_FAIL
Definition: suanPan.h:162
concept sp_i
Definition: suanPan.h:228
void suanpan_for(const IT start, const IT end, F &&FN)
Definition: utility.h:24