Hierarchically low-rank matrix representation. More...
#include <HODLRMatrix.hpp>
Public Types | |
using | real_t = typename RealType< scalar_t >::value_type |
using | mult_t = typename std::function< void(Trans, const DenseM_t &, DenseM_t &)> |
using | elem_t = typename std::function< scalar_t(std::size_t i, std::size_t j)> |
using | delem_blocks_t = typename std::function< void(VecVec_t &I, VecVec_t &J, std::vector< DistMW_t > &B, ExtractionMeta &)> |
using | elem_blocks_t = typename std::function< void(VecVec_t &I, VecVec_t &J, std::vector< DenseMW_t > &B, ExtractionMeta &)> |
Public Member Functions | |
HODLRMatrix () | |
HODLRMatrix (const MPIComm &c, kernel::Kernel< real_t > &K, const opts_t &opts) | |
HODLRMatrix (const MPIComm &c, const HSS::HSSPartitionTree &tree, const std::function< scalar_t(int i, int j)> &Aelem, const opts_t &opts) | |
HODLRMatrix (const MPIComm &c, const HSS::HSSPartitionTree &tree, const elem_blocks_t &Aelem, const opts_t &opts) | |
HODLRMatrix (const MPIComm &c, const HSS::HSSPartitionTree &tree, const std::function< void(Trans op, const DenseM_t &R, DenseM_t &S)> &Amult, const opts_t &opts) | |
HODLRMatrix (const MPIComm &c, const HSS::HSSPartitionTree &tree, const opts_t &opts) | |
template<typename integer_t > | |
HODLRMatrix (const MPIComm &c, const HSS::HSSPartitionTree &tree, const CSRGraph< integer_t > &graph, const opts_t &opts) | |
HODLRMatrix (const HODLRMatrix< scalar_t > &h)=delete | |
HODLRMatrix (HODLRMatrix< scalar_t > &&h) | |
virtual | ~HODLRMatrix () |
HODLRMatrix< scalar_t > & | operator= (const HODLRMatrix< scalar_t > &h)=delete |
HODLRMatrix< scalar_t > & | operator= (HODLRMatrix< scalar_t > &&h) |
std::size_t | rows () const |
std::size_t | cols () const |
std::size_t | lrows () const |
std::size_t | begin_row () const |
std::size_t | end_row () const |
const MPIComm & | Comm () const |
std::size_t | memory () const |
std::size_t | total_memory () const |
std::size_t | factor_memory () const |
std::size_t | total_factor_memory () const |
std::size_t | rank () const |
std::size_t | max_rank () const |
double | get_stat (const std::string &name) const |
void | print_stats () |
void | compress (const std::function< void(Trans op, const DenseM_t &R, DenseM_t &S)> &Amult) |
void | compress (const std::function< void(Trans op, const DenseM_t &R, DenseM_t &S)> &Amult, int rank_guess) |
void | compress (const delem_blocks_t &Aelem) |
void | compress (const elem_blocks_t &Aelem) |
void | mult (Trans op, const DenseM_t &X, DenseM_t &Y) const |
void | mult (Trans op, const DistM_t &X, DistM_t &Y) const |
void | factor () |
long long int | solve (const DenseM_t &B, DenseM_t &X) const |
long long int | solve (const DistM_t &B, DistM_t &X) const |
long long int | inv_mult (Trans op, const DenseM_t &B, DenseM_t &X) const |
void | extract_elements (const VecVec_t &I, const VecVec_t &J, std::vector< DistM_t > &B) |
void | extract_elements (const VecVec_t &I, const VecVec_t &J, std::vector< DenseM_t > &B) |
void | extract_elements (const Vec_t &I, const Vec_t &J, DenseM_t &B) |
DistM_t | dense (const BLACSGrid *g) const |
void | set_sampling_parameter (double sample_param) |
DenseM_t | redistribute_2D_to_1D (const DistM_t &R) const |
void | redistribute_2D_to_1D (const DistM_t &R2D, DenseM_t &R1D) const |
void | redistribute_1D_to_2D (const DenseM_t &S1D, DistM_t &S2D) const |
DenseMatrix< scalar_t > | gather_from_1D (const DenseM_t &A) const |
DenseMatrix< scalar_t > | all_gather_from_1D (const DenseM_t &A) const |
DenseMatrix< scalar_t > | scatter_to_1D (const DenseM_t &A) const |
const std::vector< int > & | perm () const |
const std::vector< int > & | iperm () const |
Friends | |
template<typename S > | |
class | ButterflyMatrix |
Hierarchically low-rank matrix representation.
This requires MPI support.
There are 3 different ways to create an HODLRMatrix
scalar_t | Can be double, or std::complex<double>. |
|
inline |
Default constructor, makes an empty 0 x 0 matrix.
strumpack::HODLR::HODLRMatrix< scalar_t >::HODLRMatrix | ( | const MPIComm & | c, |
kernel::Kernel< real_t > & | K, | ||
const opts_t & | opts | ||
) |
Construct an HODLR approximation for the kernel matrix K.
c | MPI communicator, this communicator is copied internally. |
K | Kernel matrix object. The data associated with this kernel will be permuted according to the clustering algorithm selected by the HODLROptions objects. The permutation will be stored in the kernel object. |
opts | object containing a number of HODLR options |
strumpack::HODLR::HODLRMatrix< scalar_t >::HODLRMatrix | ( | const MPIComm & | c, |
const HSS::HSSPartitionTree & | tree, | ||
const std::function< scalar_t(int i, int j)> & | Aelem, | ||
const opts_t & | opts | ||
) |
Construct an HODLR approximation using a routine to evaluate individual matrix elements. This will construct an approximation of a permuted matrix, with the permutation available through this->perm() (and it's inverse this->iperm()). This permutation is applied symmetrically to rows and columns.
c | MPI communicator, this communicator is copied internally. |
t | tree specifying the HODLR matrix partitioning |
Aelem | Routine, std::function, which can also be a lambda function or a functor (class object implementing the member "scalar_t operator()(int i, int j)"), that evaluates/returns the matrix element A(i,j) |
opts | object containing a number of HODLR options |
strumpack::HODLR::HODLRMatrix< scalar_t >::HODLRMatrix | ( | const MPIComm & | c, |
const HSS::HSSPartitionTree & | tree, | ||
const elem_blocks_t & | Aelem, | ||
const opts_t & | opts | ||
) |
Construct an HODLR approximation using a routine to evaluate multiple sub-blocks of the matrix at once. This will construct an approximation of a permuted matrix, with the permutation available through this->perm() (and it's inverse this->iperm()). This permutation is applied symmetrically to rows and columns.
c | MPI communicator, this communicator is copied internally. |
t | tree specifying the HODLR matrix partitioning |
Aelem | Routine, std::function, which can also be a lambda function or a functor. This should have the signature: void(std::vector<std::vector<std::size_t>>& I, std::vector<std::vector<std::size_t>>& J, std::vector<DenseMatrixWrapper<scalar_t>>& B, ExtractionMeta&); The ExtractionMeta object can be ignored. |
opts | object containing a number of HODLR options |
strumpack::HODLR::HODLRMatrix< scalar_t >::HODLRMatrix | ( | const MPIComm & | c, |
const HSS::HSSPartitionTree & | tree, | ||
const std::function< void(Trans op, const DenseM_t &R, DenseM_t &S)> & | Amult, | ||
const opts_t & | opts | ||
) |
Construct an HODLR matrix using a specified HODLR tree and matrix-vector multiplication routine.
c | MPI communicator, this communicator is copied internally. |
t | tree specifying the HODLR matrix partitioning |
Amult | Routine for the matrix-vector product. Trans op argument will be N, T or C for none, transpose or complex conjugate. The const DenseM_t& argument is the the random matrix R, and the final DenseM_t& argument S is what the user routine should compute as A*R, A^t*R or A^c*R. S will already be allocated. |
opts | object containing a number of options for HODLR compression |
strumpack::HODLR::HODLRMatrix< scalar_t >::HODLRMatrix | ( | const MPIComm & | c, |
const HSS::HSSPartitionTree & | tree, | ||
const opts_t & | opts | ||
) |
Construct an HODLR matrix using a specified HODLR tree. After construction, the HODLR matrix will be empty, and can be filled by calling one of the compress member routines.
c | MPI communicator, this communicator is copied internally. |
t | tree specifying the HODLR matrix partitioning |
opts | object containing a number of options for HODLR compression |
strumpack::HODLR::HODLRMatrix< scalar_t >::HODLRMatrix | ( | const MPIComm & | c, |
const HSS::HSSPartitionTree & | tree, | ||
const CSRGraph< integer_t > & | graph, | ||
const opts_t & | opts | ||
) |
Construct an HODLR matrix using a specified HODLR tree. After construction, the HODLR matrix will be empty, and can be filled by calling one of the compress member routines.
c | MPI communicator, this communicator is copied internally. |
t | tree specifying the HODLR matrix partitioning |
graph | connectivity info for the dofs |
opts | object containing a number of options for HODLR compression |
|
delete |
Copy constructor is not supported.
|
inline |
Move constructor.
h | HODLRMatrix to move from, will be emptied. |
|
virtual |
Virtual destructor.
|
inline |
Return the first row of the local rows owned by this process.
|
inline |
Return the number of columns in the matrix.
|
inline |
Return MPI communicator wrapper object.
void strumpack::HODLR::HODLRMatrix< scalar_t >::compress | ( | const delem_blocks_t & | Aelem | ) |
Construct the compressed HODLR representation of the matrix, using only a element evaluation (multiple blocks at once).
Aelem | element extraction routine, extracting multiple blocks at once. |
void strumpack::HODLR::HODLRMatrix< scalar_t >::compress | ( | const elem_blocks_t & | Aelem | ) |
Construct the compressed HODLR representation of the matrix, using only a element evaluation (multiple blocks at once).
Aelem | element extraction routine, extracting multiple blocks at once. |
void strumpack::HODLR::HODLRMatrix< scalar_t >::compress | ( | const std::function< void(Trans op, const DenseM_t &R, DenseM_t &S)> & | Amult | ) |
Construct the compressed HODLR representation of the matrix, using only a matrix-(multiple)vector multiplication routine.
Amult | Routine for the matrix-vector product. Trans op argument will be N, T or C for none, transpose or complex conjugate. The const DenseM_t& argument is the the random matrix R, and the final DenseM_t& argument S is what the user routine should compute as A*R, A^t*R or A^c*R. S will already be allocated. |
void strumpack::HODLR::HODLRMatrix< scalar_t >::compress | ( | const std::function< void(Trans op, const DenseM_t &R, DenseM_t &S)> & | Amult, |
int | rank_guess | ||
) |
Construct the compressed HODLR representation of the matrix, using only a matrix-(multiple)vector multiplication routine.
Amult | Routine for the matrix-vector product. Trans op argument will be N, T or C for none, transpose or complex conjugate. The const DenseM_t& argument is the the random matrix R, and the final DenseM_t& argument S is what the user routine should compute as A*R, A^t*R or A^c*R. S will already be allocated. |
rank_guess | Initial guess for the rank |
DistM_t strumpack::HODLR::HODLRMatrix< scalar_t >::dense | ( | const BLACSGrid * | g | ) | const |
Create a dense matrix from this HODLR compressed matrix. This is mainly for debugging.
g | BLACSFGrid to be used for the output |
|
inline |
Return last row (+1) of the local rows (begin_rows()+lrows())
void strumpack::HODLR::HODLRMatrix< scalar_t >::extract_elements | ( | const Vec_t & | I, |
const Vec_t & | J, | ||
DenseM_t & | B | ||
) |
Extract a submatrix defined by index sets I (rows) and J (columns), and put the result in matrix B (only on the root!).
I | set of row indices, needs to be specified on all ranks! |
J | set of column indices, needs to be specified on all ranks! |
B | output, the extracted elements. B should have the correct size B.rows() == I.size() and B.cols() == J.size() |
void strumpack::HODLR::HODLRMatrix< scalar_t >::factor | ( | ) |
|
inline |
double strumpack::HODLR::HODLRMatrix< scalar_t >::get_stat | ( | const std::string & | name | ) | const |
Get certain statistics about the HODLR matrix. See the HODLR code at https://github.com/liuyangzhuan/ButterflyPACK for more info.
name | fi : Rank_max, Mem_Factor, Mem_Fill, Flop_Fill, Flop_Factor, Flop_C_Mult |
long long int strumpack::HODLR::HODLRMatrix< scalar_t >::inv_mult | ( | Trans | op, |
const DenseM_t & | B, | ||
DenseM_t & | X | ||
) | const |
Solve a system of linear equations op(A)*X=B, with possibly multiple right-hand sides, where op can be none, transpose or complex conjugate.
B | Right hand side. This is the local part of the distributed matrix B. Should be B.rows() == this.lrows(). |
X | Result, should be X.cols() == B.cols(), X.rows() == this.lrows(). X should be allocated. |
|
inline |
The inverse permutation for the matrix, which is applied to both rows and columns. This is 1-based!.
|
inline |
Return the number of local rows, owned by this process.
|
inline |
|
inline |
Return the memory for this HODLR matrix, on this rank, in bytes.
void strumpack::HODLR::HODLRMatrix< scalar_t >::mult | ( | Trans | op, |
const DenseM_t & | X, | ||
DenseM_t & | Y | ||
) | const |
Multiply this HODLR matrix with a dense matrix: Y = op(this)*X, where op can be none, transpose or complex conjugate. X and Y are the local parts of block-row distributed matrices. The number of rows in X and Y should correspond to the distribution of this HODLR matrix.
op | Transpose, conjugate, or none. |
X | Right-hand side matrix. This is the local part of the distributed matrix X. Should be X.rows() == this.lrows(). |
Y | Result, should be Y.cols() == X.cols(), Y.rows() == this.lrows() |
void strumpack::HODLR::HODLRMatrix< scalar_t >::mult | ( | Trans | op, |
const DistM_t & | X, | ||
DistM_t & | Y | ||
) | const |
Multiply this HODLR matrix with a dense matrix: Y = op(this)*X, where op can be none, transpose or complex conjugate. X and Y are in 2D block cyclic distribution.
op | Transpose, conjugate, or none. |
X | Right-hand side matrix. Should be X.rows() == this.rows(). |
Y | Result, should be Y.cols() == X.cols(), Y.rows() == this.rows() |
|
delete |
Copy assignment operator is not supported.
HODLRMatrix<scalar_t>& strumpack::HODLR::HODLRMatrix< scalar_t >::operator= | ( | HODLRMatrix< scalar_t > && | h | ) |
Move assignment operator.
h | HODLRMatrix to move from, will be emptied. |
|
inline |
The permutation for the matrix, which is applied to both rows and columns. This is 1-based!.
|
inline |
Return the maximal rank encountered in this HODLR matrix.
|
inline |
Return the number of rows in the matrix.
long long int strumpack::HODLR::HODLRMatrix< scalar_t >::solve | ( | const DenseM_t & | B, |
DenseM_t & | X | ||
) | const |
Solve a system of linear equations A*X=B, with possibly multiple right-hand sides.
B | Right hand side. This is the local part of the distributed matrix B. Should be B.rows() == this.lrows(). |
X | Result, should be X.cols() == B.cols(), X.rows() == this.lrows(). X should be allocated. |
long long int strumpack::HODLR::HODLRMatrix< scalar_t >::solve | ( | const DistM_t & | B, |
DistM_t & | X | ||
) | const |
Solve a system of linear equations A*X=B, with possibly multiple right-hand sides. X and B are in 2D block cyclic distribution.
B | Right hand side. This is the local part of the distributed matrix B. Should be B.rows() == this.rows(). |
X | Result, should be X.cols() == B.cols(), X.rows() == this.rows(). X should be allocated. |
|
inline |
|
inline |
Return the total memory for this HODLR matrix, summed over all ranks, in bytes.