29#ifndef _TENSOR_IMPL_HPP_
30#define _TENSOR_IMPL_HPP_
69template <
template <
typename>
class Matrix,
typename C>
78template <
template <
typename>
class Matrix,
typename C>
80 init(shape, shape.size() / 2);
87template <
template <
typename>
class Matrix,
typename C>
89 : Mat(comm), upper_rank(0), axes_map(){};
97template <
template <
typename>
class Matrix,
typename C>
100 init(shape, shape.size() / 2);
109template <
template <
typename>
class Matrix,
typename C>
111 const size_t upper_rank)
122template <
template <
typename>
class Matrix,
typename C>
127 const size_t n = Mat.local_size();
130 for (
size_t i = 0;
i <
n; ++
i) {
142template <
template <
typename>
class Matrix,
typename C>
146 const size_t n = Mat.local_size();
148 for (
size_t i = 0;
i <
n; ++
i) {
149 idx = global_index(
i);
157template <
template <
typename>
class Matrix,
typename C>
163template <
template <
typename>
class Matrix,
typename C>
172template <
template <
typename>
class Matrix,
typename C>
178template <
template <
typename>
class Matrix,
typename C>
184template <
template <
typename>
class Matrix,
typename C>
196template <
template <
typename>
class Matrix,
typename C>
202template <
template <
typename>
class Matrix,
typename C>
208template <
template <
typename>
class Matrix,
typename C>
214template <
template <
typename>
class Matrix,
typename C>
221template <
template <
typename>
class Matrix,
typename C>
227template <
template <
typename>
class Matrix,
typename C>
236template <
template <
typename>
class Matrix,
typename C>
245template <
template <
typename>
class Matrix,
typename C>
251template <
template <
typename>
class Matrix,
typename C>
257template <
template <
typename>
class Matrix,
typename C>
263template <
template <
typename>
class Matrix,
typename C>
265 const size_t rank = Dim.size();
266 const size_t rank0 = upper_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]);
276#pragma omp parallel default(shared)
288 for (
size_t k = 0;
k <
l_row; ++
k) {
289 map = &(l2g_map_row[
k *
rank0]);
291 for (
size_t i = 0;
i <
rank0; ++
i) {
298 for (
size_t k = 0;
k <
l_col; ++
k) {
299 map = &(l2g_map_col[
k *
rank1]);
301 for (
size_t i = 0;
i <
rank1; ++
i) {
319template <
template <
typename>
class Matrix,
typename C>
322 const size_t rank = Dim.size();
323 const size_t rank0 = upper_rank;
329 const size_t *
axes_map0 = &(axes_map[0]);
331 for (
size_t i = 0;
i <
rank0; ++
i) {
334 for (
size_t i = 0;
i <
rank1; ++
i) {
352template <
template <
typename>
class Matrix,
typename C>
355 const size_t rank = Dim.size();
356 const size_t rank0 = upper_rank;
364 for (
size_t i = 0;
i <
rank0; ++
i) {
368 for (
size_t i = 0;
i <
rank1; ++
i) {
375template <
template <
typename>
class Matrix,
typename C>
389template <
template <
typename>
class Matrix,
typename C>
400template <
template <
typename>
class Matrix,
typename C>
401void Tensor<Matrix, C>::init(
const Shape &shape,
size_t urank,
405 const size_t rank = Dim.size();
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);
423template <
template <
typename>
class Matrix,
typename C>
425 const std::string &
tag)
const {
426 if (
tag.size() > 0) {
429 out <<
"Tensor: shape= " << Dim <<
" upper_rank= " << upper_rank
430 <<
" axes_map= " << axes_map <<
"\t";
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();
448 for (
int i = 0;
i < mpisize; ++
i) {
450 if (
tag.size() > 0) {
453 out <<
"mpirank: " << mpirank <<
"\t";
469template <
template <
typename>
class Matrix,
typename C>
471 const size_t rank =
gindex.size();
472 assert(rank == Dim.size());
476 for (
size_t i = 0;
i < upper_rank; ++
i) {
477 size_t j = axes_map[
i];
481 for (
size_t i = upper_rank;
i < rank; ++
i) {
482 size_t j = axes_map[
i];
494template <
template <
typename>
class Matrix,
typename C>
496 const size_t rank = Dim.size();
502 for (
size_t i = 0;
i < upper_rank; ++
i) {
503 size_t j = axes_map[
i];
508 for (
size_t i = upper_rank;
i < rank; ++
i) {
509 size_t j = axes_map[
i];
524template <
template <
typename>
class Matrix,
typename C>
526 const size_t rank = Dim.size();
530 for (
size_t i = 0;
i < upper_rank; ++
i) {
531 size_t j = axes_map[
i];
536 for (
size_t i = upper_rank;
i < rank; ++
i) {
537 size_t j = axes_map[
i];
551template <
template <
typename>
class Matrix,
typename C>
554 if (local_index(idx,
li)) {
569template <
template <
typename>
class Matrix,
typename C>
572 if (local_index(idx,
li)) {
583template <
template <
typename>
class Matrix,
typename C>
586 const size_t rank = Dim.size();
589 for (
size_t i = 0;
i < upper_rank; ++
i) {
590 size_t j = axes_map[
i];
594 for (
size_t i = upper_rank;
i < rank; ++
i) {
595 size_t j = axes_map[
i];
610template <
template <
typename>
class Matrix,
typename C>
617 const size_t rank = this->rank();
622 for (
size_t i = 0;
i < rank; ++
i) {
630 const size_t local_size =
T_old.local_size();
632 std::vector<size_t> local_position(local_size);
634 T_old.prep_local_to_global();
635 this->prep_global_to_local();
637#pragma omp parallel default(shared)
643 for (
size_t i = 0;
i < local_size; ++
i) {
644 T_old.global_index_fast(
i, index);
648 this->local_position(index,
dest,
pos);
650 local_position[
i] =
pos;
666template <
template <
typename>
class Matrix,
typename C>
668 const size_t rank = Dim.size();
669 assert(debug::check_transpose_axes(
axes, rank));
675 for (
size_t i = 0;
i < rank; ++
i) {
679 for (
size_t i = 0;
i < rank; ++
i) {
698template <
template <
typename>
class Matrix,
typename C>
703 const size_t local_size = this->local_size();
704 prep_local_to_global();
705#pragma omp parallel default(shared)
710 for (
size_t i = 0;
i < local_size; ++
i) {
711 global_index_fast(
i, idx);
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,
738 const size_t local_size = this->local_size();
739 prep_local_to_global();
740#pragma omp parallel default(shared)
745 for (
size_t i = 0;
i < local_size; ++
i) {
746 global_index_fast(
i, idx);
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,
776 const size_t local_size = this->local_size();
777 prep_local_to_global();
778#pragma omp parallel default(shared)
783 for (
size_t i = 0;
i < local_size; ++
i) {
784 global_index_fast(
i, idx);
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,
818 const size_t local_size = this->local_size();
819 prep_local_to_global();
820#pragma omp parallel default(shared)
825 for (
size_t i = 0;
i < local_size; ++
i) {
826 global_index_fast(
i, idx);
845template <
template <
typename>
class Matrix,
typename C>
849 const size_t i_end) {
859 std::vector<size_t> local_position(local_size);
862 prep_global_to_local();
864#pragma omp parallel default(shared)
869 for (
size_t i = 0;
i < local_size; ++
i) {
875 this->local_position(index,
dest,
pos);
877 local_position[
i] =
pos;
900template <
template <
typename>
class Matrix,
typename C>
904 const size_t nr = rank();
908 for (
size_t r = 0;
r <
nr; ++
r) {
916 std::vector<size_t> local_position(local_size);
919 prep_global_to_local();
921#pragma omp parallel default(shared)
926 for (
size_t i = 0;
i < local_size; ++
i) {
928 for (
size_t r = 0;
r <
nr; ++
r) {
934 this->local_position(index,
dest,
pos);
936 local_position[
i] =
pos;
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));
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));
976 return get_matrix().flatten();
980template <
template <
typename>
class Matrix,
typename C>
983 change_configuration(
rhs.get_upper_rank(),
rhs.get_axes_map());
984 Mat +=
rhs.get_matrix();
989template <
template <
typename>
class Matrix,
typename C>
992 change_configuration(
rhs.get_upper_rank(),
rhs.get_axes_map());
993 Mat -=
rhs.get_matrix();
998template <
template <
typename>
class Matrix,
typename C>
1005template <
template <
typename>
class Matrix,
typename C>
1012template <
template <
typename>
class Matrix,
typename C>
1015 for (
int i = 0;
i <
n; ++
i) Mat[
i] =
rhs;
1024template <
template <
typename>
class Matrix,
typename C>
1025template <
typename UnaryOperation>
1040template <
template <
typename>
class Matrix,
typename C>
1042 return T.transpose(
axes);
1055template <
template <
typename>
class Matrix,
typename C>
1058 const size_t rank =
T.rank();
1059 assert(debug::check_transpose_axes(
axes, rank));
1070 for (
size_t r = 0;
r < rank; ++
r) {
1078 const size_t local_size =
T.local_size();
1080 std::vector<size_t> local_position(local_size);
1083 T_new.prep_global_to_local();
1085#pragma omp parallel default(shared)
1099 Axes axes_map =
T.get_axes_map();
1102 for (
size_t r = 0;
r < rank; ++
r) {
1105 for (
size_t r = 0;
r < rank; ++
r) {
1112 for (
size_t i = 0;
i < local_size; ++
i) {
1130 replace_matrix_data(
T.get_matrix(),
dest_mpirank, local_position,
1131 T_new.get_matrix());
1148template <
template <
typename>
class Matrix,
typename C>
1152 const Shape &shape =
T.shape();
1153 const size_t rank = shape.size();
1160 const size_t local_size =
T.local_size();
1162 std::vector<size_t> local_position(local_size);
1164 T.prep_local_to_global();
1165 T_new.prep_global_to_local();
1167#pragma omp parallel default(shared)
1173 for (
size_t i = 0;
i < local_size; ++
i) {
1174 T.global_index_fast(
i, index);
1177 for (
size_t r = 0;
r < rank; ++
r) {
1191 local_position[
i] =
pos;
1197 replace_matrix_data(
T.get_matrix(),
dest_mpirank, local_position,
1198 T_new.get_matrix());
1214template <
template <
typename>
class Matrix,
typename C>
1217 const int mpisize =
T.get_comm_size();
1218 const Shape &shape =
T.shape();
1230 const size_t local_size =
T.local_size();
1232 std::vector<size_t> local_position(local_size);
1234 T.prep_local_to_global();
1235 T_new.prep_global_to_local();
1237#pragma omp parallel default(shared)
1242 for (
size_t i = 0;
i < local_size; ++
i) {
1243 T.global_index_fast(
i, index);
1251 local_position[
i] =
pos;
1255 local_position[
i] = 0;
1262 replace_matrix_data(
T.get_matrix(),
dest_mpirank, local_position,
1263 T_new.get_matrix());
1282template <
template <
typename>
class Matrix,
typename C>
1285 const int mpisize =
T.get_comm_size();
1286 const Shape &shape =
T.shape();
1287 const size_t rank =
T.rank();
1292 for (
size_t r = 0;
r < rank; ++
r) {
1303 const size_t local_size =
T.local_size();
1305 std::vector<size_t> local_position(local_size);
1307 T.prep_local_to_global();
1308 T_new.prep_global_to_local();
1310#pragma omp parallel default(shared)
1315 for (
size_t i = 0;
i < local_size; ++
i) {
1316 T.global_index_fast(
i, index);
1319 for (
size_t r = 0;
r < rank; ++
r) {
1320 size_t idx = index[
r];
1336 local_position[
i] =
pos;
1340 local_position[
i] = 0;
1347 replace_matrix_data(
T.get_matrix(),
dest_mpirank, local_position,
1348 T_new.get_matrix());
1362template <
template <
typename>
class Matrix,
typename C>
1371 const size_t local_size =
T.local_size();
1373 std::vector<size_t> local_position(local_size);
1375 T.prep_local_to_global();
1376 T_new.prep_global_to_local();
1378#pragma omp parallel default(shared)
1383 for (
size_t i = 0;
i < local_size; ++
i) {
1384 T.global_index_fast(
i, index);
1390 local_position[
i] =
pos;
1396 replace_matrix_data(
T.get_matrix(),
dest_mpirank, local_position,
1397 T_new.get_matrix());
1408template <
template <
typename>
class Matrix,
typename C>
1411 return matrix_trace(
M.get_matrix());
1428template <
template <
typename>
class Matrix,
typename C>
1434 if (
T.rank() == 2)
return trace(
T);
1436 const size_t n =
T.local_size();
1437 const size_t l =
axes_1.size();
1443 for (
size_t i = 0;
i <
n; ++
i) {
1444 T.global_index_fast(
i, index);
1446 for (
size_t k = 0;
k <
l; ++
k) {
1457 return T.get_matrix().allreduce_sum(
sum);
1474template <
template <
typename>
class Matrix,
typename C>
1482 const size_t rank =
A.rank();
1485 Axes axes_map =
A.get_axes_map();
1489 for (
size_t i = 0;
i < rank; ++
i) {
1492 for (
size_t i = 0;
i < rank; ++
i) {
1496 const size_t n =
A.local_size();
1500 for (
size_t i = 0;
i <
n; ++
i) {
1504 return A.get_matrix().allreduce_sum(
sum);
1518template <
template <
typename>
class Matrix,
typename C>
1521 const int mpisize =
T.get_comm_size();
1524 if (
axes_1.size() == 0)
return T;
1534 size_t n =
v.size();
1535 for (
size_t i = 0;
i <
T.rank(); ++
i) {
1536 if (
k <
n &&
i ==
v[
k]) {
1550 std::vector<size_t> local_position(
T.local_size());
1556 const size_t n =
axes_1.size();
1558 for (
size_t i = 0;
i <
T.local_size(); ++
i) {
1559 T.global_index_fast(
i, index);
1561 for (
size_t k = 0;
k <
n; ++
k) {
1577 local_position[
i] =
pos;
1581 local_position[
i] = 0;
1587 sum_matrix_data(
T.get_matrix(),
dest_mpirank, local_position,
1588 T_new.get_matrix());
1611template <
template <
typename>
class Matrix,
typename C>
1622 for (
size_t i = 0;
i <
shape_b.size(); ++
i) {
1647template <
template <
typename>
class Matrix,
typename C>
1655 for (
size_t i = 0;
i <
axes_a.size(); ++
i) {
1672 const size_t rank = a.
rank();
1676 for (
size_t i = 0;
i < rank; ++
i)
v[
i] =
i;
1678 std::sort(
v,
v + rank);
1688 const size_t rank =
b.rank();
1692 for (
size_t i = 0;
i < rank; ++
i)
v[
i] =
i;
1694 std::sort(
v,
v + rank);
1721template <
template <
typename>
class Matrix,
typename C>
1742template <
template <
typename>
class Matrix,
typename C>
1762template <
template <
typename>
class Matrix,
typename C>
1764 std::vector<double> &
s) {
1769 const size_t rank = a.
rank();
1785 info = matrix_svd(
a_t.get_matrix(),
s);
1810template <
template <
typename>
class Matrix,
typename C>
1817 const size_t rank = a.
rank();
1848 info = matrix_svd(
a_t.get_matrix(),
u.get_matrix(),
s,
vt.get_matrix());
1865template <
template <
typename>
class Matrix,
typename C>
1889template <
template <
typename>
class Matrix,
typename C>
1916template <
template <
typename>
class Matrix,
typename C>
1947template <
template <
typename>
class Matrix,
typename C>
1973template <
template <
typename>
class Matrix,
typename C>
2005template <
template <
typename>
class Matrix,
typename C>
2012 const size_t rank = a.
rank();
2070template <
template <
typename>
class Matrix,
typename C>
2075 assert(shape[0] == shape[1]);
2079 size_t n = shape[0];
2084 info = matrix_eigh(
a_t.get_matrix(),
w,
z.get_matrix());
2096template <
template <
typename>
class Matrix,
typename C>
2100 assert(shape[0] == shape[1]);
2104 size_t n = shape[0];
2108 info = matrix_eigh(
a_t.get_matrix(),
w);
2134template <
template <
typename>
class Matrix,
typename C>
2140 const size_t rank = a.
rank();
2164 info = matrix_eigh(
a_t.get_matrix(),
w,
z.get_matrix());
2179template <
template <
typename>
class Matrix,
typename C>
2181 std::vector<double> &
w) {
2185 const size_t rank = a.
rank();
2202 info = matrix_eigh(
a_t.get_matrix(),
w);
2239template <
template <
typename>
class Matrix,
typename C>
2250 const size_t rank_b =
b.rank();
2280 info = matrix_eigh(
a_t.get_matrix(),
b_t.get_matrix(),
w,
z.get_matrix());
2294template <
template <
typename>
class Matrix,
typename C>
2299 assert(shape[0] == shape[1]);
2303 size_t n = shape[0];
2309 info = matrix_eig(
a_t.get_matrix(),
w,
z_t.get_matrix());
2324template <
template <
typename>
class Matrix,
typename C>
2328 assert(shape[0] == shape[1]);
2332 size_t n = shape[0];
2336 info = matrix_eig(
a_t.get_matrix(),
w);
2361template <
template <
typename>
class Matrix,
typename C>
2367 const size_t rank = a.
rank();
2392 info = matrix_eig(
a_t.get_matrix(),
w,
z_t.get_matrix());
2410template <
template <
typename>
class Matrix,
typename C>
2412 std::vector<complex> &
w) {
2416 const size_t rank = a.
rank();
2433 info = matrix_eig(
a_t.get_matrix(),
w);
2446template <
template <
typename>
class Matrix,
typename C>
2448 std::vector<C> &
x) {
2470template <
template <
typename>
class Matrix,
typename C>
2474 assert(
b.rank() == 2 ||
b.rank() == 1);
2479 if (
b.rank() == 2) {
2516template <
template <
typename>
class Matrix,
typename C>
2550 info = matrix_solve(
a_t.get_matrix(),
b_t.get_matrix());
2556template <
template <
typename>
class Matrix,
typename C>
2561template <
template <
typename>
class Matrix,
typename C>
2563 return (
rhs *= -1.0);
2566template <
template <
typename>
class Matrix,
typename C>
2572template <
template <
typename>
class Matrix,
typename C>
2578template <
template <
typename>
class Matrix,
typename C,
typename D>
2583template <
template <
typename>
class Matrix,
typename C,
typename D>
2588template <
template <
typename>
class Matrix,
typename C,
typename D>
2594template <
template <
typename>
class Matrix,
typename C>
2599template <
template <
typename>
class Matrix,
typename C>
2604template <
template <
typename>
class Matrix,
typename C>
2609template <
template <
typename>
class Matrix,
typename C>
2615template <
template <
typename>
class Matrix>
2617 return t.map(
static_cast<double (*)(
double)
>(&std::sqrt));
2619template <
template <
typename>
class Matrix>
2621 return t.map(
static_cast<complex (*)(
const complex &)
>(&std::sqrt));
2623template <
template <
typename>
class Matrix>
2627template <
template <
typename>
class Matrix>
2629 return t.map(
static_cast<complex (*)(
const complex &)
>(&std::conj));
2638template <
template <
typename>
class Matrix,
typename C>
2644 std::vector<std::size_t>
accum(
dim + 1, 1);
2647 std::vector<std::size_t> idx(
dim);
2648 for (
size_t i = 0;
i < size; ++
i) {
2650 if (idx[
dim - 1] == 0) {
2652 for (
size_t d =
dim;
d != 0 && idx[
d - 1] == 0; --
d) ++
count;
2659 out <<
' ' <<
val <<
' ';
2661 if (idx[
dim - 1] == shape[
dim - 1] - 1) {
2663 for (
size_t d =
dim;
d != 0 && idx[
d - 1] == shape[
d - 1] - 1; --
d) ++
count;
2674 for (
size_t d = 0;
d <
dim; ++
d)
out <<
'[';
2675 out <<
" distributed tensor ... ";
2676 for (
size_t d = 0;
d <
dim; ++
d)
out <<
']';
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 &)
double operator-(Timer &t1, Timer &t0)
Definition timer.hpp:64