StrumpackSparseSolverMPIDist.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  */
32 #ifndef STRUMPACK_SPARSE_SOLVER_MPI_DIST_HPP
33 #define STRUMPACK_SPARSE_SOLVER_MPI_DIST_HPP
34 
35 #include "SparseSolverBase.hpp"
36 #include "dense/ScaLAPACKWrapper.hpp"
37 #include "dense/DistributedVector.hpp"
38 #include "sparse/CSRMatrixMPI.hpp"
39 
40 namespace strumpack {
41 
42  // forward declarations
43  template<typename scalar_t,typename integer_t> class EliminationTreeMPIDist;
44  template<typename scalar_t,typename integer_t> class MatrixReorderingMPI;
45 
70  template<typename scalar_t,typename integer_t>
72  public SparseSolverBase<scalar_t,integer_t> {
73 
79 
80  public:
95  SparseSolverMPIDist(MPI_Comm comm, int argc, char* argv[],
96  bool verbose=true);
97 
107  SparseSolverMPIDist(MPI_Comm comm, bool verbose=true);
108 
113 
127 
148  (integer_t N, const integer_t* row_ptr, const integer_t* col_ind,
149  const scalar_t* values, bool symmetric_pattern=false);
150 
164 
186  (integer_t local_rows, const integer_t* row_ptr,
187  const integer_t* col_ind, const scalar_t* values,
188  const integer_t* dist, bool symmetric_pattern=false);
189 
196  (integer_t local_rows,
197  const integer_t* d_ptr, const integer_t* d_ind, const scalar_t* d_val,
198  const integer_t* o_ptr, const integer_t* o_ind, const scalar_t* o_val,
199  const integer_t* garray);
200 
201 
236  (integer_t local_rows, const integer_t* row_ptr,
237  const integer_t* col_ind, const scalar_t* values,
238  const integer_t* dist, bool symmetric_pattern=false);
239 
260 
278  (integer_t local_rows, const integer_t* d_ptr, const integer_t* d_ind,
279  const scalar_t* d_val, const integer_t* o_ptr, const integer_t* o_ind,
280  const scalar_t* o_val, const integer_t* garray);
281 
286  MPI_Comm comm() const;
287 
292  const MPIComm& Comm() const { return comm_; }
293 
294  private:
297  MPIComm comm_;
298 
299  SpMat_t* matrix() override { return mat_mpi_.get(); }
300  std::unique_ptr<SpMat_t> matrix_nonzero_diag() override {
301  using real_t = typename RealType<scalar_t>::value_type;
302  return mat_mpi_->add_missing_diagonal
303  (std::sqrt(blas::lamch<real_t>('E')) * mat_mpi_->norm1());
304  }
305  Reord_t* reordering() override;
306  Tree_t* tree() override { return tree_mpi_dist_.get(); }
307  const SpMat_t* matrix() const override { return mat_mpi_.get(); }
308  const Reord_t* reordering() const override;
309  const Tree_t* tree() const override { return tree_mpi_dist_.get(); }
310 
311  void setup_tree() override;
312  void setup_reordering() override;
313  int compute_reordering(const int* p, int base, int nx, int ny, int nz,
314  int components, int width) override;
315  void separator_reordering() override;
316 
317  void perf_counters_stop(const std::string& s) override;
318  void synchronize() override { comm_.barrier(); }
319  void reduce_flop_counters() const override;
320 
321  double max_peak_memory() const override {
322  return comm_.reduce(double(params::peak_memory), MPI_MAX);
323  }
324  double min_peak_memory() const override {
325  return comm_.reduce(double(params::peak_memory), MPI_MIN);
326  }
327 
328  void redistribute_values();
329 
330  void delete_factors_internal() override;
331 
332  ReturnCode
333  solve_internal(const scalar_t* b, scalar_t* x,
334  bool use_initial_guess=false) override;
335  ReturnCode
336  solve_internal(const DenseM_t& b, DenseM_t& x,
337  bool use_initial_guess=false) override;
338  ReturnCode
339  solve_internal(int nrhs, const scalar_t* b, int ldb,
340  scalar_t* x, int ldx,
341  bool use_initial_guess=false) override;
342 
343  std::unique_ptr<CSRMatrixMPI<scalar_t,integer_t>> mat_mpi_;
344  std::unique_ptr<MatrixReorderingMPI<scalar_t,integer_t>> nd_mpi_;
345  std::unique_ptr<EliminationTreeMPIDist<scalar_t,integer_t>> tree_mpi_dist_;
346  };
347 
348  template<typename scalar_t,typename integer_t>
349  using StrumpackSparseSolverMPIDist =
350  SparseSolverMPIDist<scalar_t,integer_t>;
351 
352 } // end namespace strumpack
353 
354 #endif // STRUMPACK_SPARSE_SOLVER_MPI_DIST_HPP
Contains the CSRMatrixMPI class definition, a class representing a block-row distributed compressed s...
Contains the definition of the base (abstract/pure virtual) sparse solver class.
Block-row distributed compressed sparse row storage.
Definition: CSRMatrixMPI.hpp:72
Class for storing a compressed sparse row matrix (single node).
Definition: CSRMatrix.hpp:55
Abstract base class for compressed sparse matrix storage.
Definition: CompressedSparseMatrix.hpp:149
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: StrumpackSparseSolverMPIDist.hpp:43
Definition: SparseSolverBase.hpp:53
Wrapper class around an MPI_Comm object.
Definition: MPIWrapper.hpp:194
void barrier() const
Definition: MPIWrapper.hpp:299
T reduce(T t, MPI_Op op) const
Definition: MPIWrapper.hpp:534
Definition: StrumpackSparseSolverMPIDist.hpp:44
Definition: SparseSolverBase.hpp:52
SparseSolverBase is the virtual base for both the sequential/multithreaded and distributed sparse sol...
Definition: SparseSolverBase.hpp:78
This is the fully distributed solver.
Definition: StrumpackSparseSolverMPIDist.hpp:72
SparseSolverMPIDist(MPI_Comm comm, int argc, char *argv[], bool verbose=true)
void set_distributed_csr_matrix(integer_t local_rows, const integer_t *row_ptr, const integer_t *col_ind, const scalar_t *values, const integer_t *dist, bool symmetric_pattern=false)
void set_MPIAIJ_matrix(integer_t local_rows, const integer_t *d_ptr, const integer_t *d_ind, const scalar_t *d_val, const integer_t *o_ptr, const integer_t *o_ind, const scalar_t *o_val, const integer_t *garray)
void broadcast_matrix(const CSRMatrix< scalar_t, integer_t > &A)
void set_matrix(const CSRMatrixMPI< scalar_t, integer_t > &A)
void broadcast_csr_matrix(integer_t N, const integer_t *row_ptr, const integer_t *col_ind, const scalar_t *values, bool symmetric_pattern=false)
void update_matrix_values(const CSRMatrixMPI< scalar_t, integer_t > &A)
void update_MPIAIJ_matrix_values(integer_t local_rows, const integer_t *d_ptr, const integer_t *d_ind, const scalar_t *d_val, const integer_t *o_ptr, const integer_t *o_ind, const scalar_t *o_val, const integer_t *garray)
void update_matrix_values(integer_t local_rows, const integer_t *row_ptr, const integer_t *col_ind, const scalar_t *values, const integer_t *dist, bool symmetric_pattern=false)
const MPIComm & Comm() const
Definition: StrumpackSparseSolverMPIDist.hpp:292
SparseSolverMPIDist(MPI_Comm comm, bool verbose=true)
Definition: StrumpackOptions.hpp:43
ReturnCode
Enumeration for the possible return codes.
Definition: StrumpackParameters.hpp:50