Distributed operator

Builder

template<typename CoefficientPrecision, typename CoordinatePrecision = underlying_type<CoefficientPrecision>>
class DefaultApproximationBuilder

Public Functions

inline DefaultApproximationBuilder(const VirtualInternalGenerator<CoefficientPrecision> &generator, const Cluster<CoordinatePrecision> &target_cluster, const Cluster<CoordinatePrecision> &source_cluster, const HMatrixTreeBuilder<CoefficientPrecision, CoordinatePrecision> &hmatrix_tree_builder, MPI_Comm communicator)
inline DefaultApproximationBuilder(const VirtualGenerator<CoefficientPrecision> &generator, const Cluster<CoordinatePrecision> &target_cluster, const Cluster<CoordinatePrecision> &source_cluster, const HMatrixTreeBuilder<CoefficientPrecision, CoordinatePrecision> &hmatrix_tree_builder, MPI_Comm communicator)

Public Members

HMatrix<CoefficientPrecision, CoordinatePrecision> hmatrix
DistributedOperator<CoefficientPrecision> &distributed_operator
const HMatrix<CoefficientPrecision, CoordinatePrecision> *block_diagonal_hmatrix = {nullptr}
template<typename CoefficientPrecision, typename CoordinatePrecision = underlying_type<CoefficientPrecision>>
class DefaultLocalApproximationBuilder

Public Functions

inline DefaultLocalApproximationBuilder(const VirtualInternalGenerator<CoefficientPrecision> &generator, const Cluster<CoordinatePrecision> &target_cluster, const Cluster<CoordinatePrecision> &source_cluster, const HMatrixTreeBuilder<CoefficientPrecision, CoordinatePrecision> &hmatrix_tree_builder, MPI_Comm communicator)
inline DefaultLocalApproximationBuilder(const VirtualGenerator<CoefficientPrecision> &generator, const Cluster<CoordinatePrecision> &target_cluster, const Cluster<CoordinatePrecision> &source_cluster, const HMatrixTreeBuilder<CoefficientPrecision, CoordinatePrecision> &hmatrix_tree_builder, MPI_Comm communicator)

Public Members

HMatrix<CoefficientPrecision, CoordinatePrecision> hmatrix
DistributedOperator<CoefficientPrecision> &distributed_operator
const HMatrix<CoefficientPrecision, CoordinatePrecision> *block_diagonal_hmatrix = {nullptr}
template<typename CoefficientPrecision, typename CoordinatePrecision = underlying_type<CoefficientPrecision>>
class CustomApproximationBuilder

Public Functions

inline explicit CustomApproximationBuilder(const Cluster<CoordinatePrecision> &target_cluster, const Cluster<CoordinatePrecision> &source_cluster, char symmetry, char UPLO, MPI_Comm communicator, const VirtualLocalOperator<CoefficientPrecision> &local_operator)

Public Members

DistributedOperator<CoefficientPrecision> distributed_operator

DistributedOperator

template<typename CoefficientPrecision>
class DistributedOperator

Public Functions

DistributedOperator(const DistributedOperator&) = delete
DistributedOperator &operator=(const DistributedOperator&) = delete
DistributedOperator(DistributedOperator &&cluster) noexcept = default
DistributedOperator &operator=(DistributedOperator &&cluster) noexcept = default
virtual ~DistributedOperator() = default
inline explicit DistributedOperator(const VirtualPartition<CoefficientPrecision> &target_partition, const VirtualPartition<CoefficientPrecision> &source_partition, char symmetry, char UPLO, MPI_Comm comm)
inline void add_local_operator(const VirtualLocalOperator<CoefficientPrecision> *local_operator)
void vector_product_global_to_global(const CoefficientPrecision *const in, CoefficientPrecision *const out) const
void matrix_product_global_to_global(const CoefficientPrecision *const in, CoefficientPrecision *const out, int mu) const
void vector_product_transp_global_to_global(const CoefficientPrecision *const in, CoefficientPrecision *const out) const
void matrix_product_transp_global_to_global(const CoefficientPrecision *const in, CoefficientPrecision *const out, int mu) const
void vector_product_local_to_local(const CoefficientPrecision *const in, CoefficientPrecision *const out, CoefficientPrecision *work = nullptr) const
void matrix_product_local_to_local(const CoefficientPrecision *const in, CoefficientPrecision *const out, int mu, CoefficientPrecision *work = nullptr) const
void vector_product_transp_local_to_local(const CoefficientPrecision *const in, CoefficientPrecision *const out, CoefficientPrecision *work = nullptr) const
void matrix_product_transp_local_to_local(const CoefficientPrecision *const in, CoefficientPrecision *const out, int mu, CoefficientPrecision *work = nullptr) const
void internal_vector_product_global_to_local(const CoefficientPrecision *const in, CoefficientPrecision *const out) const
void internal_matrix_product_global_to_local(const CoefficientPrecision *const in, CoefficientPrecision *const out, int mu) const
void internal_vector_product_transp_local_to_global(const CoefficientPrecision *const in, CoefficientPrecision *const out) const
void internal_matrix_product_transp_local_to_global(const CoefficientPrecision *const in, CoefficientPrecision *const out, int mu) const
void internal_vector_product_local_to_local(const CoefficientPrecision *const in, CoefficientPrecision *const out, CoefficientPrecision *work = nullptr) const
void internal_matrix_product_local_to_local(const CoefficientPrecision *const in, CoefficientPrecision *const out, int mu, CoefficientPrecision *work = nullptr) const
void internal_vector_product_transp_local_to_local(const CoefficientPrecision *const in, CoefficientPrecision *const out, CoefficientPrecision *work = nullptr) const
void internal_matrix_product_transp_local_to_local(const CoefficientPrecision *const in, CoefficientPrecision *const out, int mu, CoefficientPrecision *work = nullptr) const
void internal_sub_matrix_product_to_local(const CoefficientPrecision *const in, CoefficientPrecision *const out, int mu, int offset, int size) const
void local_to_global_source(const CoefficientPrecision *const in, CoefficientPrecision *const out, const int &mu) const
void local_to_global_target(const CoefficientPrecision *const in, CoefficientPrecision *const out, const int &mu) const
inline bool &use_permutation()
inline const bool &use_permutation() const
inline char get_symmetry_type() const
inline char get_storage_type() const
inline MPI_Comm get_comm() const
inline const VirtualPartition<CoefficientPrecision> &get_target_partition() const
inline const VirtualPartition<CoefficientPrecision> &get_source_partition() const