strumpack::structured Namespace Reference

Classes

class  ClusterTree
 The cluster tree, or partition tree that represents the partitioning of the rows or columns of a hierarchical matrix. More...
 
class  StructuredMatrix
 Class to represent a structured matrix. This is the abstract base class for several types of structured matrices. More...
 
class  StructuredOptions
 Class containing several options for the StructuredMatrix code and data-structures. More...
 

Typedefs

template<typename scalar_t >
using extract_t = std::function< scalar_t(std::size_t i, std::size_t j)>
 
template<typename scalar_t >
using extract_block_t = std::function< void(const std::vector< std::size_t > &I, const std::vector< std::size_t > &J, DenseMatrix< scalar_t > &B)>
 
template<typename scalar_t >
using mult_t = std::function< void(Trans op, const DenseMatrix< scalar_t > &R, DenseMatrix< scalar_t > &S)>
 
template<typename scalar_t >
using extract_dist_block_t = std::function< void(const std::vector< std::size_t > &I, const std::vector< std::size_t > &J, DistributedMatrix< scalar_t > &B)>
 
template<typename scalar_t >
using mult_2d_t = std::function< void(Trans op, const DistributedMatrix< scalar_t > &R, DistributedMatrix< scalar_t > &S)>
 
template<typename scalar_t >
using mult_1d_t = std::function< void(Trans op, const DenseMatrix< scalar_t > &R, DenseMatrix< scalar_t > &S, const std::vector< int > &rdist, const std::vector< int > &cdist)>
 
using admissibility_t = DenseMatrix< bool >
 

Enumerations

enum class  Type : int {
  HSS = 0 , BLR , HODLR , HODBF ,
  BUTTERFLY , LR , LOSSY , LOSSLESS
}
 

Functions

template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > construct_from_dense (const DenseMatrix< scalar_t > &A, const StructuredOptions< scalar_t > &opts, const structured::ClusterTree *row_tree=nullptr, const structured::ClusterTree *col_tree=nullptr)
 
template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > construct_from_dense (int rows, int cols, const scalar_t *A, int ldA, const StructuredOptions< scalar_t > &opts, const structured::ClusterTree *row_tree=nullptr, const structured::ClusterTree *col_tree=nullptr)
 
template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > construct_from_elements (int rows, int cols, const extract_t< scalar_t > &A, const StructuredOptions< scalar_t > &opts, const structured::ClusterTree *row_tree=nullptr, const structured::ClusterTree *col_tree=nullptr)
 
template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > construct_from_elements (int rows, int cols, const extract_block_t< scalar_t > &A, const StructuredOptions< scalar_t > &opts, const structured::ClusterTree *row_tree=nullptr, const structured::ClusterTree *col_tree=nullptr)
 
template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > construct_matrix_free (int rows, int cols, const mult_t< scalar_t > &Amult, const StructuredOptions< scalar_t > &opts, const structured::ClusterTree *row_tree=nullptr, const structured::ClusterTree *col_tree=nullptr)
 
template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > construct_partially_matrix_free (int rows, int cols, const mult_t< scalar_t > &Amult, const extract_block_t< scalar_t > &Aelem, const StructuredOptions< scalar_t > &opts, const structured::ClusterTree *row_tree=nullptr, const structured::ClusterTree *col_tree=nullptr)
 
template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > construct_partially_matrix_free (int rows, int cols, const mult_t< scalar_t > &Amult, const extract_t< scalar_t > &Aelem, const StructuredOptions< scalar_t > &opts, const structured::ClusterTree *row_tree=nullptr, const structured::ClusterTree *col_tree=nullptr)
 
template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > construct_from_dense (const DistributedMatrix< scalar_t > &A, const StructuredOptions< scalar_t > &opts, const structured::ClusterTree *row_tree=nullptr, const structured::ClusterTree *col_tree=nullptr)
 
template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > construct_from_elements (const MPIComm &comm, const BLACSGrid *grid, int rows, int cols, const extract_t< scalar_t > &A, const StructuredOptions< scalar_t > &opts, const structured::ClusterTree *row_tree=nullptr, const structured::ClusterTree *col_tree=nullptr)
 
template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > construct_from_elements (const MPIComm &comm, const BLACSGrid *grid, int rows, int cols, const extract_dist_block_t< scalar_t > &A, const StructuredOptions< scalar_t > &opts, const structured::ClusterTree *row_tree=nullptr, const structured::ClusterTree *col_tree=nullptr)
 
template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > construct_matrix_free (const MPIComm &comm, const BLACSGrid *g, int rows, int cols, const mult_2d_t< scalar_t > &Amult, const StructuredOptions< scalar_t > &opts, const structured::ClusterTree *row_tree=nullptr, const structured::ClusterTree *col_tree=nullptr)
 
template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > construct_matrix_free (const MPIComm &comm, int rows, int cols, const mult_1d_t< scalar_t > &Amult, const StructuredOptions< scalar_t > &opts, const structured::ClusterTree *row_tree=nullptr, const structured::ClusterTree *col_tree=nullptr)
 
template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > construct_partially_matrix_free (const BLACSGrid *grid, int rows, int cols, const mult_2d_t< scalar_t > &Amult, const extract_dist_block_t< scalar_t > &Aelem, const StructuredOptions< scalar_t > &opts, const structured::ClusterTree *row_tree=nullptr, const structured::ClusterTree *col_tree=nullptr)
 
template<typename real_t >
real_t default_structured_rel_tol ()
 
template<typename real_t >
real_t default_structured_abs_tol ()
 
template<>
float default_structured_rel_tol ()
 
template<>
float default_structured_abs_tol ()
 
std::string get_name (Type a)
 

Detailed Description

Namespace containing the StructuredMatrix class, a general abstract interface to several other rank-structured matrix representations.

Typedef Documentation

◆ admissibility_t

Type to specify admissibility of individual sub-blocks.

◆ extract_block_t

template<typename scalar_t >
using strumpack::structured::extract_block_t = typedef std::function <void(const std::vector<std::size_t>& I, const std::vector<std::size_t>& J, DenseMatrix<scalar_t>& B)>

Type for block element extraction routine. This can be implemented as a lambda function, a functor or a funtion pointer.

Parameters
Ivector of row indices
Jvector of column indices
Bsubmatrix (I, J) to be computed by user code

◆ extract_dist_block_t

template<typename scalar_t >
using strumpack::structured::extract_dist_block_t = typedef std::function <void(const std::vector<std::size_t>& I, const std::vector<std::size_t>& J, DistributedMatrix<scalar_t>& B)>

Type for block element extraction routine. This can be implemented as a lambda function, a functor or a funtion pointer. The sub-matrix should be put into a 2d block cyclicly distributed matrix.

Parameters
Ivector of row indices
Jvector of column indices
B2d block cyclic matrix, submatrix (I, J) to be computed by user code

◆ extract_t

template<typename scalar_t >
using strumpack::structured::extract_t = typedef std::function <scalar_t(std::size_t i, std::size_t j)>

Type for element extraction routine. This can be implemented as a lambda function, a functor or a funtion pointer.

Parameters
ivector of row indices
jvector of column indices
Returns
element i,j of the matrix

◆ mult_1d_t

template<typename scalar_t >
using strumpack::structured::mult_1d_t = typedef std::function <void(Trans op, const DenseMatrix<scalar_t>& R, DenseMatrix<scalar_t>& S, const std::vector<int>& rdist, const std::vector<int>& cdist)>

Type for matrix multiplication routine with 1D block row distribution. The block row/column distribution of the matrix is given by the rdist and cdist vectors, with processor p owning rows [rdist[p],rdist[p+1]). cdist is only needed for non-square matrices. S is distributed according to rdist if t==Trans::N, else cdist. R is distributed according to cdist of t==Trans::N, else rdist.

Parameters
opwhether to compute the transposed/conjugate product
Rrandom matrix passed to user code
Ssample matrix to be computed by user code:
S = A*R if op is Trans::N
S = A*T*R if op is Trans::T
S = A^C*R if op is Trans::C
rdistMatrix row distribution, same on all ranks
cdistMatrix column distribution, same on all ranks

◆ mult_2d_t

template<typename scalar_t >
using strumpack::structured::mult_2d_t = typedef std::function <void(Trans op, const DistributedMatrix<scalar_t>& R, DistributedMatrix<scalar_t>& S)>

Type for matrix multiplication routine with 2D block cyclic distribution. This can be implemented as a lambda function, a functor or a funtion pointer.

Parameters
opwhether to compute the transposed/conjugate product
Rrandom matrix passed to user code
Ssample matrix to be computed by user code:
S = A*R if op is Trans::N
S = A*T*R if op is Trans::T
S = A^C*R if op is Trans::C

◆ mult_t

template<typename scalar_t >
using strumpack::structured::mult_t = typedef std::function <void(Trans op, const DenseMatrix<scalar_t>& R, DenseMatrix<scalar_t>& S)>

Type for matrix multiplication routine. This can be implemented as a lambda function, a functor or a funtion pointer.

Parameters
opwhether to compute the transposed/conjugate product
Rrandom matrix passed to user code
Ssample matrix to be computed by user code:
S = A*R if op is Trans::N
S = A*T*R if op is Trans::T
S = A^C*R if op is Trans::C

Enumeration Type Documentation

◆ Type

enum class strumpack::structured::Type : int
strong

Enumeration of possible structured matrix types.

Enumerator
HSS 

Hierarchically Semi-Separable, see HSS::HSSMatrix and HSS::HSSMatrixMPI

BLR 

Block Low Rank, see BLR::BLRMatrix and BLR::BLRMatrixMPI

HODLR 

Hierarchically Off-Diagonal Low Rank, see HODLR::HODLRMatrix. Does not support float or std::complex<float>.

HODBF 

Hierarchically Off-Diagonal Butterfly, implemented as HODLR::HODLRMatrix. Does not support float or std::complex<float>.

BUTTERFLY 

Butterfly matrix, implemented as HODLR::ButterflyMatrix. Does not support float or std::complex<float>.

LR 

Low rank matrix, implemented as HODLR::ButterflyMatrix. Does not support float or std::complex<float>.

LOSSY 

Lossy compression matrix

LOSSLESS 

Lossless compressed matrix

Function Documentation

◆ construct_from_dense() [1/3]

template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > strumpack::structured::construct_from_dense ( const DenseMatrix< scalar_t > &  A,
const StructuredOptions< scalar_t > &  opts,
const structured::ClusterTree row_tree = nullptr,
const structured::ClusterTree col_tree = nullptr 
)

Construct a StructuredMatrix from a DenseMatrix<scalar_t>. This construction, taking a dense matrix should be valid for all StructuredMatrix types (see StructuredMatrix::Type). However, faster and more memory efficient methods of constructing StructuredMatrix representations are available, see StructuredMatrix::construct_from_elements, StructuredMatrix::construct_matrix_free and StructuredMatrix::construct_partially_matrix_free.

Template Parameters
scalar_tprecision of input matrix, and of constructed StructuredMatrix. Note that not all types support every all precisions. See StructuredMatrix::Type.
Parameters
AInput dense matrix, will not be modified
optsOptions object
row_treeoptional clustertree for the rows, see also strumpack::binary_tree_clustering
col_treeoptional clustertree for the columns. If the matrix is square, this does not need to be specified.
Returns
std::unique_ptr holding a pointer to a StructuredMatrix of the requested StructuredMatrix::Type
Exceptions
std::invalid_argumentIf the operatation is not supported for the type of structured::StructuredMatrix, if the type requires a square matrix and the input is not square, if the structured::StructuredMatrix type requires MPI.
std::logic_errorIf the operation is not implemented yet
std::runtime_errorIf the operation requires a third party library which was not enabled when configuring STRUMPACK.
See also
strumpack::binary_tree_clustering, construct_from_dense construct_from_elements, construct_matrix_free and construct_partially_matrix_free

◆ construct_from_dense() [2/3]

template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > strumpack::structured::construct_from_dense ( const DistributedMatrix< scalar_t > &  A,
const StructuredOptions< scalar_t > &  opts,
const structured::ClusterTree row_tree = nullptr,
const structured::ClusterTree col_tree = nullptr 
)

Construct a StructuredMatrix from a 2D block-cyclicly distributed matrix (ScaLAPACK layout).

Template Parameters
scalar_tprecision of input matrix, and of constructed StructuredMatrix. Note that not all types support every all precisions. See StructuredMatrix::Type.
Parameters
A2D block cyclic distributed matrix (ScaLAPACK layout).
optsOptions object
row_treeoptional clustertree for the rows, see also strumpack::binary_tree_clustering
col_treeoptional clustertree for the columns. If the matrix is square, this does not need to be specified.
Returns
std::unique_ptr holding a pointer to a StructuredMatrix of the requested StructuredMatrix::Type
Exceptions
std::invalid_argumentIf the operatation is not supported for the type of structured::StructuredMatrix, if the type requires a square matrix and the input is not square, if the structured::StructuredMatrix type requires MPI.
std::logic_errorIf the operation is not implemented yet
std::runtime_errorIf the operation requires a third party library which was not enabled when configuring STRUMPACK.
See also
strumpack::binary_tree_clustering, construct_from_dense construct_from_elements, construct_matrix_free and construct_partially_matrix_free

◆ construct_from_dense() [3/3]

template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > strumpack::structured::construct_from_dense ( int  rows,
int  cols,
const scalar_t *  A,
int  ldA,
const StructuredOptions< scalar_t > &  opts,
const structured::ClusterTree row_tree = nullptr,
const structured::ClusterTree col_tree = nullptr 
)

Construct a StructuredMatrix from a DenseMatrix<scalar_t>. This construction, taking a dense matrix should be valid for all StructuredMatrix types (see StructuredMatrix::Type). However, faster and more memory efficient methods of constructing StructuredMatrix representations are available, see StructuredMatrix::construct_from_elements, StructuredMatrix::construct_matrix_free and StructuredMatrix::construct_partially_matrix_free.

Template Parameters
scalar_tprecision of input matrix, and of constructed StructuredMatrix. Note that not all types support every all precisions. See StructuredMatrix::Type.
Parameters
rowsnumber of rows in matrix to be constructed
colsnumber of columnss in matrix to be constructed
Apointer to matrix A data
ldAleading dimension of A
optsOptions object
row_treeoptional clustertree for the rows, see also strumpack::binary_tree_clustering
col_treeoptional clustertree for the columns. If the matrix is square, this does not need to be specified.
Returns
std::unique_ptr holding a pointer to a StructuredMatrix of the requested StructuredMatrix::Type
Exceptions
std::invalid_argumentIf the operatation is not supported for the type of structured::StructuredMatrix, if the type requires a square matrix and the input is not square, if the structured::StructuredMatrix type requires MPI.
std::logic_errorIf the operation is not implemented yet
std::runtime_errorIf the operation requires a third party library which was not enabled when configuring STRUMPACK.
See also
strumpack::binary_tree_clustering, construct_from_dense construct_from_elements, construct_matrix_free and construct_partially_matrix_free

◆ construct_from_elements() [1/4]

template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > strumpack::structured::construct_from_elements ( const MPIComm comm,
const BLACSGrid grid,
int  rows,
int  cols,
const extract_dist_block_t< scalar_t > &  A,
const StructuredOptions< scalar_t > &  opts,
const structured::ClusterTree row_tree = nullptr,
const structured::ClusterTree col_tree = nullptr 
)

Construct a StructuredMatrix using a routine to extract a matrix sub-block.

Template Parameters
scalar_tprecision of input matrix, and of constructed StructuredMatrix. Note that not all types support every all precisions. See StructuredMatrix::Type.
Parameters
commMPI communicator (wrapper class)
gridOnly used for HSS construction, can be NULL otherwise
rowsNumber of rows of matrix to be constructed.
colsNumber of columns of matrix to be constructed.
AMatrix sub-block extraction routine.
optsOptions object
row_treeoptional clustertree for the rows, see also strumpack::binary_tree_clustering
col_treeoptional clustertree for the columns. If the matrix is square, this does not need to be specified.
Returns
std::unique_ptr holding a pointer to a StructuredMatrix of the requested StructuredMatrix::Type
Exceptions
std::invalid_argumentIf the operatation is not supported for the type of structured::StructuredMatrix, if the type requires a square matrix and the input is not square, if the structured::StructuredMatrix type requires MPI.
std::logic_errorIf the operation is not implemented yet
std::runtime_errorIf the operation requires a third party library which was not enabled when configuring STRUMPACK.
See also
strumpack::binary_tree_clustering, construct_from_dense construct_from_elements, construct_matrix_free and construct_partially_matrix_free

◆ construct_from_elements() [2/4]

template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > strumpack::structured::construct_from_elements ( const MPIComm comm,
const BLACSGrid grid,
int  rows,
int  cols,
const extract_t< scalar_t > &  A,
const StructuredOptions< scalar_t > &  opts,
const structured::ClusterTree row_tree = nullptr,
const structured::ClusterTree col_tree = nullptr 
)

Construct a StructuredMatrix using an element extraction routine.

Template Parameters
scalar_tprecision of input matrix, and of constructed StructuredMatrix. Note that not all types support every all precisions. See StructuredMatrix::Type.
Parameters
commMPI communicator (wrapper class)
gridOnly used for HSS construction, can be NULL otherwise
rowsNumber of rows of matrix to be constructed.
colsNumber of columns of matrix to be constructed.
AMatrix element extraction routine.
optsOptions object
row_treeoptional clustertree for the rows, see also strumpack::binary_tree_clustering
col_treeoptional clustertree for the columns. If the matrix is square, this does not need to be specified.
Returns
std::unique_ptr holding a pointer to a StructuredMatrix of the requested StructuredMatrix::Type
Exceptions
std::invalid_argumentIf the operatation is not supported for the type of structured::StructuredMatrix, if the type requires a square matrix and the input is not square, if the structured::StructuredMatrix type requires MPI.
std::logic_errorIf the operation is not implemented yet
std::runtime_errorIf the operation requires a third party library which was not enabled when configuring STRUMPACK.
See also
strumpack::binary_tree_clustering, construct_from_dense construct_from_elements, construct_matrix_free and construct_partially_matrix_free

◆ construct_from_elements() [3/4]

template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > strumpack::structured::construct_from_elements ( int  rows,
int  cols,
const extract_block_t< scalar_t > &  A,
const StructuredOptions< scalar_t > &  opts,
const structured::ClusterTree row_tree = nullptr,
const structured::ClusterTree col_tree = nullptr 
)

Construct a StructuredMatrix using a routine to extract a sub-block from the matrix.

Template Parameters
scalar_tprecision of input matrix, and of constructed StructuredMatrix. Note that not all types support every all precisions. See StructuredMatrix::Type.
Parameters
rowsNumber of rows of matrix to be constructed.
colsNumber of columns of matrix to be constructed.
AMatrix block extraction routine.
optsOptions object
row_treeoptional clustertree for the rows, see also strumpack::binary_tree_clustering
col_treeoptional clustertree for the columns. If the matrix is square, this does not need to be specified.
Returns
std::unique_ptr holding a pointer to a StructuredMatrix of the requested StructuredMatrix::Type
Exceptions
std::invalid_argumentIf the operatation is not supported for the type of structured::StructuredMatrix, if the type requires a square matrix and the input is not square, if the structured::StructuredMatrix type requires MPI.
std::logic_errorIf the operation is not implemented yet
std::runtime_errorIf the operation requires a third party library which was not enabled when configuring STRUMPACK.
See also
strumpack::binary_tree_clustering, construct_from_dense construct_from_elements, construct_matrix_free and construct_partially_matrix_free

◆ construct_from_elements() [4/4]

template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > strumpack::structured::construct_from_elements ( int  rows,
int  cols,
const extract_t< scalar_t > &  A,
const StructuredOptions< scalar_t > &  opts,
const structured::ClusterTree row_tree = nullptr,
const structured::ClusterTree col_tree = nullptr 
)

Construct a StructuredMatrix using a routine, provided by the user, to extract elements from the matrix.

Template Parameters
scalar_tprecision of input matrix, and of constructed StructuredMatrix. Note that not all types support every all precisions. See StructuredMatrix::Type.
Parameters
rowsnumber of rows in matrix to be constructed
colsnumber of columnss in matrix to be constructed
Aelement extraction routine
optsOptions object
row_treeoptional clustertree for the rows, see also strumpack::binary_tree_clustering
col_treeoptional clustertree for the columns. If the matrix is square, this does not need to be specified.
Returns
std::unique_ptr holding a pointer to a StructuredMatrix of the requested StructuredMatrix::Type
Exceptions
std::invalid_argumentIf the operatation is not supported for the type of structured::StructuredMatrix, if the type requires a square matrix and the input is not square, if the structured::StructuredMatrix type requires MPI.
std::logic_errorIf the operation is not implemented yet
std::runtime_errorIf the operation requires a third party library which was not enabled when configuring STRUMPACK.
See also
strumpack::binary_tree_clustering, construct_from_dense construct_from_elements, construct_matrix_free and construct_partially_matrix_free

◆ construct_matrix_free() [1/3]

template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > strumpack::structured::construct_matrix_free ( const MPIComm comm,
const BLACSGrid g,
int  rows,
int  cols,
const mult_2d_t< scalar_t > &  Amult,
const StructuredOptions< scalar_t > &  opts,
const structured::ClusterTree row_tree = nullptr,
const structured::ClusterTree col_tree = nullptr 
)

Construct a StructuredMatrix using only a matrix-vector multiplication routine.

Parameters
commMPI communicator (wrapper class)
rowsNumber of rows of matrix to be constructed.
colsNumber of columns of matrix to be constructed.
AmultMatrix-(multi)vector multiplication routine (using 2d block cyclic layout for vectors).
optsOptions object
row_treeoptional clustertree for the rows, see also strumpack::binary_tree_clustering
col_treeoptional clustertree for the columns. If the matrix is square, this does not need to be specified.
Returns
std::unique_ptr holding a pointer to a StructuredMatrix of the requested StructuredMatrix::Type
Exceptions
std::invalid_argumentIf the operatation is not supported for the type of structured::StructuredMatrix, if the type requires a square matrix and the input is not square, if the structured::StructuredMatrix type requires MPI.
std::logic_errorIf the operation is not implemented yet
std::runtime_errorIf the operation requires a third party library which was not enabled when configuring STRUMPACK.
See also
strumpack::binary_tree_clustering, construct_from_dense construct_from_elements, construct_matrix_free and construct_partially_matrix_free

◆ construct_matrix_free() [2/3]

template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > strumpack::structured::construct_matrix_free ( const MPIComm comm,
int  rows,
int  cols,
const mult_1d_t< scalar_t > &  Amult,
const StructuredOptions< scalar_t > &  opts,
const structured::ClusterTree row_tree = nullptr,
const structured::ClusterTree col_tree = nullptr 
)

Construct a StructuredMatrix using only a matrix-vector multiplication routine.

Parameters
commMPI communicator (wrapper class)
rowsNumber of rows of matrix to be constructed.
colsNumber of columns of matrix to be constructed.
AmultMatrix-(multi)vector multiplication routine (using 1d block row distribution for vectors).
optsOptions object
row_treeoptional clustertree for the rows, see also strumpack::binary_tree_clustering
col_treeoptional clustertree for the columns. If the matrix is square, this does not need to be specified.
Returns
std::unique_ptr holding a pointer to a StructuredMatrix of the requested StructuredMatrix::Type
Exceptions
std::invalid_argumentIf the operatation is not supported for the type of structured::StructuredMatrix, if the type requires a square matrix and the input is not square, if the structured::StructuredMatrix type requires MPI.
std::logic_errorIf the operation is not implemented yet
std::runtime_errorIf the operation requires a third party library which was not enabled when configuring STRUMPACK.
See also
strumpack::binary_tree_clustering, construct_from_dense construct_from_elements, construct_matrix_free and construct_partially_matrix_free

◆ construct_matrix_free() [3/3]

template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > strumpack::structured::construct_matrix_free ( int  rows,
int  cols,
const mult_t< scalar_t > &  Amult,
const StructuredOptions< scalar_t > &  opts,
const structured::ClusterTree row_tree = nullptr,
const structured::ClusterTree col_tree = nullptr 
)

Construct a StructuredMatrix using only a matrix-vector multiplication routine.

Template Parameters
scalar_tprecision of input matrix, and of constructed StructuredMatrix. Note that not all types support every all precisions. See StructuredMatrix::Type.
Parameters
rowsNumber of rows of matrix to be constructed.
colsNumber of columns of matrix to be constructed.
AmultMatrix-(multi)vector multiplication routine.
optsOptions object
row_treeoptional clustertree for the rows, see also strumpack::binary_tree_clustering
col_treeoptional clustertree for the columns. If the matrix is square, this does not need to be specified.
Returns
std::unique_ptr holding a pointer to a StructuredMatrix of the requested StructuredMatrix::Type
Exceptions
std::invalid_argumentIf the operatation is not supported for the type of structured::StructuredMatrix, if the type requires a square matrix and the input is not square, if the structured::StructuredMatrix type requires MPI.
std::logic_errorIf the operation is not implemented yet
std::runtime_errorIf the operation requires a third party library which was not enabled when configuring STRUMPACK.
See also
strumpack::binary_tree_clustering, construct_from_dense construct_from_elements, construct_matrix_free and construct_partially_matrix_free

◆ construct_partially_matrix_free() [1/3]

template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > strumpack::structured::construct_partially_matrix_free ( const BLACSGrid grid,
int  rows,
int  cols,
const mult_2d_t< scalar_t > &  Amult,
const extract_dist_block_t< scalar_t > &  Aelem,
const StructuredOptions< scalar_t > &  opts,
const structured::ClusterTree row_tree = nullptr,
const structured::ClusterTree col_tree = nullptr 
)

Construct a StructuredMatrix using both a matrix-vector multiplication routine and a routine to extract individual matrix elements.

Template Parameters
scalar_tprecision of input matrix, and of constructed StructuredMatrix. Note that not all types support every all precisions. See StructuredMatrix::Type.
Parameters
gridBLACS grid (wrapper class)
rowsNumber of rows of matrix to be constructed.
colsNumber of columns of matrix to be constructed.
AmultMatrix-(multi)vector multiplication routine (using 2d block cyclic distribution for vectors).
AelemMatrix block extraction routine.
optsOptions object
row_treeoptional clustertree for the rows, see also strumpack::binary_tree_clustering
col_treeoptional clustertree for the columns. If the matrix is square, this does not need to be specified.
Returns
std::unique_ptr holding a pointer to a StructuredMatrix of the requested StructuredMatrix::Type
Exceptions
std::invalid_argumentIf the operatation is not supported for the type of structured::StructuredMatrix, if the type requires a square matrix and the input is not square, if the structured::StructuredMatrix type requires MPI.
std::logic_errorIf the operation is not implemented yet
std::runtime_errorIf the operation requires a third party library which was not enabled when configuring STRUMPACK.
See also
strumpack::binary_tree_clustering, construct_from_dense construct_from_elements, construct_matrix_free and construct_partially_matrix_free

◆ construct_partially_matrix_free() [2/3]

template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > strumpack::structured::construct_partially_matrix_free ( int  rows,
int  cols,
const mult_t< scalar_t > &  Amult,
const extract_block_t< scalar_t > &  Aelem,
const StructuredOptions< scalar_t > &  opts,
const structured::ClusterTree row_tree = nullptr,
const structured::ClusterTree col_tree = nullptr 
)

Construct a StructuredMatrix using both a matrix-vector multiplication routine and a routine to extract a matrix sub-block.

Template Parameters
scalar_tprecision of input matrix, and of constructed StructuredMatrix. Note that not all types support every all precisions. See StructuredMatrix::Type.
Parameters
rowsNumber of rows of matrix to be constructed.
colsNumber of columns of matrix to be constructed.
AmultMatrix-(multi)vector multiplication routine.
AelemMatrix block extraction routine.
optsOptions object
row_treeoptional clustertree for the rows, see also strumpack::binary_tree_clustering
col_treeoptional clustertree for the columns. If the matrix is square, this does not need to be specified.
Returns
std::unique_ptr holding a pointer to a StructuredMatrix of the requested StructuredMatrix::Type
Exceptions
std::invalid_argumentIf the operatation is not supported for the type of structured::StructuredMatrix, if the type requires a square matrix and the input is not square, if the structured::StructuredMatrix type requires MPI.
std::logic_errorIf the operation is not implemented yet
std::runtime_errorIf the operation requires a third party library which was not enabled when configuring STRUMPACK.
See also
strumpack::binary_tree_clustering, construct_from_dense construct_from_elements, construct_matrix_free and construct_partially_matrix_free

◆ construct_partially_matrix_free() [3/3]

template<typename scalar_t >
std::unique_ptr< StructuredMatrix< scalar_t > > strumpack::structured::construct_partially_matrix_free ( int  rows,
int  cols,
const mult_t< scalar_t > &  Amult,
const extract_t< scalar_t > &  Aelem,
const StructuredOptions< scalar_t > &  opts,
const structured::ClusterTree row_tree = nullptr,
const structured::ClusterTree col_tree = nullptr 
)

Construct a StructuredMatrix using both a matrix-vector multiplication routine and a routine to extract individual matrix elements.

Template Parameters
scalar_tprecision of input matrix, and of constructed StructuredMatrix. Note that not all types support every all precisions. See StructuredMatrix::Type.
Parameters
rowsNumber of rows of matrix to be constructed.
colsNumber of columns of matrix to be constructed.
AmultMatrix-(multi)vector multiplication routine.
AelemMatrix element extraction routine.
optsOptions object
row_treeoptional clustertree for the rows, see also strumpack::binary_tree_clustering
col_treeoptional clustertree for the columns. If the matrix is square, this does not need to be specified.
Returns
std::unique_ptr holding a pointer to a StructuredMatrix of the requested StructuredMatrix::Type
Exceptions
std::invalid_argumentIf the operatation is not supported for the type of structured::StructuredMatrix, if the type requires a square matrix and the input is not square, if the structured::StructuredMatrix type requires MPI.
std::logic_errorIf the operation is not implemented yet
std::runtime_errorIf the operation requires a third party library which was not enabled when configuring STRUMPACK.
See also
strumpack::binary_tree_clustering, construct_from_dense construct_from_elements, construct_matrix_free and construct_partially_matrix_free