StrumpackSparseSolver.hpp
1 /*
2  * STRUMPACK -- STRUctured Matrices PACKage, Copyright (c) 2014, The
3  * Regents of the University of California, through Lawrence Berkeley
4  * National Laboratory (subject to receipt of any required approvals
5  * from the U.S. Dept. of Energy). All rights reserved.
6  *
7  * If you have questions about your rights to use or distribute this
8  * software, please contact Berkeley Lab's Technology Transfer
9  * Department at TTD@lbl.gov.
10  *
11  * NOTICE. This software is owned by the U.S. Department of Energy. As
12  * such, the U.S. Government has been granted for itself and others
13  * acting on its behalf a paid-up, nonexclusive, irrevocable,
14  * worldwide license in the Software to reproduce, prepare derivative
15  * works, and perform publicly and display publicly. Beginning five
16  * (5) years after the date permission to assert copyright is obtained
17  * from the U.S. Department of Energy, and subject to any subsequent
18  * five (5) year renewals, the U.S. Government is granted for itself
19  * and others acting on its behalf a paid-up, nonexclusive,
20  * irrevocable, worldwide license in the Software to reproduce,
21  * prepare derivative works, distribute copies to the public, perform
22  * publicly and display publicly, and to permit others to do so.
23  *
24  * Developers: Pieter Ghysels, Francois-Henry Rouet, Xiaoye S. Li.
25  * (Lawrence Berkeley National Lab, Computational Research
26  * Division).
27  */
33 #ifndef STRUMPACK_SPARSE_SOLVER_HPP
34 #define STRUMPACK_SPARSE_SOLVER_HPP
35 
36 #include <new>
37 #include <memory>
38 #include <vector>
39 #include <string>
40 
41 #include "SparseSolverBase.hpp"
42 
46 namespace strumpack {
47 
48  // forward declarations
49  template<typename scalar_t,typename integer_t> class MatrixReordering;
50  template<typename scalar_t,typename integer_t> class EliminationTree;
51  class TaskTimer;
52 
73  template<typename scalar_t,typename integer_t=int>
74  class SparseSolver :
75  public SparseSolverBase<scalar_t,integer_t> {
76 
82 
83  public:
84 
97  SparseSolver(int argc, char* argv[], bool verbose=true, bool root=true);
98 
107  SparseSolver(bool verbose=true, bool root=true);
108 
113 
127 
147  void set_csr_matrix(integer_t N,
148  const integer_t* row_ptr, const integer_t* col_ind,
149  const scalar_t* values, bool symmetric_pattern=false);
150 
180  void update_matrix_values(integer_t N,
181  const integer_t* row_ptr,
182  const integer_t* col_ind,
183  const scalar_t* values,
184  bool symmetric_pattern=false);
185 
205 
206  private:
207  void setup_tree() override;
208  void setup_reordering() override;
209  int compute_reordering(const int* p, int base,
210  int nx, int ny, int nz,
211  int components, int width) override;
212  void separator_reordering() override;
213 
214  SpMat_t* matrix() override { return mat_.get(); }
215  Reord_t* reordering() override { return nd_.get(); }
216  Tree_t* tree() override { return tree_.get(); }
217  const SpMat_t* matrix() const override { return mat_.get(); }
218  const Reord_t* reordering() const override { return nd_.get(); }
219  const Tree_t* tree() const override { return tree_.get(); }
220 
221  void permute_matrix_values();
222 
223  ReturnCode solve_internal
224  (const scalar_t* b, scalar_t* x, bool use_initial_guess=false) override;
225  ReturnCode solve_internal
226  (const DenseM_t& b, DenseM_t& x, bool use_initial_guess=false) override;
227 
228  void delete_factors_internal() override;
229 
230  std::unique_ptr<CSRMatrix<scalar_t,integer_t>> mat_;
231  std::unique_ptr<MatrixReordering<scalar_t,integer_t>> nd_;
232  std::unique_ptr<EliminationTree<scalar_t,integer_t>> tree_;
233 
234  using SPBase_t = SparseSolverBase<scalar_t,integer_t>;
235  using SPBase_t::opts_;
236  using SPBase_t::is_root_;
237  using SPBase_t::matching_;
238  using SPBase_t::equil_;
239  using SPBase_t::factored_;
240  using SPBase_t::reordered_;
241  using SPBase_t::Krylov_its_;
242  };
243 
244  template<typename scalar_t,typename integer_t>
245  using StrumpackSparseSolver = SparseSolver<scalar_t,integer_t>;
246 
247 } //end namespace strumpack
248 
249 #endif // STRUMPACK_SPARSE_SOLVER_HPP
Class for storing a compressed sparse row matrix (single node).
Definition: CSRMatrix.hpp:54
Abstract base class for compressed sparse matrix storage.
Definition: CompressedSparseMatrix.hpp:163
Like DenseMatrix, this class represents a matrix, stored in column major format, to allow direct use ...
Definition: DenseMatrix.hpp:1015
This class represents a matrix, stored in column major format, to allow direct use of BLAS/LAPACK rou...
Definition: DenseMatrix.hpp:138
Definition: StrumpackSparseSolver.hpp:50
Definition: StrumpackSparseSolver.hpp:49
SparseSolver is the main sequential or multithreaded sparse solver class.
Definition: StrumpackSparseSolver.hpp:75
void update_matrix_values(const CSRMatrix< scalar_t, integer_t > &A)
void update_matrix_values(integer_t N, const integer_t *row_ptr, const integer_t *col_ind, const scalar_t *values, bool symmetric_pattern=false)
SparseSolver(bool verbose=true, bool root=true)
SparseSolver(int argc, char *argv[], bool verbose=true, bool root=true)
void set_matrix(const CSRMatrix< scalar_t, integer_t > &A)
void set_csr_matrix(integer_t N, const integer_t *row_ptr, const integer_t *col_ind, const scalar_t *values, bool symmetric_pattern=false)
Definition: StrumpackOptions.hpp:42
ReturnCode
Enumeration for the possible return codes.
Definition: StrumpackParameters.hpp:60