mptensor v0.4.0
Parallel Library for Tensor Network Methods
Loading...
Searching...
No Matches
matrix_scalapack.hpp
Go to the documentation of this file.
1/*
2 mptensor - Parallel Library for Tensor Network Methods
3
4 Copyright 2016 Satoshi Morita
5
6 mptensor is free software: you can redistribute it and/or modify it
7 under the terms of the GNU Lesser General Public License as
8 published by the Free Software Foundation, either version 3 of the
9 License, or (at your option) any later version.
10
11 mptensor is distributed in the hope that it will be useful, but
12 WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 Lesser General Public License for more details.
15
16 You should have received a copy of the GNU Lesser General Public
17 License along with mptensor. If not, see
18 <https://www.gnu.org/licenses/>.
19*/
20
28#ifndef _MATRIX_SCALAPACK_HPP_
29#define _MATRIX_SCALAPACK_HPP_
30#ifndef _NO_MPI
31
32#include <iostream>
33#include <vector>
34#include <string>
35
36#include <mpi.h>
37
38#include "blacsgrid.hpp"
39
40namespace mptensor {
41namespace scalapack {
42
45
46template <typename C>
47class Matrix {
48 public:
49 typedef C value_type;
51
52 constexpr static size_t matrix_type_tag = MATRIX_TYPE_TAG_SCALAPACK;
53 constexpr static char* matrix_type_name = (char*)"ScaLAPACK";
54
55 Matrix();
56 explicit Matrix(const MPI_Comm& comm);
57 Matrix(size_t n_row, size_t n_col);
58 Matrix(const MPI_Comm& comm, size_t n_row, size_t n_col);
59
60 void init(size_t n_row, size_t n_col);
61
62 const MPI_Comm& get_comm() const;
63 int get_comm_size() const;
64 int get_comm_rank() const;
65
66 void print_info(std::ostream&) const;
67
68 const C& operator[](size_t i) const;
69 C& operator[](size_t i);
70 const C* head() const;
71 C* head();
72
73 size_t local_size() const;
74 void global_index(size_t i, size_t& g_row, size_t& g_col) const;
75 bool local_index(size_t g_row, size_t g_col, size_t& i) const;
76 void local_position(size_t g_row, size_t g_col, int& comm_rank,
77 size_t& lindex) const;
78
79 size_t local_row_size() const;
80 size_t local_col_size() const;
81 size_t local_row_index(size_t lindex) const;
82 size_t local_col_index(size_t lindex) const;
83 size_t global_row_index(size_t lindex_row) const;
84 size_t global_col_index(size_t lindex_col) const;
85
86 Matrix& operator+=(const Matrix& rhs);
87 Matrix& operator-=(const Matrix& rhs);
90
91 template <typename UnaryOperation>
93
94 std::vector<C> flatten();
95
96 void barrier() const;
97 C allreduce_sum(C value) const;
98 template <typename D>
99 void bcast(D* buffer, int count, int root) const;
100
101 void prep_local_to_global() const;
102 void prep_global_to_local() const;
103
104 int n_row() const;
105 int n_col() const;
106 const int* descriptor() const;
107
108 const Matrix transpose();
109 void save_index(const std::string &filename) const;
110
111 private:
112 BlacsGrid grid;
113 std::vector<C> V; // local strage
114 std::vector<int> desc; // descriptor
115 static const size_t BLOCK_SIZE;
116 size_t local_row_size_;
117 size_t local_col_size_;
118
119 int get_lld() const;
120 // void set(C (*element)(size_t g_row, size_t g_col));
121 // const BlacsGrid &get_grid() const;
122 // int blacs_context() const;
123 // int nb_row() const;
124 // int nb_col() const;
125
126 mutable bool has_local_to_global;
127 mutable bool has_global_to_local;
128 mutable std::vector<size_t>
129 global_row;
130 mutable std::vector<size_t>
131 global_col;
132 mutable std::vector<size_t>
133 local_row;
134 mutable std::vector<size_t>
135 local_col;
136 mutable std::vector<int>
137 proc_row;
138 mutable std::vector<int>
139 proc_col;
140 mutable std::vector<int> lld_list;
142};
143
146template <typename C>
147void replace_matrix_data(const Matrix<C>& M, const std::vector<int>& dest_rank,
148 const std::vector<size_t>& local_position,
150template <typename C>
151void replace_matrix_data(const std::vector<C>& V, const std::vector<int>& dest_rank,
152 const std::vector<size_t>& local_position,
154template <typename C>
155void sum_matrix_data(const Matrix<C>& M, const std::vector<int>& dest_rank,
156 const std::vector<size_t>& local_position,
158
159template <typename C>
160double max_abs(const Matrix<C>& a);
161template <typename C>
162double min_abs(const Matrix<C>& a);
163template <typename C>
164C matrix_trace(const Matrix<C>& a);
165
166// defined in matrix_scalapack.cc
167template <typename C>
169template <typename C>
170int matrix_svd(Matrix<C>& a, Matrix<C>& u, std::vector<double>& s,
171 Matrix<C>& v);
172template <typename C>
173int matrix_svd(Matrix<C>& a, std::vector<double>& s);
174template <typename C>
176template <typename C>
177int matrix_eigh(Matrix<C>& a, std::vector<double>& s, Matrix<C>& u);
178template <typename C>
179int matrix_eigh(Matrix<C>& a, std::vector<double>& s);
180template <typename C>
181int matrix_eigh(Matrix<C>& a, Matrix<C>& b, std::vector<double>& s,
182 Matrix<C>& u);
183// template <typename C>
184// int matrix_eig(Matrix<C>& a, std::vector<complex>& s, Matrix<complex>& u);
185// template <typename C>
186// int matrix_eig(Matrix<C>& a, std::vector<complex>& s);
187template <typename C>
189template <typename C>
190double max(const Matrix<C>& a);
191template <typename C>
192double min(const Matrix<C>& a);
194
196
197} // namespace scalapack
198} // namespace mptensor
199
201
202#endif // _NO_MPI
203#endif // _MATRIX_SCALAPACK_HPP_
BlacsGrid class.
Definition blacsgrid.hpp:39
Distributed matrix using ScaLAPACK.
Definition matrix_scalapack.hpp:47
int n_row() const
Definition matrix_scalapack_impl.hpp:144
Matrix & operator*=(C rhs)
Definition matrix_scalapack_impl.hpp:372
const int * descriptor() const
Definition matrix_scalapack_impl.hpp:136
void bcast(D *buffer, int count, int root) const
Definition matrix_scalapack_impl.hpp:441
void print_info(std::ostream &) const
Definition matrix_scalapack_impl.hpp:351
void prep_local_to_global() const
Preprocess for fast conversion from local index to global one.
Definition matrix_scalapack_impl.hpp:240
C value_type
Definition matrix_scalapack.hpp:49
size_t local_row_index(size_t lindex) const
Definition matrix_scalapack_impl.hpp:217
bool local_index(size_t g_row, size_t g_col, size_t &i) const
Definition matrix_scalapack_impl.hpp:176
void global_index(size_t i, size_t &g_row, size_t &g_col) const
Definition matrix_scalapack_impl.hpp:165
MPI_Comm comm_type
Definition matrix_scalapack.hpp:50
int n_col() const
Definition matrix_scalapack_impl.hpp:149
size_t local_row_size() const
Definition matrix_scalapack_impl.hpp:207
Matrix & operator-=(const Matrix &rhs)
Definition matrix_scalapack_impl.hpp:365
const C * head() const
Definition matrix_scalapack_impl.hpp:103
Matrix()
Definition matrix_scalapack_impl.hpp:70
static constexpr size_t matrix_type_tag
Definition matrix_scalapack.hpp:52
void init(size_t n_row, size_t n_col)
Definition matrix_scalapack_impl.hpp:321
C allreduce_sum(C value) const
Definition matrix_scalapack_impl.hpp:435
Matrix & map(UnaryOperation op)
std::vector< C > flatten()
Definition matrix_scalapack_impl.hpp:417
const MPI_Comm & get_comm() const
Definition matrix_scalapack_impl.hpp:121
Matrix & operator/=(C rhs)
Definition matrix_scalapack_impl.hpp:378
void barrier() const
Definition matrix_scalapack_impl.hpp:430
void local_position(size_t g_row, size_t g_col, int &comm_rank, size_t &lindex) const
Definition matrix_scalapack_impl.hpp:197
const Matrix transpose()
Definition matrix_scalapack_impl.hpp:391
void save_index(const std::string &filename) const
Definition matrix_scalapack_impl.hpp:446
size_t local_col_index(size_t lindex) const
Definition matrix_scalapack_impl.hpp:222
size_t local_size() const
Definition matrix_scalapack_impl.hpp:113
size_t global_row_index(size_t lindex_row) const
Definition matrix_scalapack_impl.hpp:227
int get_comm_size() const
Definition matrix_scalapack_impl.hpp:126
int get_comm_rank() const
Definition matrix_scalapack_impl.hpp:131
void prep_global_to_local() const
Preprocess for fast conversion from global index to local one.
Definition matrix_scalapack_impl.hpp:271
Matrix & operator+=(const Matrix &rhs)
Definition matrix_scalapack_impl.hpp:358
const C & operator[](size_t i) const
Definition matrix_scalapack_impl.hpp:93
size_t local_col_size() const
Definition matrix_scalapack_impl.hpp:212
size_t global_col_index(size_t lindex_col) const
Definition matrix_scalapack_impl.hpp:233
static constexpr char * matrix_type_name
Definition matrix_scalapack.hpp:53
std::string filename(const std::string &prefix, int proc_size)
Definition common.hpp:32
std::complex< double > complex
Definition complex.hpp:38
constexpr size_t MATRIX_TYPE_TAG_SCALAPACK
Definition matrix.hpp:35
Implementation of scalapack::Matrix class.
int matrix_qr(Matrix< C > &a, Matrix< C > &r)
double min_abs(const Matrix< C > &a)
Definition matrix_scalapack_impl.hpp:716
void replace_matrix_data(const Matrix< C > &M, const std::vector< int > &dest_rank, const std::vector< size_t > &local_position, Matrix< C > &M_new)
Definition matrix_scalapack_impl.hpp:472
void sum_matrix_data(const Matrix< C > &M, const std::vector< int > &dest_rank, const std::vector< size_t > &local_position, Matrix< C > &M_new)
Definition matrix_scalapack_impl.hpp:638
double max(const Matrix< C > &a)
int matrix_svd(Matrix< C > &a, Matrix< C > &u, std::vector< double > &s, Matrix< C > &v)
double min(const Matrix< C > &a)
double max_abs(const Matrix< C > &a)
Definition matrix_scalapack_impl.hpp:704
C matrix_trace(const Matrix< C > &a)
Definition matrix_scalapack_impl.hpp:728
int matrix_solve(Matrix< C > &a, Matrix< C > &b)
int matrix_eigh(Matrix< C > &a, std::vector< double > &s, Matrix< C > &u)
void matrix_product(const Matrix< C > &a, const Matrix< C > &b, Matrix< C > &c)
Definition complex.hpp:34