mptensor  v0.3.0
Parallel Library for Tensor Network Methods
Loading...
Searching...
No Matches
matrix_lapack_impl.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_LAPACK_IMPL_HPP_
29#define _MATRIX_LAPACK_IMPL_HPP_
30
31#include <algorithm>
32#include <cassert>
33#include <fstream>
34#include <iostream>
35#include <vector>
36#include <string>
37
38#include "../complex.hpp"
39
40namespace mptensor {
41
43namespace lapack {
44
46
51/* ---------- constructors ---------- */
52
54
57template <typename C>
59 init(0, 0);
60};
61
63
66template <typename C>
68 init(0, 0);
69};
70
72
76template <typename C>
77Matrix<C>::Matrix(size_t n_row, size_t n_col) {
78 init(n_row, n_col);
79};
80
82
87template <typename C>
88Matrix<C>::Matrix(const comm_type& comm_dummy, size_t n_row, size_t n_col) {
89 init(n_row, n_col);
90};
91
92/* ---------- member functions ---------- */
93template <typename C>
94inline const C& Matrix<C>::operator[](size_t i) const {
95 return V[i];
96}
97
98template <typename C>
99inline C& Matrix<C>::operator[](size_t i) {
100 return V[i];
101}
102
103template <typename C>
104inline const C* Matrix<C>::head() const {
105 return &(V[0]);
106}
107
108template <typename C>
110 return &(V[0]);
111}
112
113template <typename C>
114inline size_t Matrix<C>::local_size() const {
115 return V.size();
116}
117
119template <typename C>
120inline const typename Matrix<C>::comm_type& Matrix<C>::get_comm() const {
121 return comm_;
122}
123
125template <typename C>
126inline int Matrix<C>::get_comm_size() const {
127 return 1;
128}
129
131template <typename C>
132inline int Matrix<C>::get_comm_rank() const {
133 return 0;
134}
135
136template <typename C>
137inline int Matrix<C>::n_row() const {
138 return n_row_;
139}
140
141template <typename C>
142inline int Matrix<C>::n_col() const {
143 return n_col_;
144}
145
146template <typename C>
147inline void Matrix<C>::global_index(size_t i, size_t& g_row,
148 size_t& g_col) const {
149 g_row = i % n_row_;
150 g_col = i / n_row_;
151};
152
153template <typename C>
154inline bool Matrix<C>::local_index(size_t g_row, size_t g_col,
155 size_t& i) const {
156 i = g_row + g_col * n_row_; // column major
157 return true;
158};
159
161
167template <typename C>
168inline void Matrix<C>::local_position(size_t g_row, size_t g_col,
169 int& comm_rank, size_t& lindex) const {
170 comm_rank = 0;
171 lindex = g_row + g_col * n_row_;
172}
173
174template <typename C>
175inline size_t Matrix<C>::local_row_size() const {
176 return n_row_;
177}
178
179template <typename C>
180inline size_t Matrix<C>::local_col_size() const {
181 return n_col_;
182}
183
184template <typename C>
185inline size_t Matrix<C>::local_row_index(size_t lindex) const {
186 return lindex % n_row_;
187}
188
189template <typename C>
190inline size_t Matrix<C>::local_col_index(size_t lindex) const {
191 return lindex / n_row_;
192}
193
194template <typename C>
195inline size_t Matrix<C>::global_row_index(size_t lindex_row) const {
196 return lindex_row;
197}
198
199template <typename C>
200inline size_t Matrix<C>::global_col_index(size_t lindex_col) const {
201 return lindex_col;
202}
203
204template <typename C>
205inline void Matrix<C>::init(size_t n_row, size_t n_col) {
206 comm_ = 0;
207 n_row_ = n_row;
208 n_col_ = n_col;
209 V.resize(n_row * n_col);
210};
211
212template <typename C>
213inline void Matrix<C>::print_info(std::ostream& out) const {
214 out << "Matrix: local_size= " << local_size() << "\n";
215};
216
217template <typename C>
219 assert(V.size() == rhs.local_size());
220 for (size_t i = 0; i < V.size(); ++i) V[i] += rhs[i];
221 return *this;
222};
223
224template <typename C>
226 assert(V.size() == rhs.local_size());
227 for (size_t i = 0; i < V.size(); ++i) V[i] -= rhs[i];
228 return *this;
229};
230
231template <typename C>
233 for (size_t i = 0; i < V.size(); ++i) V[i] *= rhs;
234 return *this;
235};
236
237template <typename C>
239 for (size_t i = 0; i < V.size(); ++i) V[i] /= rhs;
240 return *this;
241};
242
243template <typename C>
244template <typename UnaryOperation>
246 std::transform(V.begin(), V.end(), V.begin(), op);
247 return *this;
248};
249
250template <typename C>
252 /* new matrix */
253 Matrix<C> M_new(get_comm(), n_col(), n_row());
254
255 size_t g_col, g_row, i_new;
256 for (size_t i = 0; i < local_size(); i++) {
257 global_index(i, g_row, g_col);
258 M_new.local_index(g_col, g_row, i_new);
259
260 M_new[i_new] = V[i];
261 }
262
263 return M_new;
264}
265
266template <typename C>
267inline std::vector<C> Matrix<C>::flatten() {
268 return V;
269};
270
271template <typename C>
272inline void Matrix<C>::barrier() const {
273 return;
274}
275
276template <typename C>
278 return value;
279}
280
281template <typename C>
282template <typename D>
283inline void Matrix<C>::bcast(D* buffer, int count, int root) const {
284 return;
285}
286
288template <typename C>
290 return;
291}
292
294template <typename C>
296 return;
297}
298
299template <typename C>
300void Matrix<C>::save_index(const std::string &filename) const {
301 std::ofstream fout(filename);
302 fout << "local_size= " << local_size() << "\n";
303 fout << "local_n_row= " << n_row_ << "\n";
304 fout << "local_n_col= " << n_col_ << "\n";
305 fout << "[global_row]\n";
306 for (size_t i = 0; i < n_row_; ++i) {
307 fout << " " << i;
308 if (i % 10 == 9) fout << "\n";
309 }
310 fout << "\n";
311 fout << "[global_col]\n";
312 for (size_t i = 0; i < n_col_; ++i) {
313 fout << " " << i;
314 if (i % 10 == 9) fout << "\n";
315 }
316 fout << "\n";
317 fout.close();
318}
319
320/* ---------- non-member functions ---------- */
321
322template <typename C>
323void replace_matrix_data(const Matrix<C>& M, const std::vector<int>& dest_rank,
324 const std::vector<size_t>& local_position,
325 Matrix<C>& M_new) {
326 assert(dest_rank.size() == M.local_size());
327 assert(local_position.size() == M.local_size());
328
329 const C* mat = M.head();
330 C* mat_new = M_new.head();
331
332 for (size_t i = 0; i < M.local_size(); i++) {
333 if (dest_rank[i] == 0) {
334 mat_new[local_position[i]] = mat[i];
335 }
336 }
337}
338
339template <typename C>
340void replace_matrix_data(const std::vector<C>& V, const std::vector<int>& dest_rank,
341 const std::vector<size_t>& local_position,
342 Matrix<C>& M_new) {
343 assert(dest_rank.size() == V.size());
344 assert(local_position.size() == V.size());
345
346 const C* mat = &(V[0]);
347 C* mat_new = M_new.head();
348
349 for (size_t i = 0; i < V.size(); i++) {
350 if (dest_rank[i] == 0) {
351 mat_new[local_position[i]] = mat[i];
352 }
353 }
354}
355
356template <typename C>
357void sum_matrix_data(const Matrix<C>& M, const std::vector<int>& dest_rank,
358 const std::vector<size_t>& local_position,
359 Matrix<C>& M_new) {
360 assert(dest_rank.size() == M.local_size());
361 assert(local_position.size() == M.local_size());
362
363 const C* mat = M.head();
364 C* mat_new = M_new.head();
365
366 for (size_t i = 0; i < M.local_size(); i++) {
367 if (dest_rank[i] == 0) {
368 mat_new[local_position[i]] += mat[i];
369 }
370 }
371}
372
373template <typename C>
375 const size_t n = A.local_size();
376 size_t g_row, g_col;
377 C val(0.0);
378
379 for (size_t i = 0; i < n; ++i) {
380 A.global_index(i, g_row, g_col);
381 if (g_row == g_col) val += A[i];
382 }
383
384 return val;
385};
386
387template <typename C>
388double max_abs(const Matrix<C>& a) {
389 const size_t n = a.local_size();
390 double val = 0.0;
391 for (size_t i = 0; i < n; ++i) {
392 val = std::max(val, std::abs(a[i]));
393 }
394 return val;
395};
396
397template <typename C>
398double min_abs(const Matrix<C>& a) {
399 const size_t n = a.local_size();
400 if (n == 0) return 0.0;
401 double val = std::abs(a[0]);
402 for (size_t i = 0; i < n; ++i) {
403 val = std::min(val, std::abs(a[i]));
404 }
405 return val;
406};
407
408} // namespace lapack
409} // namespace mptensor
410
411#endif // _MATRIX_LAPACK_IMPL_HPP_
void print_info(std::ostream &) const
Definition matrix_lapack_impl.hpp:213
size_t global_col_index(size_t lindex_col) const
Definition matrix_lapack_impl.hpp:200
Matrix & map(UnaryOperation op)
int n_row() const
Definition matrix_lapack_impl.hpp:137
int get_comm_size() const
Always returns 1.
Definition matrix_lapack_impl.hpp:126
Matrix & operator/=(C rhs)
Definition matrix_lapack_impl.hpp:238
size_t local_row_index(size_t lindex) const
Definition matrix_lapack_impl.hpp:185
Matrix & operator*=(C rhs)
Definition matrix_lapack_impl.hpp:232
void save_index(const std::string &filename) const
Definition matrix_lapack_impl.hpp:300
Matrix & operator+=(const Matrix &rhs)
Definition matrix_lapack_impl.hpp:218
std::vector< C > flatten()
Definition matrix_lapack_impl.hpp:267
size_t local_size() const
Definition matrix_lapack_impl.hpp:114
int get_comm_rank() const
Always returns 0.
Definition matrix_lapack_impl.hpp:132
const comm_type & get_comm() const
Always returns 0.
Definition matrix_lapack_impl.hpp:120
const Matrix transpose()
Definition matrix_lapack_impl.hpp:251
Matrix & operator-=(const Matrix &rhs)
Definition matrix_lapack_impl.hpp:225
Matrix()
Default constructor.
Definition matrix_lapack_impl.hpp:58
void prep_local_to_global() const
Do nothing.
Definition matrix_lapack_impl.hpp:289
size_t global_row_index(size_t lindex_row) const
Definition matrix_lapack_impl.hpp:195
void global_index(size_t i, size_t &g_row, size_t &g_col) const
Definition matrix_lapack_impl.hpp:147
const C & operator[](size_t i) const
Definition matrix_lapack_impl.hpp:94
void init(size_t n_row, size_t n_col)
Definition matrix_lapack_impl.hpp:205
size_t local_row_size() const
Definition matrix_lapack_impl.hpp:175
size_t local_col_index(size_t lindex) const
Definition matrix_lapack_impl.hpp:190
int comm_type
Definition matrix_lapack.hpp:47
size_t local_col_size() const
Definition matrix_lapack_impl.hpp:180
void prep_global_to_local() const
Do nothing.
Definition matrix_lapack_impl.hpp:295
void local_position(size_t g_row, size_t g_col, int &comm_rank, size_t &lindex) const
Convert a global index to an index of local storage.
Definition matrix_lapack_impl.hpp:168
void bcast(D *buffer, int count, int root) const
Definition matrix_lapack_impl.hpp:283
const C * head() const
Definition matrix_lapack_impl.hpp:104
C allreduce_sum(C value) const
Definition matrix_lapack_impl.hpp:277
bool local_index(size_t g_row, size_t g_col, size_t &i) const
Definition matrix_lapack_impl.hpp:154
void barrier() const
Definition matrix_lapack_impl.hpp:272
int n_col() const
Definition matrix_lapack_impl.hpp:142
std::string filename(const std::string &prefix, int proc_size)
Definition common.hpp:32
std::complex< double > complex
Definition complex.hpp:38
C matrix_trace(const Matrix< C > &a)
Definition matrix_lapack_impl.hpp:374
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_lapack_impl.hpp:323
double max_abs(const Matrix< C > &a)
Definition matrix_lapack_impl.hpp:388
double min_abs(const Matrix< C > &a)
Definition matrix_lapack_impl.hpp:398
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_lapack_impl.hpp:357
Definition complex.hpp:34