HODLRMatrix.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  *
28  */
33 #ifndef STRUMPACK_HODLR_MATRIX_HPP
34 #define STRUMPACK_HODLR_MATRIX_HPP
35 
36 #include <cassert>
37 #include <functional>
38 
39 #include "HSS/HSSPartitionTree.hpp"
40 #include "kernel/Kernel.hpp"
41 #include "dense/DistributedMatrix.hpp"
42 #include "HODLROptions.hpp"
43 #include "sparse/CSRGraph.hpp"
44 
45 namespace strumpack {
46 
52  namespace HODLR {
53 
54  struct ExtractionMeta {
55  std::unique_ptr<int[]> iwork;
56  int Ninter, Nallrows, Nallcols, Nalldat_loc,
57  *allrows, *allcols, *rowids, *colids, *pgids, Npmap, *pmaps;
58  };
59 
78  template<typename scalar_t> class HODLRMatrix {
84  using Vec_t = std::vector<std::size_t>;
85  using VecVec_t = std::vector<std::vector<std::size_t>>;
86  using F2Cptr = void*;
87 
88  public:
89  using real_t = typename RealType<scalar_t>::value_type;
90  using mult_t = typename std::function
91  <void(Trans, const DenseM_t&, DenseM_t&)>;
92  using elem_t = typename std::function
93  <scalar_t(std::size_t i, std::size_t j)>;
94  using delem_blocks_t = typename std::function
95  <void(VecVec_t& I, VecVec_t& J, std::vector<DistMW_t>& B,
96  ExtractionMeta&)>;
97  using elem_blocks_t = typename std::function
98  <void(VecVec_t& I, VecVec_t& J, std::vector<DenseMW_t>& B,
99  ExtractionMeta&)>;
100 
105 
118  const opts_t& opts);
119 
137  HODLRMatrix(const MPIComm& c, const HSS::HSSPartitionTree& tree,
138  const std::function<scalar_t(int i, int j)>& Aelem,
139  const opts_t& opts);
140 
141 
162  HODLRMatrix(const MPIComm& c, const HSS::HSSPartitionTree& tree,
163  const elem_blocks_t& Aelem, const opts_t& opts);
164 
165 
183  HODLRMatrix(const MPIComm& c, const HSS::HSSPartitionTree& tree,
184  const std::function<
185  void(Trans op,const DenseM_t& R,DenseM_t& S)>& Amult,
186  const opts_t& opts);
187 
200  HODLRMatrix(const MPIComm& c, const HSS::HSSPartitionTree& tree,
201  const opts_t& opts);
202 
216  template<typename integer_t>
217  HODLRMatrix(const MPIComm& c, const HSS::HSSPartitionTree& tree,
218  const CSRGraph<integer_t>& graph, const opts_t& opts);
219 
223  HODLRMatrix(const HODLRMatrix<scalar_t>& h) = delete;
224 
229  HODLRMatrix(HODLRMatrix<scalar_t>&& h) { *this = std::move(h); }
230 
234  virtual ~HODLRMatrix();
235 
240 
246 
251  std::size_t rows() const { return rows_; }
256  std::size_t cols() const { return cols_; }
261  std::size_t lrows() const { return lrows_; }
266  std::size_t begin_row() const { return dist_[c_->rank()]; }
271  std::size_t end_row() const { return dist_[c_->rank()+1]; }
275  const MPIComm& Comm() const { return *c_; }
276 
281  std::size_t memory() const {
282  return get_stat("Mem_Fill") * 1024 * 1024;
283  }
284 
289  std::size_t total_memory() const {
290  return c_->all_reduce(memory(), MPI_SUM);
291  }
292 
299  std::size_t factor_memory() const {
300  return get_stat("Mem_Factor") * 1.e6;
301  }
302 
310  std::size_t total_factor_memory() const {
311  return c_->all_reduce(factor_memory(), MPI_SUM);
312  }
313 
317  std::size_t rank() const { return get_stat("Rank_max"); }
318 
324  std::size_t max_rank() const {
325  return c_->all_reduce(rank(), MPI_MAX);
326  }
327 
328 
337  double get_stat(const std::string& name) const;
338 
339  void print_stats();
340 
352  void compress
353  (const std::function<void(Trans op,const DenseM_t& R,DenseM_t& S)>& Amult);
354 
367  void compress
368  (const std::function<void(Trans op,const DenseM_t& R,DenseM_t& S)>& Amult,
369  int rank_guess);
370 
378  void compress(const delem_blocks_t& Aelem);
379 
387  void compress(const elem_blocks_t& Aelem);
388 
403  void mult(Trans op, const DenseM_t& X, DenseM_t& Y) const;
404 
417  void mult(Trans op, const DistM_t& X, DistM_t& Y) const;
418 
425  void factor();
426 
438  long long int solve(const DenseM_t& B, DenseM_t& X) const;
439 
452  long long int solve(const DistM_t& B, DistM_t& X) const;
453 
466  long long int inv_mult(Trans op, const DenseM_t& B, DenseM_t& X) const;
467 
468  void extract_elements
469  (const VecVec_t& I, const VecVec_t& J, std::vector<DistM_t>& B);
470  void extract_elements
471  (const VecVec_t& I, const VecVec_t& J, std::vector<DenseM_t>& B);
472 
485  void extract_elements(const Vec_t& I, const Vec_t& J, DenseM_t& B);
486 
495  DistM_t dense(const BLACSGrid* g) const;
496 
497  void set_sampling_parameter(double sample_param);
498 
499  DenseM_t redistribute_2D_to_1D(const DistM_t& R) const;
500  void redistribute_2D_to_1D(const DistM_t& R2D, DenseM_t& R1D) const;
501  void redistribute_1D_to_2D(const DenseM_t& S1D, DistM_t& S2D) const;
502 
503  DenseMatrix<scalar_t> gather_from_1D(const DenseM_t& A) const;
504  DenseMatrix<scalar_t> all_gather_from_1D(const DenseM_t& A) const;
505  DenseMatrix<scalar_t> scatter_to_1D(const DenseM_t& A) const;
506 
511  const std::vector<int>& perm() const { return perm_; }
516  const std::vector<int>& iperm() const { return iperm_; }
517 
518  private:
519  F2Cptr ho_bf_ = nullptr; // HODLR handle returned by Fortran code
520  F2Cptr options_ = nullptr; // options structure returned by Fortran code
521  F2Cptr stats_ = nullptr; // statistics structure returned by Fortran code
522  F2Cptr msh_ = nullptr; // mesh structure returned by Fortran code
523  F2Cptr kerquant_ = nullptr; // kernel quantities structure returned by Fortran code
524  F2Cptr ptree_ = nullptr; // process tree returned by Fortran code
525  MPI_Fint Fcomm_; // the fortran MPI communicator
526  const MPIComm* c_;
527  int rows_ = 0, cols_ = 0, lrows_ = 0, lvls_ = 0;
528  std::vector<int> perm_, iperm_; // permutation used by the HODLR code
529  std::vector<int> dist_; // begin rows of each rank
530  std::vector<int> leafs_; // leaf sizes of the tree
531 
532  void options_init(const opts_t& opts);
533  void perm_init();
534  void dist_init();
535 
536  template<typename S> friend class ButterflyMatrix;
537  };
538 
539 
540  template<typename scalar_t> struct AelemCommPtrs {
541  const typename HODLRMatrix<scalar_t>::delem_blocks_t* Aelem;
542  const MPIComm* c;
543  const BLACSGrid* gl;
544  const BLACSGrid* g0;
545  };
546 
547  template<typename scalar_t> void HODLR_block_evaluation
548  (int* Ninter, int* Nallrows, int* Nallcols, int* Nalldat_loc,
549  int* allrows, int* allcols, scalar_t* alldat_loc,
550  int* rowids, int* colids, int* pgids, int* Npmap, int* pmaps,
551  void* AC);
552 
553  template<typename scalar_t> void HODLR_block_evaluation_seq
554  (int* Ninter, int* Nallrows, int* Nallcols, int* Nalldat_loc,
555  int* allrows, int* allcols, scalar_t* alldat_loc,
556  int* rowids, int* colids, int* pgids, int* Npmap, int* pmaps,
557  void* f);
558 
559  } // end namespace HODLR
560 } // end namespace strumpack
561 
562 #endif // STRUMPACK_HODLR_MATRIX_HPP
strumpack::BLACSGrid
This is a small wrapper class around a BLACS grid and a BLACS context.
Definition: BLACSGrid.hpp:66
strumpack::HODLR::HODLRMatrix::get_stat
double get_stat(const std::string &name) const
strumpack::HODLR::AelemCommPtrs
Definition: HODLRMatrix.hpp:540
strumpack::HODLR::HODLRMatrix::begin_row
std::size_t begin_row() const
Definition: HODLRMatrix.hpp:266
strumpack::HSS::HSSPartitionTree
The cluster tree, or partition tree that represents the matrix partitioning of an HSS matrix.
Definition: HSSPartitionTree.hpp:65
strumpack::HODLR::HODLRMatrix::end_row
std::size_t end_row() const
Definition: HODLRMatrix.hpp:271
strumpack::MPIComm
Wrapper class around an MPI_Comm object.
Definition: MPIWrapper.hpp:190
strumpack::HODLR::HODLRMatrix::total_memory
std::size_t total_memory() const
Definition: HODLRMatrix.hpp:289
strumpack::HODLR::HODLRMatrix::memory
std::size_t memory() const
Definition: HODLRMatrix.hpp:281
strumpack::HODLR::HODLRMatrix::HODLRMatrix
HODLRMatrix(HODLRMatrix< scalar_t > &&h)
Definition: HODLRMatrix.hpp:229
strumpack
Definition: StrumpackOptions.hpp:42
strumpack::HODLR::HODLRMatrix::rows
std::size_t rows() const
Definition: HODLRMatrix.hpp:251
strumpack::HODLR::HODLRMatrix::perm
const std::vector< int > & perm() const
Definition: HODLRMatrix.hpp:511
strumpack::DistributedMatrix
Definition: CompressedSparseMatrix.hpp:58
strumpack::HODLR::HODLRMatrix::Comm
const MPIComm & Comm() const
Definition: HODLRMatrix.hpp:275
strumpack::DenseMatrixWrapper
Like DenseMatrix, this class represents a matrix, stored in column major format, to allow direct use ...
Definition: DenseMatrix.hpp:990
Kernel.hpp
Definitions of several kernel functions, and helper routines. Also provides driver routines for kerne...
strumpack::HODLR::ExtractionMeta
Definition: HODLRMatrix.hpp:54
strumpack::HODLR::HODLRMatrix::~HODLRMatrix
virtual ~HODLRMatrix()
strumpack::CompressionType::HODLR
@ HODLR
strumpack::DenseMatrix< scalar_t >
strumpack::DistributedMatrixWrapper
Definition: DistributedMatrix.hpp:271
strumpack::MPIComm::all_reduce
T all_reduce(T t, MPI_Op op) const
Definition: MPIWrapper.hpp:483
strumpack::HODLR::HODLRMatrix::factor
void factor()
strumpack::HODLR::HODLRMatrix::total_factor_memory
std::size_t total_factor_memory() const
Definition: HODLRMatrix.hpp:310
strumpack::HODLR::HODLRMatrix::HODLRMatrix
HODLRMatrix()
Definition: HODLRMatrix.hpp:104
strumpack::HODLR::HODLROptions
Class containing several options for the HODLR code and data-structures.
Definition: HODLROptions.hpp:115
strumpack::HODLR::HODLRMatrix::iperm
const std::vector< int > & iperm() const
Definition: HODLRMatrix.hpp:516
strumpack::HODLR::HODLRMatrix::inv_mult
long long int inv_mult(Trans op, const DenseM_t &B, DenseM_t &X) const
strumpack::HODLR::HODLRMatrix::max_rank
std::size_t max_rank() const
Definition: HODLRMatrix.hpp:324
strumpack::HODLR::HODLRMatrix::solve
long long int solve(const DenseM_t &B, DenseM_t &X) const
strumpack::HODLR::HODLRMatrix::factor_memory
std::size_t factor_memory() const
Definition: HODLRMatrix.hpp:299
strumpack::MPIComm::rank
int rank() const
Definition: MPIWrapper.hpp:267
strumpack::HODLR::HODLRMatrix
Hierarchically low-rank matrix representation.
Definition: HODLRMatrix.hpp:78
strumpack::HODLR::HODLRMatrix::dense
DistM_t dense(const BLACSGrid *g) const
strumpack::HODLR::HODLRMatrix::operator=
HODLRMatrix< scalar_t > & operator=(const HODLRMatrix< scalar_t > &h)=delete
strumpack::HODLR::HODLRMatrix::mult
void mult(Trans op, const DenseM_t &X, DenseM_t &Y) const
HODLROptions.hpp
Contains the class holding HODLR matrix options.
strumpack::HODLR::HODLRMatrix::rank
std::size_t rank() const
Definition: HODLRMatrix.hpp:317
strumpack::HODLR::HODLRMatrix::cols
std::size_t cols() const
Definition: HODLRMatrix.hpp:256
strumpack::HODLR::HODLRMatrix::compress
void compress(const std::function< void(Trans op, const DenseM_t &R, DenseM_t &S)> &Amult)
strumpack::kernel::Kernel
Representation of a kernel matrix.
Definition: Kernel.hpp:73
strumpack::HODLR::HODLRMatrix::lrows
std::size_t lrows() const
Definition: HODLRMatrix.hpp:261
strumpack::Trans
Trans
Definition: DenseMatrix.hpp:50
HSSPartitionTree.hpp
This file contains the HSSPartitionTree class definition.
strumpack::CSRGraph
Definition: CompressedSparseMatrix.hpp:56