StrumpackSparseSolverBase.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_BASE_HPP
34 #define STRUMPACK_SPARSE_SOLVER_BASE_HPP
35 
36 #include <new>
37 #include <memory>
38 #include <vector>
39 #include <string>
40 
41 #include "StrumpackConfig.hpp"
42 #include "StrumpackOptions.hpp"
43 #include "sparse/CSRMatrix.hpp"
44 #include "dense/DenseMatrix.hpp"
45 
49 namespace strumpack {
50 
51  // forward declarations
52  template<typename scalar_t,typename integer_t> class MatrixReordering;
53  template<typename scalar_t,typename integer_t> class EliminationTree;
54  class TaskTimer;
55 
76  template<typename scalar_t,typename integer_t=int>
78 
84 
85  public:
86 
99  StrumpackSparseSolverBase(int argc, char* argv[],
100  bool verbose=true, bool root=true);
101 
110  StrumpackSparseSolverBase(bool verbose=true, bool root=true);
111 
115  virtual ~StrumpackSparseSolverBase();
116 
146  ReturnCode reorder(int nx=1, int ny=1, int nz=1,
147  int components=1, int width=1);
148 
158  ReturnCode reorder(const int* p, int base=0);
159 
173  ReturnCode factor();
174 
175  void move_to_gpu();
176  void remove_from_gpu();
177 
202  ReturnCode solve(const scalar_t* b, scalar_t* x,
203  bool use_initial_guess=false);
204 
231  ReturnCode solve(const DenseM_t& b, DenseM_t& x,
232  bool use_initial_guess=false);
233 
238 
242  const SPOptions<scalar_t>& options() const;
243 
250  void set_from_options();
251 
260  void set_from_options(int argc, char* argv[]);
261 
270  int maximum_rank() const;
271 
279  std::size_t factor_nonzeros() const;
280 
291  std::size_t factor_memory() const;
292 
297  int Krylov_iterations() const;
298 
307  void draw(const std::string& name) const;
308 
309  protected:
310  virtual void setup_tree() = 0;
311  virtual void setup_reordering() = 0;
312  virtual
313  int compute_reordering(const int* p, int base,
314  int nx, int ny, int nz,
315  int components, int width) = 0;
316  virtual void separator_reordering() = 0;
317 
318  virtual SpMat_t* matrix() = 0;
319  virtual Reord_t* reordering() = 0;
320  virtual Tree_t* tree() = 0;
321  virtual const SpMat_t* matrix() const = 0;
322  virtual const Reord_t* reordering() const = 0;
323  virtual const Tree_t* tree() const = 0;
324 
325  virtual void perf_counters_start();
326  virtual void perf_counters_stop(const std::string& s);
327 
328  virtual void synchronize() {}
329  virtual void communicate_ordering() {}
330 
331  void papi_initialize();
332  long long dense_factor_nonzeros() const;
333  void print_solve_stats(TaskTimer& t) const;
334 
335  virtual void reduce_flop_counters() const {}
336  void print_flop_breakdown_HSS() const;
337  void print_flop_breakdown_HODLR() const;
338  void flop_breakdown_reset() const;
339 
340  void print_wrong_sparsity_error();
341 
342  SPOptions<scalar_t> opts_;
343  bool is_root_;
344 
347 
348  std::new_handler old_handler_;
349  std::ostream* rank_out_ = nullptr;
350  bool factored_ = false;
351  bool reordered_ = false;
352  int Krylov_its_ = 0;
353 
354 #if defined(STRUMPACK_USE_PAPI)
355  float rtime_ = 0., ptime_ = 0.;
356  long_long _flpops = 0;
357 #endif
358 #if defined(STRUMPACK_COUNT_FLOPS)
359  long long int f0_ = 0, ftot_ = 0, fmin_ = 0, fmax_ = 0;
360  long long int b0_ = 0, btot_ = 0, bmin_ = 0, bmax_ = 0;
361  long long int m0_ = 0, mtot_ = 0, mmin_ = 0, mmax_ = 0;
362  long long int ptot_ = 0, pmin_ = 0, pmax_ = 0;
363  long long int dm0_ = 0, dmtot_ = 0, dmmin_ = 0, dmmax_ = 0;
364  long long int dptot_ = 0, dpmin_ = 0, dpmax_ = 0;
365 #endif
366 
367  private:
368  ReturnCode reorder_internal(const int* p, int base,
369  int nx, int ny, int nz,
370  int components, int width);
371 
372  virtual
373  ReturnCode solve_internal(const scalar_t* b, scalar_t* x,
374  bool use_initial_guess=false) = 0;
375  virtual
376  ReturnCode solve_internal(const DenseM_t& b, DenseM_t& x,
377  bool use_initial_guess=false) = 0;
378  };
379 
380 } //end namespace strumpack
381 
382 #endif // STRUMPACK_SPARSE_SOLVER_HPP
strumpack::StrumpackSparseSolverBase::draw
void draw(const std::string &name) const
strumpack::Equilibration< scalar_t >
strumpack::StrumpackSparseSolverBase::factor_memory
std::size_t factor_memory() const
DenseMatrix.hpp
Contains the DenseMatrix and DenseMatrixWrapper classes, simple wrappers around BLAS/LAPACK style den...
strumpack
Definition: StrumpackOptions.hpp:42
strumpack::StrumpackSparseSolverBase::solve
ReturnCode solve(const scalar_t *b, scalar_t *x, bool use_initial_guess=false)
strumpack::DenseMatrixWrapper
Like DenseMatrix, this class represents a matrix, stored in column major format, to allow direct use ...
Definition: DenseMatrix.hpp:990
strumpack::StrumpackSparseSolverBase::reorder
ReturnCode reorder(int nx=1, int ny=1, int nz=1, int components=1, int width=1)
strumpack::SPOptions
Options for the sparse solver.
Definition: StrumpackOptions.hpp:192
strumpack::StrumpackSparseSolverBase::StrumpackSparseSolverBase
StrumpackSparseSolverBase(int argc, char *argv[], bool verbose=true, bool root=true)
strumpack::DenseMatrix< scalar_t >
StrumpackOptions.hpp
Holds options for the sparse solver.
strumpack::StrumpackSparseSolverBase::maximum_rank
int maximum_rank() const
strumpack::StrumpackSparseSolverBase::factor
ReturnCode factor()
strumpack::ReturnCode
ReturnCode
Enumeration for the possible return codes.
Definition: StrumpackParameters.hpp:60
strumpack::StrumpackSparseSolverBase::~StrumpackSparseSolverBase
virtual ~StrumpackSparseSolverBase()
strumpack::StrumpackSparseSolverBase::options
SPOptions< scalar_t > & options()
strumpack::StrumpackSparseSolverBase::factor_nonzeros
std::size_t factor_nonzeros() const
strumpack::StrumpackSparseSolverBase::set_from_options
void set_from_options()
strumpack::CompressedSparseMatrix
Abstract base class for compressed sparse matrix storage.
Definition: CompressedSparseMatrix.hpp:163
strumpack::StrumpackSparseSolverBase::Krylov_iterations
int Krylov_iterations() const
strumpack::EliminationTree
Definition: StrumpackSparseSolverBase.hpp:53
strumpack::MatrixReordering
Definition: StrumpackSparseSolverBase.hpp:52
CSRMatrix.hpp
Contains the compressed sparse row matrix storage class.
strumpack::MatchingData
Definition: CompressedSparseMatrix.hpp:69
strumpack::StrumpackSparseSolverBase
StrumpackSparseSolver is the main sequential or multithreaded sparse solver class.
Definition: StrumpackSparseSolverBase.hpp:77