StrumpackSparseSolver.hpp
Go to the documentation of this file.
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 
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>
75  public StrumpackSparseSolverBase<scalar_t,integer_t> {
76 
82 
83  public:
84 
98  (int argc, char* argv[], bool verbose=true, bool root=true);
99 
108  StrumpackSparseSolver(bool verbose=true, bool root=true);
109 
114 
128 
148  void set_csr_matrix(integer_t N,
149  const integer_t* row_ptr, const integer_t* col_ind,
150  const scalar_t* values, bool symmetric_pattern=false);
151 
181  void update_matrix_values(integer_t N,
182  const integer_t* row_ptr,
183  const integer_t* col_ind,
184  const scalar_t* values,
185  bool symmetric_pattern=false);
186 
206 
207  private:
208  void setup_tree() override;
209  void setup_reordering() override;
210  int compute_reordering(const int* p, int base,
211  int nx, int ny, int nz,
212  int components, int width) override;
213  void separator_reordering() override;
214 
215  SpMat_t* matrix() override { return mat_.get(); }
216  Reord_t* reordering() override { return nd_.get(); }
217  Tree_t* tree() { return tree_.get(); }
218  const SpMat_t* matrix() const override { return mat_.get(); }
219  const Reord_t* reordering() const override { return nd_.get(); }
220  const Tree_t* tree() const override { return tree_.get(); }
221 
222  void permute_matrix_values();
223 
224  ReturnCode solve_internal
225  (const scalar_t* b, scalar_t* x, bool use_initial_guess=false) override;
226  ReturnCode solve_internal
227  (const DenseM_t& b, DenseM_t& x, bool use_initial_guess=false) override;
228 
229  std::unique_ptr<CSRMatrix<scalar_t,integer_t>> mat_;
230  std::unique_ptr<MatrixReordering<scalar_t,integer_t>> nd_;
231  std::unique_ptr<EliminationTree<scalar_t,integer_t>> tree_;
232 
234  using SPBase_t::opts_;
235  using SPBase_t::is_root_;
236  using SPBase_t::matching_;
237  using SPBase_t::equil_;
238  using SPBase_t::factored_;
239  using SPBase_t::reordered_;
240  using SPBase_t::Krylov_its_;
241  };
242 
243 } //end namespace strumpack
244 
245 #endif // STRUMPACK_SPARSE_SOLVER_HPP
strumpack::StrumpackSparseSolver::~StrumpackSparseSolver
~StrumpackSparseSolver()
strumpack::StrumpackSparseSolver
StrumpackSparseSolver is the main sequential or multithreaded sparse solver class.
Definition: StrumpackSparseSolver.hpp:74
strumpack::Trans::N
@ N
strumpack::StrumpackSparseSolver::StrumpackSparseSolver
StrumpackSparseSolver(int argc, char *argv[], bool verbose=true, bool root=true)
strumpack::CSRMatrix
Class for storing a compressed sparse row matrix (single node).
Definition: CSRMatrix.hpp:53
strumpack::StrumpackSparseSolver::update_matrix_values
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)
strumpack
Definition: StrumpackOptions.hpp:42
strumpack::DenseMatrixWrapper
Like DenseMatrix, this class represents a matrix, stored in column major format, to allow direct use ...
Definition: DenseMatrix.hpp:946
StrumpackSparseSolverBase.hpp
Contains the definition of the base (abstract/pure virtual) sparse solver class.
strumpack::DenseMatrix< scalar_t >
strumpack::ReturnCode
ReturnCode
Enumeration for the possible return codes.
Definition: StrumpackParameters.hpp:60
strumpack::StrumpackSparseSolver::set_matrix
void set_matrix(const CSRMatrix< scalar_t, integer_t > &A)
strumpack::StrumpackSparseSolver::set_csr_matrix
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)
strumpack::CompressedSparseMatrix
Abstract base class for compressed sparse matrix storage.
Definition: CompressedSparseMatrix.hpp:163
strumpack::EliminationTree
Definition: StrumpackSparseSolverBase.hpp:53
strumpack::MatrixReordering
Definition: StrumpackSparseSolverBase.hpp:52
strumpack::StrumpackSparseSolverBase
StrumpackSparseSolver is the main sequential or multithreaded sparse solver class.
Definition: StrumpackSparseSolverBase.hpp:77