Hierarchical matrix¶
Builder¶
-
template<typename CoefficientsPrecision, typename CoordinatesPrecision = htool::underlying_type<CoefficientsPrecision>>
class HMatrixBuilder¶ Public Functions
-
inline HMatrixBuilder(int target_number_of_points, int target_spatial_dimension, const CoordinatesPrecision *target_coordinates, const ClusterTreeBuilder<CoordinatesPrecision> *target_cluster_tree_builder, int source_number_of_points, int source_spatial_dimension, const CoordinatesPrecision *source_coordinates, const ClusterTreeBuilder<CoordinatesPrecision> *source_cluster_tree_builder)¶
-
inline HMatrixBuilder(int target_number_of_points, int target_spatial_dimension, const CoordinatesPrecision *target_coordinates, int source_number_of_points, int source_spatial_dimension, const CoordinatesPrecision *source_coordinates)¶
-
inline HMatrixBuilder(int number_of_points, int spatial_dimension, const CoordinatesPrecision *coordinates, const ClusterTreeBuilder<CoordinatesPrecision> *cluster_tree_builder = nullptr)¶
-
template<typename ExecutionPolicy>
inline HMatrix<CoefficientsPrecision, CoordinatesPrecision> build(ExecutionPolicy &&execution_policy, const VirtualGenerator<CoefficientsPrecision> &generator, const HMatrixTreeBuilder<CoefficientsPrecision, CoordinatesPrecision> &hmatrix_tree_builder)¶
-
inline HMatrix<CoefficientsPrecision, CoordinatesPrecision> build(const VirtualGenerator<CoefficientsPrecision> &generator, const HMatrixTreeBuilder<CoefficientsPrecision, CoordinatesPrecision> &hmatrix_tree_builder)¶
Public Members
-
Cluster<CoordinatesPrecision> target_cluster¶
-
std::unique_ptr<Cluster<CoordinatesPrecision>> source_cluster_ptr = nullptr¶
-
inline HMatrixBuilder(int target_number_of_points, int target_spatial_dimension, const CoordinatesPrecision *target_coordinates, const ClusterTreeBuilder<CoordinatesPrecision> *target_cluster_tree_builder, int source_number_of_points, int source_spatial_dimension, const CoordinatesPrecision *source_coordinates, const ClusterTreeBuilder<CoordinatesPrecision> *source_cluster_tree_builder)¶
-
template<typename CoefficientPrecision, typename CoordinatePrecision = underlying_type<CoefficientPrecision>>
class HMatrixTreeBuilder¶ Public Functions
-
HMatrixTreeBuilder(const HMatrixTreeBuilder&) = delete¶
-
HMatrixTreeBuilder &operator=(const HMatrixTreeBuilder&) = delete¶
-
HMatrixTreeBuilder(HMatrixTreeBuilder&&) noexcept = default¶
-
HMatrixTreeBuilder &operator=(HMatrixTreeBuilder&&) noexcept = default¶
-
virtual ~HMatrixTreeBuilder() = default¶
-
template<typename ExecutionPolicy>
HMatrixType build(ExecutionPolicy&&, const VirtualInternalGenerator<CoefficientPrecision> &generator, const ClusterType &target_root_cluster_tree, const ClusterType &source_root_cluster_tree, int target_partition_number = -1, int partition_number_for_symmetry = -1) const¶
-
template<typename ExecutionPolicy>
inline HMatrixType build(ExecutionPolicy &&execution_policy, const VirtualGenerator<CoefficientPrecision> &generator, const ClusterType &target_root_cluster_tree, const ClusterType &source_root_cluster_tree, int target_partition_number = -1, int partition_number_for_symmetry = -1) const¶
-
inline HMatrixType build(const VirtualInternalGenerator<CoefficientPrecision> &generator, const ClusterType &target_root_cluster_tree, const ClusterType &source_root_cluster_tree, int target_partition_number = -1, int partition_number_for_symmetry = -1) const¶
-
inline HMatrixType build(const VirtualGenerator<CoefficientPrecision> &generator, const ClusterType &target_root_cluster_tree, const ClusterType &source_root_cluster_tree, int target_partition_number = -1, int partition_number_for_symmetry = -1) const¶
-
HMatrixType sequential_build(const VirtualInternalGenerator<CoefficientPrecision> &generator, const ClusterType &target_root_cluster_tree, const ClusterType &source_root_cluster_tree, int target_partition_number = -1, int partition_number_for_symmetry = -1) const¶
-
inline HMatrixType sequential_build(const VirtualGenerator<CoefficientPrecision> &generator, const ClusterType &target_root_cluster_tree, const ClusterType &source_root_cluster_tree, int target_partition_number = -1, int partition_number_for_symmetry = -1) const¶
-
HMatrixType openmp_build(const VirtualInternalGenerator<CoefficientPrecision> &generator, const ClusterType &target_root_cluster_tree, const ClusterType &source_root_cluster_tree, int target_partition_number = -1, int partition_number_for_symmetry = -1) const¶
-
inline HMatrixType openmp_build(const VirtualGenerator<CoefficientPrecision> &generator, const ClusterType &target_root_cluster_tree, const ClusterType &source_root_cluster_tree, int target_partition_number = -1, int partition_number_for_symmetry = -1) const¶
-
HMatrixType task_based_build(const VirtualInternalGenerator<CoefficientPrecision> &generator, const ClusterType &target_root_cluster_tree, const ClusterType &source_root_cluster_tree, std::vector<HMatrixType*> &L0, int max_nb_nodes, int target_partition_number = -1, int partition_number_for_symmetry = -1) const¶
-
inline HMatrixType task_based_build(const VirtualGenerator<CoefficientPrecision> &generator, const ClusterType &target_root_cluster_tree, const ClusterType &source_root_cluster_tree, std::vector<HMatrixType*> &L0, int max_nb_nodes, int target_partition_number = -1, int partition_number_for_symmetry = -1) const¶
-
inline void set_symmetry(char symmetry_type)¶
-
inline void set_UPLO(char UPLO_type)¶
-
inline void set_minimal_source_depth(int minimal_source_depth)¶
-
inline void set_minimal_target_depth(int minimal_target_depth)¶
-
inline void set_block_tree_consistency(bool consistency)¶
-
inline int get_false_positive() const¶
-
inline char get_symmetry() const¶
-
inline char get_UPLO() const¶
-
inline double get_epsilon() const¶
-
inline double get_eta() const¶
-
inline std::shared_ptr<VirtualInternalLowRankGenerator<CoefficientPrecision>> get_internal_low_rank_generator() const¶
-
inline std::shared_ptr<VirtualLowRankGenerator<CoefficientPrecision>> get_low_rank_generator() const¶
-
template<typename ExecutionPolicy>
HMatrix<CoefficientPrecision, CoordinatePrecision> build(ExecutionPolicy &&execution_policy, const VirtualInternalGenerator<CoefficientPrecision> &generator, const ClusterType &root_target_cluster_tree, const ClusterType &root_source_cluster_tree, int target_partition_number, int partition_number_for_symmetry) const¶
-
HMatrixTreeBuilder(const HMatrixTreeBuilder&) = delete¶
HMatrix¶
-
template<typename CoefficientPrecision, typename CoordinatePrecision = underlying_type<CoefficientPrecision>>
class HMatrix : public htool::TreeNode<HMatrix<CoefficientPrecision, underlying_type<CoefficientPrecision>>, HMatrixTreeData<CoefficientPrecision, underlying_type<CoefficientPrecision>>>¶ Public Types
Public Functions
-
inline HMatrix(const Cluster<CoordinatePrecision> &target_cluster, const Cluster<CoordinatePrecision> &source_cluster)¶
-
inline HMatrix(const HMatrix &parent, const Cluster<CoordinatePrecision> *target_cluster, const Cluster<CoordinatePrecision> *source_cluster)¶
-
virtual ~HMatrix() = default¶
-
inline const Cluster<CoordinatePrecision> &get_target_cluster() const¶
-
inline const Cluster<CoordinatePrecision> &get_source_cluster() const¶
-
inline int nb_cols() const¶
-
inline int nb_rows() const¶
-
inline htool::underlying_type<CoefficientPrecision> get_epsilon() const¶
-
inline HMatrix<CoefficientPrecision, CoordinatePrecision> *get_child_or_this(const Cluster<CoordinatePrecision> &required_target_cluster, const Cluster<CoordinatePrecision> &required_source_cluster)¶
-
inline const HMatrix<CoefficientPrecision, CoordinatePrecision> *get_child_or_this(const Cluster<CoordinatePrecision> &required_target_cluster, const Cluster<CoordinatePrecision> &required_source_cluster) const¶
-
inline int get_rank() const¶
-
inline const Matrix<CoefficientPrecision> *get_dense_data() const¶
-
inline Matrix<CoefficientPrecision> *get_dense_data()¶
-
inline const LowRankMatrix<CoefficientPrecision> *get_low_rank_data() const¶
-
inline LowRankMatrix<CoefficientPrecision> *get_low_rank_data()¶
-
inline char get_symmetry() const¶
-
inline char get_UPLO() const¶
-
inline const HMatrixTreeData<CoefficientPrecision, CoordinatePrecision> *get_hmatrix_tree_data() const¶
-
inline const HMatrix<CoefficientPrecision> *get_sub_hmatrix(const Cluster<CoordinatePrecision> &target_cluster, const Cluster<CoordinatePrecision> &source_cluster) const¶
-
inline HMatrix<CoefficientPrecision> *get_sub_hmatrix(const Cluster<CoordinatePrecision> &target_cluster, const Cluster<CoordinatePrecision> &source_cluster)¶
-
inline StorageType get_storage_type() const¶
-
inline void set_symmetry(char symmetry)¶
-
inline void set_UPLO(char UPLO)¶
-
inline void set_symmetry_for_leaves(char symmetry)¶
-
inline void set_UPLO_for_leaves(char UPLO)¶
-
inline void set_target_cluster(const Cluster<CoordinatePrecision> *new_target_cluster)¶
-
inline bool is_dense() const¶
-
inline bool is_low_rank() const¶
-
inline bool is_hierarchical() const¶
-
inline void set_eta(CoordinatePrecision eta)¶
-
inline void set_epsilon(underlying_type<CoefficientPrecision> epsilon)¶
-
inline void set_minimal_target_depth(unsigned int minimal_target_depth)¶
-
inline void set_minimal_source_depth(unsigned int minimal_source_depth)¶
-
inline void set_block_tree_consistency(bool consistency)¶
-
inline char get_symmetry_for_leaves() const¶
-
inline char get_UPLO_for_leaves() const¶
-
inline bool is_block_tree_consistent() const¶
-
inline void compute_dense_data(const VirtualInternalGenerator<CoefficientPrecision> &generator)¶
-
inline bool compute_low_rank_data(const VirtualInternalLowRankGenerator<CoefficientPrecision> &low_rank_generator, int reqrank, underlying_type<CoefficientPrecision> epsilon)¶
-
inline void clear_low_rank_data()¶
-
inline void set_dense_data(std::unique_ptr<Matrix<CoefficientPrecision>> dense_matrix_ptr)¶
-
inline HMatrix(const Cluster<CoordinatePrecision> &target_cluster, const Cluster<CoordinatePrecision> &source_cluster)¶
Admissibility conditions¶
-
template<typename CoordinatePrecision>
class VirtualAdmissibilityCondition¶ Subclassed by htool::RjasanowSteinbach< CoordinatePrecision >
Public Functions
-
virtual bool ComputeAdmissibility(const Cluster<CoordinatePrecision> &target, const Cluster<CoordinatePrecision> &source, double eta) const = 0¶
-
inline virtual ~VirtualAdmissibilityCondition()¶
-
virtual bool ComputeAdmissibility(const Cluster<CoordinatePrecision> &target, const Cluster<CoordinatePrecision> &source, double eta) const = 0¶
-
template<typename CoordinatePrecision>
class RjasanowSteinbach : public htool::VirtualAdmissibilityCondition<CoordinatePrecision>¶ Public Functions
-
inline virtual bool ComputeAdmissibility(const Cluster<CoordinatePrecision> &target, const Cluster<CoordinatePrecision> &source, double eta) const override¶
-
inline virtual bool ComputeAdmissibility(const Cluster<CoordinatePrecision> &target, const Cluster<CoordinatePrecision> &source, double eta) const override¶
Low-rank compression¶
-
template<typename CoefficientPrecision, typename CoordinatesPrecision = underlying_type<CoefficientPrecision>>
class VirtualLowRankGenerator¶ Public Functions
-
inline VirtualLowRankGenerator()¶
-
virtual bool copy_low_rank_approximation(int M, int N, const int *rows, const int *cols, LowRankMatrix<CoefficientPrecision> &lrmat) const = 0¶
-
virtual bool copy_low_rank_approximation(int M, int N, const int *rows, const int *cols, int reqrank, LowRankMatrix<CoefficientPrecision> &lrmat) const = 0¶
-
inline virtual ~VirtualLowRankGenerator()¶
-
inline VirtualLowRankGenerator()¶
-
template<typename CoefficientPrecision>
class SVD : public htool::VirtualInternalLowRankGenerator<CoefficientPrecision>¶ Public Functions
-
inline SVD(const VirtualInternalGenerator<CoefficientPrecision> &A)¶
-
inline SVD(const VirtualGenerator<CoefficientPrecision> &A, const int *target_permutation, const int *source_permutation)¶
-
inline virtual bool copy_low_rank_approximation(int M, int N, int row_offset, int col_offset, LowRankMatrix<CoefficientPrecision> &lrmat) const override¶
-
inline virtual bool copy_low_rank_approximation(int M, int N, int row_offset, int col_offset, int rank, LowRankMatrix<CoefficientPrecision> &lrmat) const override¶
-
inline SVD(const VirtualInternalGenerator<CoefficientPrecision> &A)¶
-
template<typename CoefficientPrecision>
class fullACA : public htool::VirtualInternalLowRankGenerator<CoefficientPrecision>¶ Public Functions
-
inline fullACA(const VirtualInternalGenerator<CoefficientPrecision> &A)¶
-
inline fullACA(const VirtualGenerator<CoefficientPrecision> &A, const int *target_permutation, const int *source_permutation)¶
-
inline virtual bool copy_low_rank_approximation(int M, int N, int row_offset, int col_offset, LowRankMatrix<CoefficientPrecision> &lrmat) const override¶
-
inline virtual bool copy_low_rank_approximation(int M, int N, int row_offset, int col_offset, int reqrank, LowRankMatrix<CoefficientPrecision> &lrmat) const override¶
-
inline fullACA(const VirtualInternalGenerator<CoefficientPrecision> &A)¶
-
template<typename CoefficientPrecision>
class partialACA : public htool::VirtualInternalLowRankGenerator<CoefficientPrecision>¶ Public Functions
-
inline partialACA(const VirtualInternalGenerator<CoefficientPrecision> &A)¶
-
inline partialACA(const VirtualGenerator<CoefficientPrecision> &A, const int *target_permutation, const int *source_permutation)¶
-
inline virtual bool copy_low_rank_approximation(int M, int N, int row_offset, int col_offset, LowRankMatrix<CoefficientPrecision> &lrmat) const override¶
-
inline virtual bool copy_low_rank_approximation(int M, int N, int row_offset, int col_offset, int reqrank, LowRankMatrix<CoefficientPrecision> &lrmat) const override¶
-
inline partialACA(const VirtualInternalGenerator<CoefficientPrecision> &A)¶
-
template<typename CoefficientPrecision>
class sympartialACA : public htool::VirtualInternalLowRankGenerator<CoefficientPrecision>¶ Public Functions
-
inline sympartialACA(const VirtualInternalGenerator<CoefficientPrecision> &A)¶
-
inline sympartialACA(const VirtualGenerator<CoefficientPrecision> &A, const int *target_permutation, const int *source_permutation)¶
-
inline virtual bool copy_low_rank_approximation(int M, int N, int row_offset, int col_offset, LowRankMatrix<CoefficientPrecision> &lrmat) const override¶
-
inline virtual bool copy_low_rank_approximation(int M, int N, int row_offset, int col_offset, int reqrank, LowRankMatrix<CoefficientPrecision> &lrmat) const override¶
-
inline sympartialACA(const VirtualInternalGenerator<CoefficientPrecision> &A)¶
-
template<typename CoefficientPrecision>
class RecompressedLowRankGenerator : public htool::VirtualInternalLowRankGenerator<CoefficientPrecision>¶ Public Functions
-
inline RecompressedLowRankGenerator(const VirtualInternalLowRankGenerator<CoefficientPrecision> &low_rank_generator, std::function<void(LowRankMatrix<CoefficientPrecision>&)> recompression)¶
-
inline virtual bool copy_low_rank_approximation(int M, int N, int row_offset, int col_offset, LowRankMatrix<CoefficientPrecision> &lrmat) const¶
-
inline virtual bool copy_low_rank_approximation(int M, int N, int row_offset, int col_offset, int reqrank, LowRankMatrix<CoefficientPrecision> &lrmat) const¶
-
inline RecompressedLowRankGenerator(const VirtualInternalLowRankGenerator<CoefficientPrecision> &low_rank_generator, std::function<void(LowRankMatrix<CoefficientPrecision>&)> recompression)¶
Visualisation¶
-
template<typename CoefficientPrecision, typename CoordinatePrecision = underlying_type<CoefficientPrecision>>
void htool::save_leaves_with_rank(const HMatrix<CoefficientPrecision, CoordinatePrecision> &hmatrix, std::string filename)¶
-
template<typename CoefficientPrecision, typename CoordinatePrecision = underlying_type<CoefficientPrecision>>
std::map<std::string, std::string> htool::get_tree_parameters(const HMatrix<CoefficientPrecision, CoordinatePrecision> &hmatrix)¶
-
template<typename CoefficientPrecision, typename CoordinatePrecision = underlying_type<CoefficientPrecision>>
void htool::print_tree_parameters(const HMatrix<CoefficientPrecision, CoordinatePrecision> &hmatrix, std::ostream &os)¶
-
template<typename CoefficientPrecision, typename CoordinatePrecision = underlying_type<CoefficientPrecision>>
std::map<std::string, std::string> htool::get_hmatrix_information(const HMatrix<CoefficientPrecision, CoordinatePrecision> &hmatrix)¶
-
template<typename CoefficientPrecision, typename CoordinatePrecision = underlying_type<CoefficientPrecision>>
void htool::print_hmatrix_information(const HMatrix<CoefficientPrecision, CoordinatePrecision> &hmatrix, std::ostream &os)¶
Linear algebra¶
-
template<typename ExecutionPolicy, typename CoefficientPrecision, typename CoordinatePrecision = underlying_type<CoefficientPrecision>>
void htool::add_hmatrix_vector_product(ExecutionPolicy&&, char trans, CoefficientPrecision alpha, const HMatrix<CoefficientPrecision, CoordinatePrecision> &A, const CoefficientPrecision *x, CoefficientPrecision beta, CoefficientPrecision *y, CoefficientPrecision *buffer = nullptr)¶ It performs one of the \(\mathcal{H}\) matrix-vector operations.
where \(\alpha\) and \(\beta\) are scalars, x and y are vectors and A \(\mathcal{H}\) matrix of dimension \((m,n)\).\[\begin{split} \begin{array}{rl} y:=& \alpha*A*x+\beta*y, \\ y:=& \alpha*A^T*x+\beta*y, \\ y:=& \alpha*A^H*x+\beta*y, \end{array} \end{split}\]- Template Parameters:
ExecutionPolicy –
CoefficientPrecision –
CoordinatePrecision –
- Parameters:
execution_policy – [in]
trans – [in] is a character. It specifies the operation to be performed as follows:
\[\begin{split} \begin{array}{rl} \text{trans}=\text{'N'},\quad y:=& \alpha*A*x+\beta*y, \\ \text{trans}=\text{'T'},\quad y:=& \alpha*A^T*x+\beta*y, \\ \text{trans}=\text{'C'},\quad y:=& \alpha*A^H*x+\beta*y, \end{array} \end{split}\]alpha – [in] is a
Tvalue specifying the scalar \(\alpha\).A – [in] is a \(\mathcal{H}\) matrix of dimension \((m,n)\).
x – [in] is a
Tarray of dimension \(n\) whentransis ‘N’, and \(m\) otherwise.beta – [in] is a
Tvalue specifying the scalar \(\beta\).y – [inout] is a
Tarray of dimension \(m\) whentransis ‘N’, and \(n\) otherwise.buffer – [in] is a
Tarray. It defaults to nullptr, and it needs to be of dimension \(m+n\) to avoid internal allocation.
-
template<typename ExecutionPolicy, typename MatB, typename MatC, typename CoordinatePrecision = underlying_type<typename MatB::value_type>, typename = std::enable_if_t<std::is_same<typename MatB::value_type, typename MatC::value_type>::value>>
void htool::add_hmatrix_matrix_product(ExecutionPolicy &&execution_policy, char transa, char transb, typename MatB::value_type alpha, const HMatrix<typename MatB::value_type, CoordinatePrecision> &A, const MatB &B, typename MatB::value_type beta, MatC &C, typename MatB::value_type *buffer = nullptr)¶ It performs one of the \(\mathcal{H}\) matrix-matrix operations.
where \(\operatorname{op}(X)\) is one of\[ C:= \alpha*\operatorname{op}(A)*\operatorname{op}(B)+\beta*C \]\(\alpha\) and \(\beta\) are scalars, x and y are vectors, \(\operatorname{op}(A)\) is a \(m\) by \(k\) \(\mathcal{H}\)-matrix, \(\operatorname{op}(B)\) is a \(k\) by \(n\) matrix and C is a \(m\) by \(n\) matrix.\[ \operatorname{op}(X) = X, \quad \text{or} \quad \operatorname{op}(X) = X^T, \quad \text{or} \quad \operatorname{op}(X) = X^H, \]- Template Parameters:
ExecutionPolicy –
MatB –
MatC –
CoordinatePrecision –
- Parameters:
execution_policy – [in]
transa – [in] is a character. It specifies the form of \(\operatorname{op}(A)\) as follows:
\[\begin{split} \begin{array}{rlrl} \text{trans}=&\text{'N'},\quad &\operatorname{op}(A)=&A\\ \text{trans}=&\text{'T'},\quad &\operatorname{op}(A)=&A^T\\ \text{trans}=&\text{'C'},\quad &\operatorname{op}(A)=&A^H. \end{array} \end{split}\]transb – [in] is a character. It specifies the form of \(\operatorname{op}(B)\) as follows:
\[\begin{split} \begin{array}{rlrl} \text{trans}=&\text{'N'},\quad &\operatorname{op}(B)=&B\\ \text{trans}=&\text{'T'},\quad &\operatorname{op}(B)=&B^T\\ \text{trans}=&\text{'C'},\quad &\operatorname{op}(B)=&B^H. \end{array} \end{split}\]alpha – [in] is a
Tvalue specifying the scalar \(\alpha\).A – [in] is a \(\mathcal{H}\) matrix of dimension \((m,k)\) when
transais ‘N’, \((k,n)\) otherwise.B – [in] is a
Tarray of dimension \((k,n)\) whentransais ‘N’, \((n,k)\) otherwise.beta – [in] is a
Tvalue specifying the scalar \(\beta\).C – [inout] is a
Tarray of dimension \((m,n)\).buffer – [in] is a
Tarray. It defaults to nullptr, and it needs to be of dimension \((m+n)*k\) to avoid internal allocation.
-
template<typename ExecutionPolicy, typename CoefficientPrecision, typename CoordinatePrecision = underlying_type<CoefficientPrecision>>
void htool::lu_factorization(ExecutionPolicy&&, HMatrix<CoefficientPrecision, CoordinatePrecision> &hmatrix)¶
-
template<typename Mat, typename CoordinatePrecision = underlying_type<typename Mat::value_type>>
void htool::lu_solve(char trans, const HMatrix<typename Mat::value_type, CoordinatePrecision> &A, Mat &X)¶
-
template<typename ExecutionPolicy, typename CoefficientPrecision, typename CoordinatePrecision = underlying_type<CoefficientPrecision>>
void htool::cholesky_factorization(ExecutionPolicy&&, char UPLO, HMatrix<CoefficientPrecision, CoordinatePrecision> &hmatrix)¶