suanPan
MetaMat.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 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 
79  virtual const T& operator()(uword, uword) const = 0;
84  virtual T& unsafe_at(uword, uword);
89  virtual T& at(uword, uword) = 0;
90 
91  [[nodiscard]] virtual const T* memptr() const = 0;
92  virtual T* memptr() = 0;
93 
94  virtual void operator+=(const shared_ptr<MetaMat>&) = 0;
95  virtual void operator-=(const shared_ptr<MetaMat>&) = 0;
96 
97  virtual void operator+=(const triplet_form<T, uword>&) = 0;
98  virtual void operator-=(const triplet_form<T, uword>&) = 0;
99 
100  virtual Mat<T> operator*(const Mat<T>&) const = 0;
101 
102  virtual void operator*=(T) = 0;
103 
104  template<ArmaContainer<T> C> Mat<T> solve(const C& B) {
105  if(IterativeSolver::NONE == this->setting.iterative_solver) return this->direct_solve(B);
106  return this->iterative_solve(B);
107  }
108 
109  template<ArmaContainer<T> C> Mat<T> solve(C&& B) {
110  if(IterativeSolver::NONE == this->setting.iterative_solver) return this->direct_solve(std::forward<C>(B));
111  return this->iterative_solve(std::forward<C>(B));
112  }
113 
114  template<ArmaContainer<T> C> int solve(Mat<T>& X, const C& B) {
115  if(IterativeSolver::NONE == this->setting.iterative_solver) return this->direct_solve(X, B);
116  return this->iterative_solve(X, B);
117  }
118 
119  template<ArmaContainer<T> C> int solve(Mat<T>& X, C&& B) {
120  if(IterativeSolver::NONE == this->setting.iterative_solver) return this->direct_solve(X, std::forward<C>(B));
121  return this->iterative_solve(X, std::forward<C>(B));
122  }
123 
124  template<ArmaContainer<T> C> Mat<T> direct_solve(const C& B) {
125  Mat<T> X;
126  if(0 != this->direct_solve(X, B)) X.reset();
127  return X;
128  }
129 
130  template<ArmaContainer<T> C> Mat<T> direct_solve(C&& B) {
131  Mat<T> X;
132  if(0 != this->direct_solve(X, std::forward<C>(B))) X.reset();
133  return X;
134  }
135 
136  virtual int direct_solve(Mat<T>&, const Mat<T>&) = 0;
137  virtual int direct_solve(Mat<T>&, const SpMat<T>&);
138  virtual int direct_solve(Mat<T>&, Mat<T>&&);
139  virtual int direct_solve(Mat<T>&, SpMat<T>&&);
140 
141  [[nodiscard]] virtual int sign_det() const = 0;
142 
143  void save(const char*);
144 
145  virtual void csc_condense();
146  virtual void csr_condense();
147 
148  Mat<T> iterative_solve(const Mat<T>&);
149  Mat<T> iterative_solve(const SpMat<T>&);
150 
151  virtual int iterative_solve(Mat<T>&, const Mat<T>&);
152  int iterative_solve(Mat<T>&, const SpMat<T>&);
153 
154  [[nodiscard]] Col<T> evaluate(const Col<T>&) const;
155 };
156 
157 template<sp_d T> MetaMat<T>::MetaMat(const uword in_rows, const uword in_cols, const uword in_elem)
158  : triplet_mat(in_rows, in_cols)
159  , n_rows(in_rows)
160  , n_cols(in_cols)
161  , n_elem(in_elem) {}
162 
163 template<sp_d T> void MetaMat<T>::set_solver_setting(const SolverSetting<T>& SS) { setting = SS; }
164 
165 template<sp_d T> SolverSetting<T>& MetaMat<T>::get_solver_setting() { return setting; }
166 
167 template<sp_d T> void MetaMat<T>::set_factored(const bool F) { factored = F; }
168 
169 template<sp_d T> T& MetaMat<T>::unsafe_at(const uword I, const uword J) { return this->at(I, J); }
170 
171 template<sp_d T> int MetaMat<T>::direct_solve(Mat<T>& X, const SpMat<T>& B) { return this->direct_solve(X, Mat<T>(B)); }
172 
173 template<sp_d T> int MetaMat<T>::direct_solve(Mat<T>& X, Mat<T>&& B) { return this->direct_solve(X, B); }
174 
175 template<sp_d T> int MetaMat<T>::direct_solve(Mat<T>& X, SpMat<T>&& B) { return this->direct_solve(X, B); }
176 
177 template<sp_d T> void MetaMat<T>::save(const char* name) {
178  if(!to_mat(*this).save(name))
179  suanpan_error("Cannot save to file \"{}\".\n", name);
180 }
181 
182 template<sp_d T> void MetaMat<T>::csc_condense() {}
183 
184 template<sp_d T> void MetaMat<T>::csr_condense() {}
185 
186 template<sp_d T> Mat<T> MetaMat<T>::iterative_solve(const Mat<T>& B) {
187  Mat<T> X;
188  if(SUANPAN_SUCCESS != this->iterative_solve(X, B)) X.reset();
189  return X;
190 }
191 
192 template<sp_d T> Mat<T> MetaMat<T>::iterative_solve(const SpMat<T>& B) { return this->iterative_solve(mat(B)); }
193 
194 template<sp_d T> int MetaMat<T>::iterative_solve(Mat<T>& X, const Mat<T>& B) {
195  X.zeros(arma::size(B));
196 
197  unique_ptr<Preconditioner<T>> preconditioner;
198  if(PreconditionerType::JACOBI == this->setting.preconditioner_type) preconditioner = std::make_unique<Jacobi<T>>(this->diag());
199 #ifndef SUANPAN_SUPERLUMT
200  else if(PreconditionerType::ILU == this->setting.preconditioner_type) {
201  if(this->triplet_mat.is_empty()) preconditioner = std::make_unique<ILU<T>>(to_triplet_form<T, int>(this));
202  else preconditioner = std::make_unique<ILU<T>>(this->triplet_mat);
203  }
204 #endif
205  else if(PreconditionerType::NONE == this->setting.preconditioner_type) preconditioner = std::make_unique<UnityPreconditioner<T>>();
206 
207  if(SUANPAN_SUCCESS != preconditioner->init()) return SUANPAN_FAIL;
208 
209  this->setting.preconditioner = preconditioner.get();
210 
211  std::atomic_int code = 0;
212 
213  if(IterativeSolver::GMRES == setting.iterative_solver)
214  suanpan_for(0llu, B.n_cols, [&](const uword I) {
215  Col<T> sub_x(X.colptr(I), X.n_rows, false, true);
216  const Col<T> sub_b(B.colptr(I), B.n_rows);
217  auto col_setting = setting;
218  code += GMRES(this, sub_x, sub_b, col_setting);
219  });
220  else if(IterativeSolver::BICGSTAB == setting.iterative_solver)
221  suanpan_for(0llu, B.n_cols, [&](const uword I) {
222  Col<T> sub_x(X.colptr(I), X.n_rows, false, true);
223  const Col<T> sub_b(B.colptr(I), B.n_rows);
224  auto col_setting = setting;
225  code += BiCGSTAB(this, sub_x, sub_b, col_setting);
226  });
227  else throw invalid_argument("no proper iterative solver assigned but somehow iterative solving is called");
228 
229  return 0 == code ? SUANPAN_SUCCESS : SUANPAN_FAIL;
230 }
231 
232 template<sp_d T> int MetaMat<T>::iterative_solve(Mat<T>& X, const SpMat<T>& B) { return this->iterative_solve(X, mat(B)); }
233 
234 template<sp_d T> Col<T> MetaMat<T>::evaluate(const Col<T>& X) const { return this->operator*(X); }
235 
236 template<sp_d T> Mat<T> to_mat(const MetaMat<T>& in_mat) {
237  Mat<T> out_mat(in_mat.n_rows, in_mat.n_cols);
238  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);
239  return out_mat;
240 }
241 
242 template<sp_d T> Mat<T> to_mat(const shared_ptr<MetaMat<T>>& in_mat) { return to_mat(*in_mat); }
243 
244 template<sp_d data_t, sp_i index_t> Mat<data_t> to_mat(const triplet_form<data_t, index_t>& in_mat) {
245  Mat<data_t> out_mat(in_mat.n_rows, in_mat.n_cols, fill::zeros);
246  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);
247  return out_mat;
248 }
249 
250 template<sp_d data_t, sp_i index_t> Mat<data_t> to_mat(const csr_form<data_t, index_t>& in_mat) {
251  Mat<data_t> out_mat(in_mat.n_rows, in_mat.n_cols, fill::zeros);
252 
253  index_t c_idx = 1;
254  for(index_t I = 0; I < in_mat.n_elem; ++I) {
255  if(I >= in_mat.row_mem()[c_idx]) ++c_idx;
256  out_mat(c_idx - 1, in_mat.col_mem()[I]) += in_mat.val_mem()[I];
257  }
258 
259  return out_mat;
260 }
261 
262 template<sp_d data_t, sp_i index_t> Mat<data_t> to_mat(const csc_form<data_t, index_t>& in_mat) {
263  Mat<data_t> out_mat(in_mat.n_rows, in_mat.n_cols, fill::zeros);
264 
265  index_t c_idx = 1;
266  for(index_t I = 0; I < in_mat.n_elem; ++I) {
267  if(I >= in_mat.col_mem()[c_idx]) ++c_idx;
268  out_mat(in_mat.row_mem()[I], c_idx - 1) += in_mat.val_mem()[I];
269  }
270 
271  return out_mat;
272 }
273 
274 template<sp_d data_t, sp_i index_t> triplet_form<data_t, index_t> to_triplet_form(MetaMat<data_t>* in_mat) {
275  if(!in_mat->triplet_mat.is_empty()) return triplet_form<data_t, index_t>(in_mat->triplet_mat);
276 
277  const sp_i auto n_rows = index_t(in_mat->n_rows);
278  const sp_i auto n_cols = index_t(in_mat->n_cols);
279  const sp_i auto n_elem = index_t(in_mat->n_elem);
280 
281  triplet_form<data_t, index_t> out_mat(n_rows, n_cols, n_elem);
282  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);
283 
284  return out_mat;
285 }
286 
287 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()); }
288 
289 #endif
290 
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:119
Mat< T > solve(C &&B)
Definition: MetaMat.hpp:109
virtual T & at(uword, uword)=0
Access element with bound check.
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:130
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:104
int solve(Mat< T > &X, const C &B)
Definition: MetaMat.hpp:114
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:124
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:163
void set_factored(bool)
Definition: MetaMat.hpp:167
Col< T > evaluate(const Col< T > &) const
Definition: MetaMat.hpp:234
void save(const char *)
Definition: MetaMat.hpp:177
Mat< T > iterative_solve(const Mat< T > &)
Definition: MetaMat.hpp:186
SolverSetting< T > & get_solver_setting()
Definition: MetaMat.hpp:165
virtual void csr_condense()
Definition: MetaMat.hpp:184
virtual void csc_condense()
Definition: MetaMat.hpp:182
virtual T & unsafe_at(uword, uword)
Access element without bound check.
Definition: MetaMat.hpp:169
MetaMat(uword, uword, uword)
Definition: MetaMat.hpp:157
triplet_form< data_t, index_t > to_triplet_form(MetaMat< data_t > *in_mat)
Definition: MetaMat.hpp:274
concept ArmaContainer
Definition: MetaMat.hpp:37
Mat< T > to_mat(const MetaMat< T > &in_mat)
Definition: MetaMat.hpp:236
unique_ptr< MetaMat< T > > operator*(const T value, const unique_ptr< MetaMat< T >> &M)
Definition: operator_times.hpp:24
IterativeSolver iterative_solver
Definition: SolverSetting.hpp:46
constexpr auto SUANPAN_SUCCESS
Definition: suanPan.h:162
constexpr auto SUANPAN_FAIL
Definition: suanPan.h:163
#define suanpan_error(...)
Definition: suanPan.h:297
concept sp_i
Definition: suanPan.h:319
void suanpan_for(const IT start, const IT end, F &&FN)
Definition: utility.h:27