mptensor v0.4.0
Parallel Library for Tensor Network Methods
Loading...
Searching...
No Matches
tensor_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
29#ifndef _TENSOR_IMPL_HPP_
30#define _TENSOR_IMPL_HPP_
31
32#include <algorithm>
33#include <cassert>
34#include <cstdio>
35#include <cstdlib>
36#include <cstring>
37#include <fstream>
38#include <iostream>
39#include <vector>
40
41#include "complex.hpp"
42#include "index.hpp"
43#include "matrix.hpp"
44#include "tensor.hpp"
45
46namespace mptensor {
47
48/* Utilities */
49bool is_no_transpose(const Axes &axes, const Axes &axes_map, size_t rank);
50namespace debug {
51bool check_total_size(const Shape &s1, const Shape &s2);
52bool check_extend(const Shape &s_old, const Shape &s_new);
53bool check_transpose_axes(const Axes &axes, size_t rank);
54bool check_svd_axes(const Axes &axes_row, const Axes &axes_col, size_t rank);
55bool check_trace_axes(const Axes &axes_1, const Axes &axes_2, size_t rank);
56bool check_trace_axes(const Axes &axes_a, const Axes &axes_b,
57 const Shape &shape_a, const Shape &shape_b);
58bool check_contract_axes(const Axes &axes_1, const Axes &axes_2, size_t rank);
59bool check_square(const Shape &shape, size_t urank);
60} // namespace debug
61
62/* ---------- constructors ---------- */
69template <template <typename> class Matrix, typename C>
70Tensor<Matrix, C>::Tensor() : Mat(), upper_rank(0), axes_map(){};
71
78template <template <typename> class Matrix, typename C>
79Tensor<Matrix, C>::Tensor(const Shape &shape) : Mat() {
80 init(shape, shape.size() / 2);
81};
87template <template <typename> class Matrix, typename C>
89 : Mat(comm), upper_rank(0), axes_map(){};
90
97template <template <typename> class Matrix, typename C>
98Tensor<Matrix, C>::Tensor(const comm_type &comm, const Shape &shape)
99 : Mat(comm) {
100 init(shape, shape.size() / 2);
101};
102
104
109template <template <typename> class Matrix, typename C>
110Tensor<Matrix, C>::Tensor(const comm_type &comm, const Shape &shape,
111 const size_t upper_rank)
112 : Mat(comm) {
113 init(shape, upper_rank);
114};
115
122template <template <typename> class Matrix, typename C>
125 : Mat(comm) {
126 init(t.shape(), t.get_upper_rank());
127 const size_t n = Mat.local_size();
128 size_t idx;
129 int dummy;
130 for (size_t i = 0; i < n; ++i) {
131 t.local_position(global_index(i), dummy, idx);
132 Mat[i] = t[idx];
135
142template <template <typename> class Matrix, typename C>
143Tensor<Matrix, C>::Tensor(const comm_type &comm, const std::vector<C> &v)
144 : Mat(comm) {
145 init(Shape(v.size()), 0);
146 const size_t n = Mat.local_size();
147 Index idx;
148 for (size_t i = 0; i < n; ++i) {
149 idx = global_index(i);
150 Mat[i] = v[idx[0]];
151 }
152};
153
154/* ---------- inline member functions ---------- */
155
157template <template <typename> class Matrix, typename C>
158inline const Shape &Tensor<Matrix, C>::shape() const {
159 return Dim;
160};
161
163template <template <typename> class Matrix, typename C>
164inline size_t Tensor<Matrix, C>::rank() const {
165 return Dim.size();
166};
167
169
172template <template <typename> class Matrix, typename C>
173inline size_t Tensor<Matrix, C>::ndim() const {
174 return Dim.size();
175};
176
178template <template <typename> class Matrix, typename C>
179inline size_t Tensor<Matrix, C>::local_size() const {
180 return Mat.local_size();
181};
182
184template <template <typename> class Matrix, typename C>
186 return upper_rank;
187};
188
190
196template <template <typename> class Matrix, typename C>
198 return axes_map;
199};
200
202template <template <typename> class Matrix, typename C>
204 return Mat;
205};
206
208template <template <typename> class Matrix, typename C>
210 return Mat;
211};
212
214template <template <typename> class Matrix, typename C>
215inline const typename Tensor<Matrix, C>::comm_type &
217 return Mat.get_comm();
218};
219
221template <template <typename> class Matrix, typename C>
223 return Mat.get_comm_size();
224};
225
227template <template <typename> class Matrix, typename C>
229 return Mat.get_comm_rank();
230};
231
233
236template <template <typename> class Matrix, typename C>
237inline const C &Tensor<Matrix, C>::operator[](size_t local_idx) const {
238 return Mat[local_idx];
239}
240
242
245template <template <typename> class Matrix, typename C>
247 return Mat[local_idx];
248}
249
251template <template <typename> class Matrix, typename C>
255
257template <template <typename> class Matrix, typename C>
261
263template <template <typename> class Matrix, typename C>
265 const size_t rank = Dim.size();
266 const size_t rank0 = upper_rank;
267 const size_t rank1 = rank - rank0; // lower_rank
268 const size_t l_row = Mat.local_row_size();
269 const size_t l_col = Mat.local_col_size();
270 const size_t *axes_map0 = &(axes_map[0]);
271 const size_t *axes_map1 = &(axes_map[rank0]);
272
273 l2g_map_row.resize(l_row * rank0);
274 l2g_map_col.resize(l_col * rank1);
275
276#pragma omp parallel default(shared)
277 {
278 // Type int instead of size_t is used because of std::div()
279 int dim0[rank0];
280 int dim1[rank1];
281 for (size_t i = 0; i < rank0; ++i) dim0[i] = int(Dim[axes_map0[i]]);
282 for (size_t i = 0; i < rank1; ++i) dim1[i] = int(Dim[axes_map1[i]]);
283
284 int g_row, g_col; // global matrix row- and column-index
285 std::div_t divresult;
286 size_t *map;
287#pragma omp for
288 for (size_t k = 0; k < l_row; ++k) {
289 map = &(l2g_map_row[k * rank0]);
290 g_row = int(Mat.global_row_index(k));
291 for (size_t i = 0; i < rank0; ++i) {
292 divresult = std::div(g_row, dim0[i]);
293 map[i] = divresult.rem;
294 g_row = divresult.quot;
295 }
296 }
297#pragma omp for
298 for (size_t k = 0; k < l_col; ++k) {
299 map = &(l2g_map_col[k * rank1]);
300 g_col = int(Mat.global_col_index(k));
301 for (size_t i = 0; i < rank1; ++i) {
302 divresult = std::div(g_col, dim1[i]);
303 map[i] = divresult.rem;
304 g_col = divresult.quot;
305 }
306 }
307 }
308
309 return;
310};
311
313
319template <template <typename> class Matrix, typename C>
321 size_t gindex[]) const {
322 const size_t rank = Dim.size();
323 const size_t rank0 = upper_rank;
324 const size_t rank1 = rank - rank0; // lower_rank
325 const size_t lindex_row = Mat.local_row_index(lindex);
326 const size_t lindex_col = Mat.local_col_index(lindex);
327 const size_t *map_row = &(l2g_map_row[lindex_row * rank0]);
328 const size_t *map_col = &(l2g_map_col[lindex_col * rank1]);
329 const size_t *axes_map0 = &(axes_map[0]);
330 const size_t *axes_map1 = &(axes_map[rank0]);
331 for (size_t i = 0; i < rank0; ++i) {
333 }
334 for (size_t i = 0; i < rank1; ++i) {
336 }
337 return;
338};
339
341
352template <template <typename> class Matrix, typename C>
354 size_t lindex, const size_t axes_trans[], size_t index_new[]) const {
355 const size_t rank = Dim.size();
356 const size_t rank0 = upper_rank;
357 const size_t rank1 = rank - rank0; // lower_rank
358 const size_t lindex_row = Mat.local_row_index(lindex);
359 const size_t lindex_col = Mat.local_col_index(lindex);
360 const size_t *map_row = &(l2g_map_row[lindex_row * rank0]);
361 const size_t *map_col = &(l2g_map_col[lindex_col * rank1]);
362 // const size_t* axes_map0 = &(axes_map[0]);
363 // const size_t* axes_map1 = &(axes_map[rank0]);
364 for (size_t i = 0; i < rank0; ++i) {
366 // index_new[axes_inv[axes_map0[i]]] = map_row[i];
367 }
368 for (size_t i = 0; i < rank1; ++i) {
370 // index_new[axes_inv[axes_map1[i]]] = map_col[i];
371 }
372 return;
373};
374
375template <template <typename> class Matrix, typename C>
377 int &comm_rank,
378 size_t &local_idx) const {
380}
381
382/* ---------- member functions ---------- */
383
385
389template <template <typename> class Matrix, typename C>
390inline void Tensor<Matrix, C>::init(const Shape &shape, size_t urank) {
391 init(shape, urank, range(shape.size()));
392}
393
395
400template <template <typename> class Matrix, typename C>
401void Tensor<Matrix, C>::init(const Shape &shape, size_t urank,
402 const Axes &map) {
403 Dim = shape;
404 axes_map = map;
405 const size_t rank = Dim.size();
406 assert(urank <= rank);
407 upper_rank = urank;
408
409 size_t n_row = 1;
410 size_t n_col = 1;
411 for (size_t i = 0; i < upper_rank; ++i) n_row *= shape[map[i]];
412 for (size_t i = upper_rank; i < rank; ++i) n_col *= shape[map[i]];
413 Mat.init(n_row, n_col);
414}
415
417
423template <template <typename> class Matrix, typename C>
425 const std::string &tag) const {
426 if (tag.size() > 0) {
427 out << tag << ":\t";
428 }
429 out << "Tensor: shape= " << Dim << " upper_rank= " << upper_rank
430 << " axes_map= " << axes_map << "\t";
431 Mat.print_info(out);
432};
433
435
441template <template <typename> class Matrix, typename C>
443 const std::string &tag) const {
444 const int mpisize = get_comm_size();
445 const int mpirank = get_comm_rank();
446 const typename Tensor<Matrix, C>::comm_type &comm = get_comm();
447 Mat.barrier();
448 for (int i = 0; i < mpisize; ++i) {
449 if (i == mpirank) {
450 if (tag.size() > 0) {
451 out << tag << ":\t";
452 }
453 out << "mpirank: " << mpirank << "\t";
454 print_info(out, "");
455 }
456 Mat.barrier();
457 }
458}
459
461
469template <template <typename> class Matrix, typename C>
470bool Tensor<Matrix, C>::local_index(const Index &gindex, size_t &lindex) const {
471 const size_t rank = gindex.size();
472 assert(rank == Dim.size());
473 assert(rank > 0);
474 size_t g_row(0), g_col(0);
475 size_t d_row(1), d_col(1);
476 for (size_t i = 0; i < upper_rank; ++i) {
477 size_t j = axes_map[i];
478 g_row += gindex[j] * d_row;
479 d_row *= Dim[j];
480 }
481 for (size_t i = upper_rank; i < rank; ++i) {
482 size_t j = axes_map[i];
483 g_col += gindex[j] * d_col;
484 d_col *= Dim[j];
485 }
486 return Mat.local_index(g_row, g_col, lindex);
487};
488
490
494template <template <typename> class Matrix, typename C>
496 const size_t rank = Dim.size();
498 size_t g_row, g_col;
499 Mat.global_index(lindex, g_row, g_col);
500 gindex.resize(rank);
501 std::div_t divresult;
502 for (size_t i = 0; i < upper_rank; ++i) {
503 size_t j = axes_map[i];
504 divresult = std::div(int(g_row), int(Dim[j]));
505 gindex[j] = divresult.rem;
506 g_row = divresult.quot;
507 }
508 for (size_t i = upper_rank; i < rank; ++i) {
509 size_t j = axes_map[i];
510 divresult = std::div(int(g_col), int(Dim[j]));
511 gindex[j] = divresult.rem;
512 g_col = divresult.quot;
513 }
514 return gindex;
515};
516
518
524template <template <typename> class Matrix, typename C>
526 const size_t rank = Dim.size();
527 size_t g_row, g_col;
528 Mat.global_index(lindex, g_row, g_col);
529 std::div_t divresult;
530 for (size_t i = 0; i < upper_rank; ++i) {
531 size_t j = axes_map[i];
532 divresult = std::div(int(g_row), int(Dim[j]));
533 gindex[j] = divresult.rem;
534 g_row = divresult.quot;
535 }
536 for (size_t i = upper_rank; i < rank; ++i) {
537 size_t j = axes_map[i];
538 divresult = std::div(int(g_col), int(Dim[j]));
539 gindex[j] = divresult.rem;
540 g_col = divresult.quot;
541 }
542 return;
543};
544
546
551template <template <typename> class Matrix, typename C>
552bool Tensor<Matrix, C>::get_value(const Index &idx, C &val) const {
553 size_t li;
554 if (local_index(idx, li)) {
555 val = Mat[li];
556 return true;
557 } else {
558 return false;
559 }
560}
561
563
569template <template <typename> class Matrix, typename C>
571 size_t li;
572 if (local_index(idx, li)) {
573 Mat[li] = val;
574 }
575}
576
578
583template <template <typename> class Matrix, typename C>
585 size_t &local_idx) const {
586 const size_t rank = Dim.size();
587 size_t g_row(0), g_col(0);
588 size_t d_row(1), d_col(1);
589 for (size_t i = 0; i < upper_rank; ++i) {
590 size_t j = axes_map[i];
591 g_row += index[j] * d_row;
592 d_row *= Dim[j];
593 }
594 for (size_t i = upper_rank; i < rank; ++i) {
595 size_t j = axes_map[i];
596 g_col += index[j] * d_col;
597 d_col *= Dim[j];
598 }
599 Mat.local_position(g_row, g_col, comm_rank, local_idx);
600 return;
601}
602
604
610template <template <typename> class Matrix, typename C>
612 const Axes &new_axes_map) {
613 if ((upper_rank == new_upper_rank) && (axes_map == new_axes_map)) return;
614
616
617 const size_t rank = this->rank();
618 Shape dim;
619 Axes axes;
620 dim.resize(rank);
621 axes.resize(rank);
622 for (size_t i = 0; i < rank; ++i) {
623 dim[i] = Dim[new_axes_map[i]];
624 axes[new_axes_map[i]] = i;
625 }
626 init(dim, new_upper_rank);
627 transpose(axes);
628
629 /* create lists of local position and destination rank */
630 const size_t local_size = T_old.local_size();
631 std::vector<int> dest_mpirank(local_size);
632 std::vector<size_t> local_position(local_size);
633
634 T_old.prep_local_to_global();
635 this->prep_global_to_local();
636
637#pragma omp parallel default(shared)
638 {
639 Index index;
640 index.resize(rank);
641
642#pragma omp for
643 for (size_t i = 0; i < local_size; ++i) {
644 T_old.global_index_fast(i, index);
645
646 int dest;
647 size_t pos;
648 this->local_position(index, dest, pos);
649
650 local_position[i] = pos;
652 }
653 }
654
655 /* exchange data */
656 replace_matrix_data(T_old.get_matrix(), dest_mpirank, local_position, Mat);
657
658 return;
659}
660
662
666template <template <typename> class Matrix, typename C>
668 const size_t rank = Dim.size();
669 assert(debug::check_transpose_axes(axes, rank));
670
671 Shape dim_now = Dim;
672 Axes map_now = axes_map;
674 axes_inv.resize(rank);
675 for (size_t i = 0; i < rank; ++i) {
676 axes_inv[axes[i]] = i;
677 }
678
679 for (size_t i = 0; i < rank; ++i) {
680 Dim[i] = dim_now[axes[i]];
681 axes_map[i] = axes_inv[map_now[i]];
682 }
683
684 return (*this);
685}
686
688
698template <template <typename> class Matrix, typename C>
699template <typename D>
701 size_t n_axes) {
702 assert(Dim[n_axes] <= vec.size());
703 const size_t local_size = this->local_size();
704 prep_local_to_global();
705#pragma omp parallel default(shared)
706 {
707 Index idx;
708 idx.resize(rank());
709#pragma omp for
710 for (size_t i = 0; i < local_size; ++i) {
711 global_index_fast(i, idx);
712 Mat[i] *= vec[idx[n_axes]];
713 }
714 }
715 return (*this);
716};
717
719
731template <template <typename> class Matrix, typename C>
732template <typename D0, typename D1>
734 const std::vector<D0> &vec0, size_t n_axes0, const std::vector<D1> &vec1,
735 size_t n_axes1) {
736 assert(Dim[n_axes0] <= vec0.size());
737 assert(Dim[n_axes1] <= vec1.size());
738 const size_t local_size = this->local_size();
739 prep_local_to_global();
740#pragma omp parallel default(shared)
741 {
742 Index idx;
743 idx.resize(rank());
744#pragma omp for
745 for (size_t i = 0; i < local_size; ++i) {
746 global_index_fast(i, idx);
747 Mat[i] *= vec0[idx[n_axes0]] * vec1[idx[n_axes1]];
748 }
749 }
750 return (*this);
751};
752
754
768template <template <typename> class Matrix, typename C>
769template <typename D0, typename D1, typename D2>
771 const std::vector<D0> &vec0, size_t n_axes0, const std::vector<D1> &vec1,
772 size_t n_axes1, const std::vector<D2> &vec2, size_t n_axes2) {
773 assert(Dim[n_axes0] <= vec0.size());
774 assert(Dim[n_axes1] <= vec1.size());
775 assert(Dim[n_axes2] <= vec2.size());
776 const size_t local_size = this->local_size();
777 prep_local_to_global();
778#pragma omp parallel default(shared)
779 {
780 Index idx;
781 idx.resize(rank());
782#pragma omp for
783 for (size_t i = 0; i < local_size; ++i) {
784 global_index_fast(i, idx);
785 Mat[i] *= vec0[idx[n_axes0]] * vec1[idx[n_axes1]] * vec2[idx[n_axes2]];
786 }
787 }
788 return (*this);
789};
790
792
808template <template <typename> class Matrix, typename C>
809template <typename D0, typename D1, typename D2, typename D3>
811 const std::vector<D0> &vec0, size_t n_axes0, const std::vector<D1> &vec1,
812 size_t n_axes1, const std::vector<D2> &vec2, size_t n_axes2,
813 const std::vector<D3> &vec3, size_t n_axes3) {
814 assert(Dim[n_axes0] <= vec0.size());
815 assert(Dim[n_axes1] <= vec1.size());
816 assert(Dim[n_axes2] <= vec2.size());
817 assert(Dim[n_axes3] <= vec3.size());
818 const size_t local_size = this->local_size();
819 prep_local_to_global();
820#pragma omp parallel default(shared)
821 {
822 Index idx;
823 idx.resize(rank());
824#pragma omp for
825 for (size_t i = 0; i < local_size; ++i) {
826 global_index_fast(i, idx);
827 Mat[i] *= vec0[idx[n_axes0]] * vec1[idx[n_axes1]] * vec2[idx[n_axes2]] *
828 vec3[idx[n_axes3]];
829 }
830 }
831 return (*this);
832};
833
835
845template <template <typename> class Matrix, typename C>
847 const size_t n_axes,
848 const size_t i_begin,
849 const size_t i_end) {
850 assert(rank() == a.rank());
851 assert(n_axes < rank());
853 assert(i_end <= shape()[n_axes]);
854 assert(i_end - i_begin == a.shape()[n_axes]);
855
856 /* create lists of local position and destination rank */
857 const size_t local_size = a.local_size();
858 std::vector<int> dest_mpirank(local_size);
859 std::vector<size_t> local_position(local_size);
860
862 prep_global_to_local();
863
864#pragma omp parallel default(shared)
865 {
866 Index index;
867 index.resize(rank());
868#pragma omp for
869 for (size_t i = 0; i < local_size; ++i) {
870 a.global_index_fast(i, index);
871 index[n_axes] += i_begin;
872
873 int dest;
874 size_t pos;
875 this->local_position(index, dest, pos);
876
877 local_position[i] = pos;
879 }
880 }
881
882 /* exchange data */
883 replace_matrix_data(a.get_matrix(), dest_mpirank, local_position,
884 get_matrix());
885
886 return (*this);
887}
888
890
900template <template <typename> class Matrix, typename C>
902 const Index &index_begin,
903 const Index &index_end) {
904 const size_t nr = rank();
905 assert(nr == a.rank());
906 assert(index_begin.size() == nr);
907 assert(index_end.size() == nr);
908 for (size_t r = 0; r < nr; ++r) {
909 assert(index_end[r] <= shape()[r]);
911 }
912
913 /* create lists of local position and destination rank */
914 const size_t local_size = a.local_size();
915 std::vector<int> dest_mpirank(local_size);
916 std::vector<size_t> local_position(local_size);
917
919 prep_global_to_local();
920
921#pragma omp parallel default(shared)
922 {
923 Index index;
924 index.resize(nr);
925#pragma omp for
926 for (size_t i = 0; i < local_size; ++i) {
927 a.global_index_fast(i, index);
928 for (size_t r = 0; r < nr; ++r) {
929 if (index_begin[r] != index_end[r]) index[r] += index_begin[r];
930 }
931
932 int dest;
933 size_t pos;
934 this->local_position(index, dest, pos);
935
936 local_position[i] = pos;
938 }
939 }
940
941 /* exchange data */
942 replace_matrix_data(a.get_matrix(), dest_mpirank, local_position,
943 get_matrix());
944
945 return (*this);
946}
947
948
950
954template <template <typename> class Matrix, typename C>
956 const size_t n = rank();
957 if (!(axes_map == range(n))) {
958 change_configuration(n / 2, range(n));
959 }
960 Tensor<lapack::Matrix, C> T(get_comm(), get_matrix().flatten());
961 return reshape(T, Dim);
962}
963
964
966
970template <template <typename> class Matrix, typename C>
972 const size_t n = rank();
973 if (!(axes_map == range(n))) {
974 change_configuration(n / 2, range(n));
975 }
976 return get_matrix().flatten();
977}
978
980template <template <typename> class Matrix, typename C>
982 assert(Dim == rhs.shape());
983 change_configuration(rhs.get_upper_rank(), rhs.get_axes_map());
984 Mat += rhs.get_matrix();
985 return (*this);
986}
987
989template <template <typename> class Matrix, typename C>
991 assert(Dim == rhs.shape());
992 change_configuration(rhs.get_upper_rank(), rhs.get_axes_map());
993 Mat -= rhs.get_matrix();
994 return (*this);
995}
996
998template <template <typename> class Matrix, typename C>
1000 Mat *= rhs;
1001 return (*this);
1002}
1003
1005template <template <typename> class Matrix, typename C>
1007 Mat /= rhs;
1008 return (*this);
1009}
1010
1012template <template <typename> class Matrix, typename C>
1014 const size_t n = Mat.local_size();
1015 for (int i = 0; i < n; ++i) Mat[i] = rhs;
1016 return (*this);
1017}
1018
1020
1024template <template <typename> class Matrix, typename C>
1025template <typename UnaryOperation>
1027 Mat.map(op);
1028 return (*this);
1029}
1030
1031/* ---------- non-member functions (Opeations) ---------- */
1032
1034
1040template <template <typename> class Matrix, typename C>
1042 return T.transpose(axes);
1043}
1044
1046
1055template <template <typename> class Matrix, typename C>
1057 size_t urank_new) {
1058 const size_t rank = T.rank();
1059 assert(debug::check_transpose_axes(axes, rank));
1060 assert(urank_new <= rank);
1061
1062 if (urank_new == T.get_upper_rank()) {
1063 if (is_no_transpose(axes, T.get_axes_map(), rank)) return T;
1064 }
1065
1066 /* new index dimension */
1067 Shape dim_old = T.shape();
1068 Shape dim_new;
1069 dim_new.resize(rank);
1070 for (size_t r = 0; r < rank; ++r) {
1071 dim_new[r] = dim_old[axes[r]];
1072 }
1073
1074 /* new tensor */
1076
1077 /* create lists of local position and destination rank */
1078 const size_t local_size = T.local_size();
1079 std::vector<int> dest_mpirank(local_size);
1080 std::vector<size_t> local_position(local_size);
1081
1082 T.make_l2g_map();
1083 T_new.prep_global_to_local();
1084
1085#pragma omp parallel default(shared)
1086 {
1087 // Assuming axes_map of T_new is [0,1,2,...]
1088 size_t d_offset[rank];
1089 size_t d_row(1), d_col(1);
1090 for (size_t r = 0; r < urank_new; ++r) {
1091 d_offset[r] = d_row;
1092 d_row *= dim_new[r];
1093 }
1094 for (size_t r = urank_new; r < rank; ++r) {
1095 d_offset[r] = d_col;
1096 d_col *= dim_new[r];
1097 }
1098
1099 Axes axes_map = T.get_axes_map();
1100 size_t axes_inv[rank];
1101 size_t axes_trans[rank];
1102 for (size_t r = 0; r < rank; ++r) {
1103 axes_inv[axes[r]] = r;
1104 }
1105 for (size_t r = 0; r < rank; ++r) {
1106 axes_trans[r] = axes_inv[axes_map[r]];
1107 }
1108
1109 size_t index_new[rank];
1110
1111#pragma omp for
1112 for (size_t i = 0; i < local_size; ++i) {
1113 T.global_index_l2g_map_transpose(i, axes_trans, index_new);
1114
1115 size_t g_row(0),
1116 g_col(0); // indices of global matrix for transposed tensor
1117 for (size_t r = 0; r < urank_new; ++r) {
1118 g_row += index_new[r] * d_offset[r];
1119 }
1120 for (size_t r = urank_new; r < rank; ++r) {
1121 g_col += index_new[r] * d_offset[r];
1122 }
1123
1124 T_new.local_position_fast(g_row, g_col, dest_mpirank[i],
1125 local_position[i]);
1126 }
1127 }
1128
1129 /* exchange data */
1130 replace_matrix_data(T.get_matrix(), dest_mpirank, local_position,
1131 T_new.get_matrix());
1132
1133 return T_new;
1134};
1135
1137
1148template <template <typename> class Matrix, typename C>
1150 assert(debug::check_total_size(shape_new, T.shape()));
1151
1152 const Shape &shape = T.shape();
1153 const size_t rank = shape.size();
1154 const size_t rank_new = shape_new.size();
1155
1156 /* initialize new tensor */
1157 Tensor<Matrix, C> T_new(T.get_comm(), shape_new);
1158
1159 /* create lists of local position and destination rank */
1160 const size_t local_size = T.local_size();
1161 std::vector<int> dest_mpirank(local_size);
1162 std::vector<size_t> local_position(local_size);
1163
1164 T.prep_local_to_global();
1165 T_new.prep_global_to_local();
1166
1167#pragma omp parallel default(shared)
1168 {
1169 Index index, index_new;
1170 index.resize(rank);
1171 index_new.resize(rank_new);
1172#pragma omp for
1173 for (size_t i = 0; i < local_size; ++i) {
1174 T.global_index_fast(i, index);
1175 size_t global_position(0);
1176 size_t dim(1);
1177 for (size_t r = 0; r < rank; ++r) {
1178 global_position += index[r] * dim;
1179 dim *= shape[r];
1180 }
1181
1182 for (size_t r = 0; r < rank_new; ++r) {
1185 }
1186
1187 int dest;
1188 size_t pos;
1189 T_new.local_position(index_new, dest, pos);
1190
1191 local_position[i] = pos;
1192 dest_mpirank[i] = dest;
1193 }
1194 }
1195
1196 /* exchange data */
1197 replace_matrix_data(T.get_matrix(), dest_mpirank, local_position,
1198 T_new.get_matrix());
1199
1200 return T_new;
1201};
1202
1204
1214template <template <typename> class Matrix, typename C>
1216 size_t i_begin, size_t i_end) {
1217 const int mpisize = T.get_comm_size();
1218 const Shape &shape = T.shape();
1219 assert(n_axes < T.rank());
1220 assert(i_begin < i_end);
1221 assert(i_end <= shape[n_axes]);
1222
1223 Shape shape_new = shape;
1225
1226 /* initialize new tensor */
1227 Tensor<Matrix, C> T_new(T.get_comm(), shape_new);
1228
1229 /* create lists of local position and destination rank */
1230 const size_t local_size = T.local_size();
1231 std::vector<int> dest_mpirank(local_size);
1232 std::vector<size_t> local_position(local_size);
1233
1234 T.prep_local_to_global();
1235 T_new.prep_global_to_local();
1236
1237#pragma omp parallel default(shared)
1238 {
1239 Index index;
1240 index.resize(T.rank());
1241#pragma omp for
1242 for (size_t i = 0; i < local_size; ++i) {
1243 T.global_index_fast(i, index);
1244 if (index[n_axes] >= i_begin && index[n_axes] < i_end) {
1245 index[n_axes] -= i_begin;
1246
1247 int dest;
1248 size_t pos;
1249 T_new.local_position(index, dest, pos);
1250
1251 local_position[i] = pos;
1252 dest_mpirank[i] = dest;
1253 } else {
1254 // not send
1255 local_position[i] = 0;
1256 dest_mpirank[i] = mpisize;
1257 }
1258 }
1259 }
1260
1261 /* exchange data */
1262 replace_matrix_data(T.get_matrix(), dest_mpirank, local_position,
1263 T_new.get_matrix());
1264
1265 return T_new;
1266}
1267
1269
1282template <template <typename> class Matrix, typename C>
1284 const Index &index_end) {
1285 const int mpisize = T.get_comm_size();
1286 const Shape &shape = T.shape();
1287 const size_t rank = T.rank();
1288 assert(index_begin.size() == rank);
1289 assert(index_end.size() == rank);
1290
1291 Shape shape_new = shape;
1292 for (size_t r = 0; r < rank; ++r) {
1293 assert(index_end[r] <= shape[r]);
1296 if (shape_new[r] == 0) shape_new[r] = shape[r];
1297 }
1298
1299 /* initialize new tensor */
1300 Tensor<Matrix, C> T_new(T.get_comm(), shape_new);
1301
1302 /* create lists of local position and destination rank */
1303 const size_t local_size = T.local_size();
1304 std::vector<int> dest_mpirank(local_size);
1305 std::vector<size_t> local_position(local_size);
1306
1307 T.prep_local_to_global();
1308 T_new.prep_global_to_local();
1309
1310#pragma omp parallel default(shared)
1311 {
1312 Index index;
1313 index.resize(T.rank());
1314#pragma omp for
1315 for (size_t i = 0; i < local_size; ++i) {
1316 T.global_index_fast(i, index);
1317
1318 bool is_send = true;
1319 for (size_t r = 0; r < rank; ++r) {
1320 size_t idx = index[r];
1321 if (index_begin[r] == index_end[r]) {
1322 continue;
1323 } else if (idx >= index_begin[r] && idx < index_end[r]) {
1324 index[r] -= index_begin[r];
1325 } else {
1326 is_send = false;
1327 break;
1328 }
1329 }
1330
1331 if (is_send) {
1332 int dest;
1333 size_t pos;
1334 T_new.local_position(index, dest, pos);
1335
1336 local_position[i] = pos;
1337 dest_mpirank[i] = dest;
1338 } else {
1339 // not send
1340 local_position[i] = 0;
1341 dest_mpirank[i] = mpisize;
1342 }
1343 }
1344 }
1345
1346 /* exchange data */
1347 replace_matrix_data(T.get_matrix(), dest_mpirank, local_position,
1348 T_new.get_matrix());
1349
1350 return T_new;
1351};
1352
1354
1362template <template <typename> class Matrix, typename C>
1364 assert(T.rank() == shape_new.size());
1365 assert(debug::check_extend(T.shape(), shape_new));
1366
1367 /* initialize new tensor */
1368 Tensor<Matrix, C> T_new(T.get_comm(), shape_new);
1369
1370 /* create lists of local position and destination rank */
1371 const size_t local_size = T.local_size();
1372 std::vector<int> dest_mpirank(local_size);
1373 std::vector<size_t> local_position(local_size);
1374
1375 T.prep_local_to_global();
1376 T_new.prep_global_to_local();
1377
1378#pragma omp parallel default(shared)
1379 {
1380 Index index;
1381 index.resize(T.rank());
1382#pragma omp for
1383 for (size_t i = 0; i < local_size; ++i) {
1384 T.global_index_fast(i, index);
1385
1386 int dest;
1387 size_t pos;
1388 T_new.local_position(index, dest, pos);
1389
1390 local_position[i] = pos;
1391 dest_mpirank[i] = dest;
1392 }
1393 }
1394
1395 /* exchange data */
1396 replace_matrix_data(T.get_matrix(), dest_mpirank, local_position,
1397 T_new.get_matrix());
1398
1399 return T_new;
1400}
1401
1403
1408template <template <typename> class Matrix, typename C>
1410 assert(M.rank() == 2);
1411 return matrix_trace(M.get_matrix());
1412};
1413
1415
1428template <template <typename> class Matrix, typename C>
1429C trace(const Tensor<Matrix, C> &T, const Axes &axes_1, const Axes &axes_2) {
1430 assert(axes_1.size() == axes_2.size());
1431 assert(axes_1.size() + axes_2.size() == T.rank());
1432 assert(debug::check_trace_axes(axes_1, axes_2, T.rank()));
1433
1434 if (T.rank() == 2) return trace(T);
1435
1436 const size_t n = T.local_size();
1437 const size_t l = axes_1.size();
1438 Index index;
1439 bool check;
1440 C sum(0.0);
1441
1442 index.resize(T.rank());
1443 for (size_t i = 0; i < n; ++i) {
1444 T.global_index_fast(i, index);
1445 check = true;
1446 for (size_t k = 0; k < l; ++k) {
1447 if (index[axes_1[k]] != index[axes_2[k]]) {
1448 check = false;
1449 break;
1450 }
1451 }
1452 if (check) {
1453 sum += T[i];
1454 }
1455 }
1456
1457 return T.get_matrix().allreduce_sum(sum);
1458};
1459
1461
1474template <template <typename> class Matrix, typename C>
1476 const Axes &axes_a, const Axes &axes_b) {
1477 assert(A.rank() == B.rank());
1478 assert(A.rank() == axes_a.size());
1479 assert(B.rank() == axes_b.size());
1480 assert(debug::check_trace_axes(axes_a, axes_b, A.shape(), B.shape()));
1481
1482 const size_t rank = A.rank();
1483 Axes axes;
1485 Axes axes_map = A.get_axes_map();
1486 axes.resize(rank);
1487 axes_a_inv.resize(rank);
1488
1489 for (size_t i = 0; i < rank; ++i) {
1490 axes_a_inv[axes_a[i]] = i;
1491 }
1492 for (size_t i = 0; i < rank; ++i) {
1493 axes[i] = axes_b[axes_a_inv[axes_map[i]]];
1494 }
1495
1496 const size_t n = A.local_size();
1497 Tensor<Matrix, C> B_t = transpose(B, axes, A.get_upper_rank());
1498 C sum(0.0);
1499
1500 for (size_t i = 0; i < n; ++i) {
1501 sum += A[i] * B_t[i];
1502 }
1503
1504 return A.get_matrix().allreduce_sum(sum);
1505}
1506
1508
1518template <template <typename> class Matrix, typename C>
1520 const Axes &axes_2) {
1521 const int mpisize = T.get_comm_size();
1522 assert(axes_1.size() == axes_2.size());
1523 assert(axes_1.size() + axes_2.size() < T.rank());
1524 if (axes_1.size() == 0) return T;
1525 assert(debug::check_contract_axes(axes_1, axes_2, T.rank()));
1526
1527 Shape shape = T.shape();
1529 Axes axes_new;
1530 {
1531 Axes v = axes_1 + axes_2;
1532 v.sort();
1533 size_t k = 0;
1534 size_t n = v.size();
1535 for (size_t i = 0; i < T.rank(); ++i) {
1536 if (k < n && i == v[k]) {
1537 ++k;
1538 } else {
1539 shape_new.push(shape[i]);
1540 axes_new.push(i);
1541 }
1542 }
1543 }
1544
1545 /* initialize new tensor */
1546 Tensor<Matrix, C> T_new(T.get_comm(), shape_new);
1547
1548 /* create lists of local position and destination rank */
1549 std::vector<int> dest_mpirank(T.local_size());
1550 std::vector<size_t> local_position(T.local_size());
1551
1552 Index index, index_new;
1553 index.resize(T.rank());
1554 index_new.resize(T_new.rank());
1555 bool check;
1556 const size_t n = axes_1.size();
1557
1558 for (size_t i = 0; i < T.local_size(); ++i) {
1559 T.global_index_fast(i, index);
1560 check = true;
1561 for (size_t k = 0; k < n; ++k) {
1562 if (index[axes_1[k]] != index[axes_2[k]]) {
1563 check = false;
1564 continue;
1565 }
1566 }
1567
1568 if (check) {
1569 for (size_t k = 0; k < index_new.size(); ++k) {
1570 index_new[k] = index[axes_new[k]];
1571 }
1572
1573 int dest;
1574 size_t pos;
1575 T_new.local_position(index_new, dest, pos);
1576
1577 local_position[i] = pos;
1578 dest_mpirank[i] = dest;
1579 } else {
1580 // not send
1581 local_position[i] = 0;
1582 dest_mpirank[i] = mpisize;
1583 }
1584 }
1585
1586 /* sum data */
1587 sum_matrix_data(T.get_matrix(), dest_mpirank, local_position,
1588 T_new.get_matrix());
1589
1590 return T_new;
1591}
1592
1594
1611template <template <typename> class Matrix, typename C>
1613 assert(a.rank() == b.rank());
1614 assert(a.get_comm() == b.get_comm());
1615
1616 Shape shape_a = a.shape();
1617 Shape shape_b = b.shape();
1619 size_t n = shape_a.size();
1621 axes_trans.resize(2 * n);
1622 for (size_t i = 0; i < shape_b.size(); ++i) {
1623 shape_c[i] *= shape_b[i];
1624 axes_trans[2 * i] = i;
1625 axes_trans[2 * i + 1] = i + n;
1626 }
1627
1628 return reshape(
1630 reshape(b, Shape(1) + shape_b), Axes(n), Axes(0)),
1631 axes_trans),
1632 shape_c);
1633};
1634
1636
1647template <template <typename> class Matrix, typename C>
1649 const Tensor<Matrix, C> &b, const Axes &axes_a,
1650 const Axes &axes_b) {
1651 assert(axes_a.size() == axes_b.size());
1652 assert(a.get_comm() == b.get_comm());
1653 Shape shape_a = a.shape();
1654 Shape shape_b = b.shape();
1655 for (size_t i = 0; i < axes_a.size(); ++i) {
1657 }
1658
1659 const typename Tensor<Matrix, C>::comm_type &comm = a.get_comm();
1660
1661 Shape shape_c;
1662 const size_t rank_row_c = a.rank() - axes_a.size();
1663 const size_t rank_col_c = b.rank() - axes_b.size();
1664 shape_c.resize(rank_row_c + rank_col_c);
1665
1668 size_t urank_a;
1669 size_t urank_b;
1670
1671 {
1672 const size_t rank = a.rank();
1673 const size_t rank_row = rank - axes_a.size();
1674 const size_t rank_col = axes_a.size();
1675 size_t v[rank];
1676 for (size_t i = 0; i < rank; ++i) v[i] = i;
1677 for (size_t i = 0; i < rank_col; ++i) v[axes_a[i]] = rank;
1678 std::sort(v, v + rank);
1679 for (size_t i = 0; i < rank_col; ++i) v[rank_row + i] = axes_a[i];
1680
1681 trans_axes_a.assign(rank, v);
1682 urank_a = rank_row;
1683
1684 for (size_t i = 0; i < rank_row; ++i) shape_c[i] = shape_a[v[i]];
1685 }
1686
1687 {
1688 const size_t rank = b.rank();
1689 const size_t rank_row = axes_b.size();
1690 const size_t rank_col = rank - axes_b.size();
1691 size_t v[rank];
1692 for (size_t i = 0; i < rank; ++i) v[i] = i;
1693 for (size_t i = 0; i < rank_row; ++i) v[axes_b[i]] = 0;
1694 std::sort(v, v + rank);
1695 for (size_t i = 0; i < rank_row; ++i) v[i] = axes_b[i];
1696
1697 trans_axes_b.assign(rank, v);
1698 urank_b = rank_row;
1699
1700 for (size_t i = 0; i < rank_col; ++i)
1702 }
1703
1705 matrix_product(transpose(a, trans_axes_a, urank_a).get_matrix(),
1706 transpose(b, trans_axes_b, urank_b).get_matrix(),
1707 c.get_matrix());
1708
1709 return c;
1710};
1711
1713
1721template <template <typename> class Matrix, typename C>
1722int svd(const Tensor<Matrix, C> &a, std::vector<double> &s) {
1723 assert(a.rank() == 2);
1724 int info;
1725 info = svd(a, Axes(0), Axes(1), s);
1726 return info;
1727}
1728
1730
1742template <template <typename> class Matrix, typename C>
1744 std::vector<double> &s, Tensor<Matrix, C> &vt) {
1745 assert(a.rank() == 2);
1746 int info;
1747 info = svd(a, Axes(0), Axes(1), u, s, vt);
1748 return info;
1749}
1750
1752
1762template <template <typename> class Matrix, typename C>
1763int svd(const Tensor<Matrix, C> &a, const Axes &axes_row, const Axes &axes_col,
1764 std::vector<double> &s) {
1765 assert(axes_row.size() > 0);
1766 assert(axes_col.size() > 0);
1767 assert(debug::check_svd_axes(axes_row, axes_col, a.rank()));
1768
1769 const size_t rank = a.rank();
1770 const size_t urank = axes_row.size();
1771
1773
1775 const Shape &shape = a_t.shape();
1776
1777 size_t d_row(1), d_col(1);
1778 for (size_t i = 0; i < urank; ++i) d_row *= shape[i];
1779 for (size_t i = urank; i < rank; ++i) d_col *= shape[i];
1780 size_t size = (d_row < d_col) ? d_row : d_col;
1781
1782 s.resize(size);
1783
1784 int info;
1785 info = matrix_svd(a_t.get_matrix(), s);
1786
1787 return info;
1788}
1789
1791
1810template <template <typename> class Matrix, typename C>
1811int svd(const Tensor<Matrix, C> &a, const Axes &axes_row, const Axes &axes_col,
1812 Tensor<Matrix, C> &u, std::vector<double> &s, Tensor<Matrix, C> &vt) {
1813 assert(axes_row.size() > 0);
1814 assert(axes_col.size() > 0);
1815 assert(debug::check_svd_axes(axes_row, axes_col, a.rank()));
1816
1817 const size_t rank = a.rank();
1818 const size_t urank = axes_row.size();
1819
1821
1823 const Shape &shape = a_t.shape();
1824
1825 size_t d_row(1), d_col(1);
1826 for (size_t i = 0; i < urank; ++i) d_row *= shape[i];
1827 for (size_t i = urank; i < rank; ++i) d_col *= shape[i];
1828 size_t size = (d_row < d_col) ? d_row : d_col;
1829
1830 Shape shape_u;
1831 shape_u.resize(urank + 1);
1832 for (size_t i = 0; i < urank; ++i) shape_u[i] = shape[i];
1833 shape_u[urank] = size;
1834
1836 shape_vt.resize(rank - urank + 1);
1837 shape_vt[0] = size;
1838 for (size_t i = urank; i < rank; ++i) shape_vt[i - urank + 1] = shape[i];
1839
1840 size_t urank_u = urank;
1841 size_t urank_vt = 1;
1842
1845 s.resize(size);
1846
1847 int info;
1848 info = matrix_svd(a_t.get_matrix(), u.get_matrix(), s, vt.get_matrix());
1849
1850 return info;
1851}
1852
1854
1865template <template <typename> class Matrix, typename C>
1866int psvd(const Tensor<Matrix, C> &a, std::vector<double> &s,
1867 const size_t target_rank) {
1868 int info = svd(a, s);
1869 if (s.size() > target_rank) s.resize(target_rank);
1870 return info;
1871}
1872
1874
1889template <template <typename> class Matrix, typename C>
1891 std::vector<double> &s, Tensor<Matrix, C> &vt,
1892 const size_t target_rank) {
1893 int info = svd(a, u, s, vt);
1894 if (s.size() > target_rank) {
1895 s.resize(target_rank);
1896 u = slice(u, 1, 0, target_rank);
1897 vt = slice(vt, 0, 0, target_rank);
1898 }
1899 return info;
1900}
1901
1903
1916template <template <typename> class Matrix, typename C>
1917int psvd(const Tensor<Matrix, C> &a, const Axes &axes_row, const Axes &axes_col,
1918 std::vector<double> &s, const size_t target_rank) {
1919 int info = svd(a, axes_row, axes_col, s);
1920 if (s.size() > target_rank) s.resize(target_rank);
1921 return info;
1922}
1923
1925
1947template <template <typename> class Matrix, typename C>
1948int psvd(const Tensor<Matrix, C> &a, const Axes &axes_row, const Axes &axes_col,
1949 Tensor<Matrix, C> &u, std::vector<double> &s, Tensor<Matrix, C> &vt,
1950 const size_t target_rank) {
1951 int info = svd(a, axes_row, axes_col, u, s, vt);
1952 if (s.size() > target_rank) {
1953 s.resize(target_rank);
1954 u = slice(u, axes_row.size(), 0, target_rank);
1955 vt = slice(vt, 0, 0, target_rank);
1956 }
1957 return info;
1958}
1959
1961
1973template <template <typename> class Matrix, typename C>
1975 assert(a.rank() == 2);
1976 int info;
1977 info = qr(a, Axes(0), Axes(1), q, r);
1978 return info;
1979}
1980
1982
2005template <template <typename> class Matrix, typename C>
2006int qr(const Tensor<Matrix, C> &a, const Axes &axes_row, const Axes &axes_col,
2008 assert(axes_row.size() > 0);
2009 assert(axes_col.size() > 0);
2010 assert(debug::check_svd_axes(axes_row, axes_col, a.rank()));
2011
2012 const size_t rank = a.rank();
2013 const size_t urank = axes_row.size();
2014 const Shape shape_a = a.shape();
2015
2017 Shape shape;
2018 shape.resize(rank);
2019 for (size_t i = 0; i < rank; ++i) shape[i] = shape_a[axes[i]];
2020
2021 size_t d_row(1), d_col(1);
2022 for (size_t i = 0; i < urank; ++i) d_row *= shape[i];
2023 for (size_t i = urank; i < rank; ++i) d_col *= shape[i];
2024 size_t size = (d_row < d_col) ? d_row : d_col;
2025
2026 // rank-2 tensors (matrices)
2029 Tensor<Matrix, C> mat_r(mat_q.get_comm(), mat_q.shape(), 1);
2030
2031 // QR decomposition. Elements of mat_q change from a to q.
2032 int info;
2033 info = matrix_qr(mat_q.get_matrix(), mat_r.get_matrix());
2034
2035 // Get shape of q and r.
2036 Shape shape_q;
2037 shape_q.resize(urank + 1);
2038 for (size_t i = 0; i < urank; ++i) shape_q[i] = shape[i];
2039 shape_q[urank] = size;
2040
2041 Shape shape_r;
2042 shape_r.resize(rank - urank + 1);
2043 shape_r[0] = size;
2044 for (size_t i = urank; i < rank; ++i) shape_r[i - urank + 1] = shape[i];
2045
2046 // Reshape q and r.
2047 if (d_row > d_col) {
2048 q = reshape(mat_q, shape_q);
2049 r = reshape(slice(mat_r, 0, 0, size), shape_r);
2050 } else if (d_row < d_col) {
2051 q = reshape(slice(mat_q, 1, 0, size), shape_q);
2052 r = reshape(mat_r, shape_r);
2053 } else { // d_row==d_col
2054 q = reshape(mat_q, shape_q);
2055 r = reshape(mat_r, shape_r);
2056 }
2057
2058 return info;
2059}
2060
2063
2070template <template <typename> class Matrix, typename C>
2071int eigh(const Tensor<Matrix, C> &a, std::vector<double> &w,
2073 assert(a.rank() == 2);
2074 Shape shape = a.shape();
2075 assert(shape[0] == shape[1]);
2076
2077 Tensor<Matrix, C> a_t = transpose(a, Axes(0, 1), 1);
2078
2079 size_t n = shape[0];
2080 w.resize(n);
2081 z = Tensor<Matrix, C>(a.get_comm(), Shape(n, n), 1);
2082
2083 int info;
2084 info = matrix_eigh(a_t.get_matrix(), w, z.get_matrix());
2085 return info;
2086};
2087
2090
2096template <template <typename> class Matrix, typename C>
2097int eigh(const Tensor<Matrix, C> &a, std::vector<double> &w) {
2098 assert(a.rank() == 2);
2099 Shape shape = a.shape();
2100 assert(shape[0] == shape[1]);
2101
2102 Tensor<Matrix, C> a_t = transpose(a, Axes(0, 1), 1);
2103
2104 size_t n = shape[0];
2105 w.resize(n);
2106
2107 int info;
2108 info = matrix_eigh(a_t.get_matrix(), w);
2109 return info;
2110};
2111
2114
2134template <template <typename> class Matrix, typename C>
2135int eigh(const Tensor<Matrix, C> &a, const Axes &axes_row, const Axes &axes_col,
2136 std::vector<double> &w, Tensor<Matrix, C> &z) {
2137 assert(axes_row.size() > 0);
2138 assert(axes_col.size() > 0);
2139
2140 const size_t rank = a.rank();
2141 const size_t urank = axes_row.size();
2142
2144
2146 const Shape &shape = a_t.shape();
2147
2148 size_t d_row(1), d_col(1);
2149 for (size_t i = 0; i < urank; ++i) d_row *= shape[i];
2150 for (size_t i = urank; i < rank; ++i) d_col *= shape[i];
2151 size_t size = (d_row < d_col) ? d_row : d_col;
2152
2153 assert(d_row == d_col);
2154
2155 Shape shape_z;
2156 shape_z.resize(urank + 1);
2157 for (size_t i = 0; i < urank; ++i) shape_z[i] = shape[i];
2158 shape_z[urank] = size;
2159
2161 w.resize(size);
2162
2163 int info;
2164 info = matrix_eigh(a_t.get_matrix(), w, z.get_matrix());
2165 return info;
2166};
2167
2169
2179template <template <typename> class Matrix, typename C>
2180int eigh(const Tensor<Matrix, C> &a, const Axes &axes_row, const Axes &axes_col,
2181 std::vector<double> &w) {
2182 assert(axes_row.size() > 0);
2183 assert(axes_col.size() > 0);
2184
2185 const size_t rank = a.rank();
2186 const size_t urank = axes_row.size();
2187
2189
2191 const Shape &shape = a_t.shape();
2192
2193 size_t d_row(1), d_col(1);
2194 for (size_t i = 0; i < urank; ++i) d_row *= shape[i];
2195 for (size_t i = urank; i < rank; ++i) d_col *= shape[i];
2196 size_t size = (d_row < d_col) ? d_row : d_col;
2197 w.resize(size);
2198
2199 assert(d_row == d_col);
2200
2201 int info;
2202 info = matrix_eigh(a_t.get_matrix(), w);
2203 return info;
2204};
2205
2208
2239template <template <typename> class Matrix, typename C>
2241 const Axes &axes_col_a, const Tensor<Matrix, C> &b,
2242 const Axes &axes_row_b, const Axes &axes_col_b, std::vector<double> &w,
2244 assert(axes_row_a.size() > 0);
2245 assert(axes_col_a.size() > 0);
2246 assert(axes_row_b.size() > 0);
2247 assert(axes_col_b.size() > 0);
2248
2249 const size_t rank_a = a.rank();
2250 const size_t rank_b = b.rank();
2251 const size_t urank_a = axes_row_a.size();
2252 const size_t urank_b = axes_row_b.size();
2257 const Shape &shape_a = a_t.shape();
2258 const Shape &shape_b = b_t.shape();
2259 size_t d_row_a(1), d_col_a(1);
2260 size_t d_row_b(1), d_col_b(1);
2261 for (size_t i = 0; i < urank_a; ++i) d_row_a *= shape_a[i];
2262 for (size_t i = 0; i < urank_b; ++i) d_row_b *= shape_b[i];
2263 for (size_t i = urank_a; i < rank_a; ++i) d_col_a *= shape_a[i];
2264 for (size_t i = urank_b; i < rank_b; ++i) d_col_b *= shape_b[i];
2265
2269 size_t size = d_row_a;
2270
2271 Shape shape_z;
2273 for (size_t i = urank_a; i < rank_a; ++i) shape_z[i - urank_a] = shape_a[i];
2274 shape_z[rank_a - urank_a] = size;
2275
2277 w.resize(size);
2278
2279 int info;
2280 info = matrix_eigh(a_t.get_matrix(), b_t.get_matrix(), w, z.get_matrix());
2281
2282 return info;
2283};
2284
2286
2294template <template <typename> class Matrix, typename C>
2295int eig(const Tensor<Matrix, C> &a, std::vector<complex> &w,
2297 assert(a.rank() == 2);
2298 Shape shape = a.shape();
2299 assert(shape[0] == shape[1]);
2300
2301 Tensor<lapack::Matrix, C> a_t = transpose(a, Axes(0, 1), 1).gather();
2302
2303 size_t n = shape[0];
2304 w.resize(n);
2305
2307
2308 int info;
2309 info = matrix_eig(a_t.get_matrix(), w, z_t.get_matrix());
2310
2312
2313 return info;
2314};
2315
2317
2324template <template <typename> class Matrix, typename C>
2325int eig(const Tensor<Matrix, C> &a, std::vector<complex> &w) {
2326 assert(a.rank() == 2);
2327 Shape shape = a.shape();
2328 assert(shape[0] == shape[1]);
2329
2330 Tensor<lapack::Matrix, C> a_t = transpose(a, Axes(0, 1), 1).gather();
2331
2332 size_t n = shape[0];
2333 w.resize(n);
2334
2335 int info;
2336 info = matrix_eig(a_t.get_matrix(), w);
2337 return info;
2338};
2339
2341
2361template <template <typename> class Matrix, typename C>
2362int eig(const Tensor<Matrix, C> &a, const Axes &axes_row, const Axes &axes_col,
2363 std::vector<complex> &w, Tensor<Matrix, complex> &z) {
2364 assert(axes_row.size() > 0);
2365 assert(axes_col.size() > 0);
2366
2367 const size_t rank = a.rank();
2368 const size_t urank = axes_row.size();
2369
2371
2373 const Shape &shape = a_t.shape();
2374
2375 size_t d_row(1), d_col(1);
2376 for (size_t i = 0; i < urank; ++i) d_row *= shape[i];
2377 for (size_t i = urank; i < rank; ++i) d_col *= shape[i];
2378 size_t size = (d_row < d_col) ? d_row : d_col;
2379
2380 assert(d_row == d_col);
2381
2382 Shape shape_z;
2383 shape_z.resize(urank + 1);
2384 for (size_t i = 0; i < urank; ++i) shape_z[i] = shape[i];
2385 shape_z[urank] = size;
2386
2388
2389 w.resize(size);
2390
2391 int info;
2392 info = matrix_eig(a_t.get_matrix(), w, z_t.get_matrix());
2393
2395
2396 return info;
2397};
2398
2400
2410template <template <typename> class Matrix, typename C>
2411int eig(const Tensor<Matrix, C> &a, const Axes &axes_row, const Axes &axes_col,
2412 std::vector<complex> &w) {
2413 assert(axes_row.size() > 0);
2414 assert(axes_col.size() > 0);
2415
2416 const size_t rank = a.rank();
2417 const size_t urank = axes_row.size();
2418
2420
2422 const Shape &shape = a_t.shape();
2423
2424 size_t d_row(1), d_col(1);
2425 for (size_t i = 0; i < urank; ++i) d_row *= shape[i];
2426 for (size_t i = urank; i < rank; ++i) d_col *= shape[i];
2427 size_t size = (d_row < d_col) ? d_row : d_col;
2428 w.resize(size);
2429
2430 assert(d_row == d_col);
2431
2432 int info;
2433 info = matrix_eig(a_t.get_matrix(), w);
2434 return info;
2435};
2436
2437
2439
2446template <template <typename> class Matrix, typename C>
2447int solve(const Tensor<Matrix, C> &a, const std::vector<C> &b,
2448 std::vector<C> &x) {
2449 assert(a.rank() == 2);
2450 assert(a.shape()[0] == a.shape()[1]);
2451 assert(a.shape()[0] == b.size());
2452
2454 int info;
2455 info = solve(a, Tensor<Matrix, C>(a.get_comm(), b), x_t, Axes(0), Axes(1),
2456 Axes(0), Axes());
2457 x = x_t.flatten();
2458
2459 return info;
2460};
2461
2463
2470template <template <typename> class Matrix, typename C>
2473 assert(a.rank() == 2);
2474 assert(b.rank() == 2 || b.rank() == 1);
2475 assert(a.shape()[0] == a.shape()[1]);
2476 assert(a.shape()[0] == b.shape()[0]);
2477 int info;
2478
2479 if (b.rank() == 2) {
2480 info = solve(a, b, x, Axes(0), Axes(1), Axes(0), Axes(1));
2481 } else {
2482 info = solve(a, b, x, Axes(0), Axes(1), Axes(0), Axes());
2483 }
2484
2485 return info;
2486};
2487
2489
2516template <template <typename> class Matrix, typename C>
2519 const Axes &axes_row_b, const Axes &axes_col_b) {
2520 assert(a.rank() == axes_row_a.size() + axes_col_a.size());
2521 assert(b.rank() == axes_row_b.size() + axes_col_b.size());
2522 assert(axes_row_a.size() > 0);
2523 assert(axes_col_a.size() > 0);
2524 assert(axes_row_b.size() > 0);
2525 assert(axes_col_b.size() >= 0);
2526
2527 // size_t rank_a = a.rank();
2528 // size_t rank_b = b.rank();
2529 size_t rank_row_a = axes_row_a.size();
2530 size_t rank_col_a = axes_col_a.size();
2531 size_t rank_row_b = axes_row_b.size();
2532 size_t rank_col_b = axes_col_b.size();
2535
2538 const Shape &shape_a = a_t.shape();
2539 const Shape &shape_b = b_t.shape();
2540
2541 assert(debug::check_square(shape_a, rank_row_a));
2542
2543 Shape shape_x;
2545 for (size_t i = 0; i < rank_col_a; ++i) shape_x[i] = shape_a[i + rank_row_a];
2546 for (size_t i = 0; i < rank_col_b; ++i)
2548
2549 int info;
2550 info = matrix_solve(a_t.get_matrix(), b_t.get_matrix());
2551 x = reshape(b_t, shape_x);
2552 return info;
2553};
2554
2556template <template <typename> class Matrix, typename C>
2561template <template <typename> class Matrix, typename C>
2563 return (rhs *= -1.0);
2564}
2566template <template <typename> class Matrix, typename C>
2568 const Tensor<Matrix, C> &rhs) {
2569 return (lhs += rhs);
2570}
2572template <template <typename> class Matrix, typename C>
2574 const Tensor<Matrix, C> &rhs) {
2575 return (lhs -= rhs);
2576}
2578template <template <typename> class Matrix, typename C, typename D>
2580 return (lhs *= rhs);
2581}
2583template <template <typename> class Matrix, typename C, typename D>
2585 return (lhs /= rhs);
2586}
2588template <template <typename> class Matrix, typename C, typename D>
2590 return (rhs *= lhs);
2591}
2592
2594template <template <typename> class Matrix, typename C>
2595inline double max(const Tensor<Matrix, C> &t) {
2596 return max(t.get_matrix());
2597}
2599template <template <typename> class Matrix, typename C>
2600inline double min(const Tensor<Matrix, C> &t) {
2601 return min(t.get_matrix());
2602}
2604template <template <typename> class Matrix, typename C>
2605inline double max_abs(const Tensor<Matrix, C> &t) {
2606 return max_abs(t.get_matrix());
2607}
2609template <template <typename> class Matrix, typename C>
2610inline double min_abs(const Tensor<Matrix, C> &t) {
2611 return min_abs(t.get_matrix());
2612}
2613
2615template <template <typename> class Matrix>
2617 return t.map(static_cast<double (*)(double)>(&std::sqrt));
2618}
2619template <template <typename> class Matrix>
2621 return t.map(static_cast<complex (*)(const complex &)>(&std::sqrt));
2622}
2623template <template <typename> class Matrix>
2625 return t;
2626}
2627template <template <typename> class Matrix>
2629 return t.map(static_cast<complex (*)(const complex &)>(&std::conj));
2630}
2632
2634
2638template <template <typename> class Matrix, typename C>
2639std::ostream &operator<<(std::ostream &out, const Tensor<Matrix, C> &t) {
2640 std::size_t dim = t.ndim();
2641 if (t.get_comm_size() == 1) {
2642 std::size_t size = t.local_size();
2643 Shape shape = t.shape();
2644 std::vector<std::size_t> accum(dim + 1, 1);
2645 for (size_t d = dim; d != 0; --d) accum[d - 1] = accum[d] * shape[d - 1];
2646
2647 std::vector<std::size_t> idx(dim);
2648 for (size_t i = 0; i < size; ++i) {
2649 for (size_t d = 0; d < dim; ++d) idx[d] = (i % accum[d]) / accum[d + 1];
2650 if (idx[dim - 1] == 0) {
2651 size_t count = 0;
2652 for (size_t d = dim; d != 0 && idx[d - 1] == 0; --d) ++count;
2653 for (size_t i = 0; i < dim - count; ++i) out << ' ';
2654 for (size_t i = 0; i < count; ++i) out << '[';
2655 }
2656
2657 C val;
2658 t.get_value(idx, val);
2659 out << ' ' << val << ' ';
2660
2661 if (idx[dim - 1] == shape[dim - 1] - 1) {
2662 size_t count = 0;
2663 for (size_t d = dim; d != 0 && idx[d - 1] == shape[d - 1] - 1; --d) ++count;
2664 for (size_t i = 0; i < count; ++i) out << ']';
2665 if (count < dim) {
2666 for (size_t i = 0; i < count; ++i) out << '\n';
2667 } else {
2668 out << '\n';
2669 }
2670 }
2671 }
2672 } else {
2673 if (t.get_comm_rank() == 0) {
2674 for (size_t d = 0; d < dim; ++d) out << '[';
2675 out << " distributed tensor ... ";
2676 for (size_t d = 0; d < dim; ++d) out << ']';
2677 out << '\n';
2678 }
2679 }
2680 return out;
2681}
2682
2683} // namespace mptensor
2684
2685#endif // _TENSOR_IMPL_HPP_
Definition index.hpp:39
void resize(size_t n)
Definition index.hpp:81
void sort()
Definition index.cc:42
size_t size() const
Definition index.hpp:79
Tensor class. The main object of mptensor.
Definition tensor.hpp:54
void print_info(std::ostream &out, const std::string &tag="") const
Output information of tensor.
Definition tensor_impl.hpp:424
size_t rank() const
Rank of tensor.
Definition tensor_impl.hpp:164
void global_index_fast(size_t i, Index &idx) const
Convert local index to global index fast.
Definition tensor_impl.hpp:525
size_t ndim() const
Rank of tensor.
Definition tensor_impl.hpp:173
void prep_local_to_global() const
Preprocess for fast conversion from local index to global one.
Definition tensor_impl.hpp:258
Tensor< Matrix, C > & map(UnaryOperation op)
Apply a unary operation to each element.
Definition tensor_impl.hpp:1026
void prep_global_to_local() const
Preprocess for fast conversion from global index to local one.
Definition tensor_impl.hpp:252
const Matrix< C > & get_matrix() const
distributed Matrix
Definition tensor_impl.hpp:203
void local_position(const Index &idx, int &comm_rank, size_t &local_idx) const
Calculate rank which has the given global index and its local index.
Definition tensor_impl.hpp:584
size_t local_size() const
Size of local storage.
Definition tensor_impl.hpp:179
int get_comm_rank() const
Rank of process.
Definition tensor_impl.hpp:228
int get_comm_size() const
Size of communicator.
Definition tensor_impl.hpp:222
const comm_type & get_comm() const
Communicator.
Definition tensor_impl.hpp:216
scalapack::Matrix< double >::comm_type comm_type
type of communicator.
Definition tensor.hpp:58
const Shape & shape() const
Shape of tensor.
Definition tensor_impl.hpp:158
size_t get_upper_rank() const
Upper rank for matrix representation.
Definition tensor_impl.hpp:185
bool get_value(const Index &idx, C &val) const
Get an element.
Definition tensor_impl.hpp:552
Define the type of a complex number.
Tensor< Matrix, C > operator*(Tensor< Matrix, C > lhs, D rhs)
Tensor-scalar multiplication.
Definition tensor_impl.hpp:2579
Tensor< Matrix, C > operator/(Tensor< Matrix, C > lhs, D rhs)
Scalar division.
Definition tensor_impl.hpp:2584
std::complex< double > complex
Definition complex.hpp:38
int eig(const Tensor< Matrix, C > &a, std::vector< complex > &eigval, Tensor< Matrix, complex > &eigvec)
Compute the eigenvalues and eigenvectors of a general matrix (rank-2 tensor)
Definition tensor_impl.hpp:2295
int svd(const Tensor< Matrix, C > &a, std::vector< double > &s)
Singular value decomposition for rank-2 tensor (matrix)
Definition tensor_impl.hpp:1722
int eigh(const Tensor< Matrix, C > &a, std::vector< double > &eigval, Tensor< Matrix, C > &eigvec)
Definition tensor_impl.hpp:2071
int qr(const Tensor< Matrix, C > &a, Tensor< Matrix, C > &q, Tensor< Matrix, C > &r)
QR decomposition of matrix (rank-2 tensor)
Definition tensor_impl.hpp:1974
int psvd(const Tensor< Matrix, C > &a, std::vector< double > &s, const size_t target_rank)
Partial SVD for rank-2 tensor (matrix)
Definition tensor_impl.hpp:1866
Tensor< Matrix, C > tensordot(const Tensor< Matrix, C > &a, const Tensor< Matrix, C > &b, const Axes &axes_a, const Axes &axes_b)
Compute tensor dot product.
Definition tensor_impl.hpp:1648
Tensor< Matrix, C > contract(const Tensor< Matrix, C > &a, const Axes &axes_1, const Axes &axes_2)
Partial trace of tensor.
Definition tensor_impl.hpp:1519
Tensor< Matrix, C > kron(const Tensor< Matrix, C > &a, const Tensor< Matrix, C > &b)
Compute the Kronecker product.
Definition tensor_impl.hpp:1612
C trace(const Tensor< Matrix, C > &a)
Trace of matrix (rank-2 tensor)
Definition tensor_impl.hpp:1409
int solve(const Tensor< Matrix, C > &a, const std::vector< C > &b, std::vector< C > &x)
Solve linear equation .
Definition tensor_impl.hpp:2447
Tensor< Matrix, C > sqrt(Tensor< Matrix, C > t)
Take square-root of each element.
Tensor< Matrix, C > conj(Tensor< Matrix, C > t)
Take conjugate of each element.
Tensor< Matrix, C > transpose(Tensor< Matrix, C > a, const Axes &axes)
Transposition of tensor with lazy evaluation.
Definition tensor_impl.hpp:1041
Tensor< Matrix, C > slice(const Tensor< Matrix, C > &a, size_t n_axes, size_t i_begin, size_t i_end)
Slice a tensor.
Definition tensor_impl.hpp:1215
Tensor< Matrix, C > reshape(const Tensor< Matrix, C > &a, const Shape &shape_new)
Change the shape of tensor.
Definition tensor_impl.hpp:1149
Tensor< Matrix, C > extend(const Tensor< Matrix, C > &a, const Shape &shape_new)
Extend the size of a tensor.
Definition tensor_impl.hpp:1363
Tensor()
Default constructor of tensor.
Definition tensor_impl.hpp:70
header file of Index class
Index(size_t j0)
Definition index_constructor.hpp:11
List of header files for matrix classes.
mptensor::complex complex
Definition matrix_lapack.cc:36
bool check_total_size(const Shape &s1, const Shape &s2)
Definition tensor.cc:49
bool check_trace_axes(const Axes &axes_1, const Axes &axes_2, size_t rank)
Definition tensor.cc:86
bool check_extend(const Shape &s_old, const Shape &s_new)
Definition tensor.cc:57
bool check_svd_axes(const Axes &axes_row, const Axes &axes_col, size_t rank)
Definition tensor.cc:76
bool check_contract_axes(const Axes &axes_1, const Axes &axes_2, size_t rank)
Definition tensor.cc:114
bool check_transpose_axes(const Axes &axes, size_t rank)
Definition tensor.cc:66
bool check_square(const Shape &shape, size_t urank)
Definition tensor.cc:125
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
Definition complex.hpp:34
Index range(const size_t start, const size_t stop)
Create an increasing sequence. it is similar to range() in python.
Definition index.cc:91
Index Axes
Definition tensor.hpp:45
bool is_no_transpose(const Axes &axes, const Axes &axes_map, size_t rank)
Definition tensor.cc:35
Index operator+(const Index &, const Index &)
Tensor class.
double operator-(Timer &t1, Timer &t0)
Definition timer.hpp:64