19#include <initializer_list>
21#include <unordered_set>
25#include <boost/archive/text_iarchive.hpp>
26#include <boost/archive/text_oarchive.hpp>
27#include <boost/serialization/valarray.hpp>
28#include <boost/serialization/vector.hpp>
30#include "QubitRegisterCalculator.h"
39template <
class T = std::complex<
double>,
class Storage = std::valarray<T>>
50 return QC::QubitRegisterCalculator<>::GetNumberOfThreads();
60 Tensor(
const std::vector<size_t> &
dims,
bool dummy =
false)
63 if constexpr (std::is_convertible<int, T>::value)
73 if constexpr (std::is_convertible<int, T>::value)
80 template <
class indT =
size_t>
81 Tensor(std::initializer_list<indT>
dims,
bool dummy =
false)
84 if constexpr (std::is_convertible<int, T>::value)
95 dims.swap(other.dims);
97 std::swap(
sz, other.sz);
102 template <
class Archive>
104 if (
typename Archive::is_loading())
values.clear();
111 *
this = std::move(temp);
117 if (
this != &other) {
118 dims.swap(other.dims);
119 values.swap(other.values);
120 std::swap(
sz, other.sz);
129 std::swap(
sz, other.
sz);
134 sz = std::accumulate(
dims.begin(),
dims.end(), 1,
135 std::multiplies<size_t>());
141 assert(index <
dims.size());
153 long long pos =
dims.size() - 1;
156 if (indices[pos] <
dims[pos] - 1) {
169 long long pos =
dims.size() - 1;
172 const bool notOnPos = pos != skip;
173 if (notOnPos && indices[pos] <
dims[pos] - 1) {
177 if (notOnPos) indices[pos] = 0;
186 const std::vector<int> &skipv)
const {
187 long long int pos =
dims.size() - 1;
191 if (!skipv.empty()) {
193 skipPos =
static_cast<int>(skipv.size()) - 1;
197 const bool notOnPos = pos != skip;
198 if (notOnPos && indices[pos] <
dims[pos] - 1) {
202 if (notOnPos) indices[pos] = 0;
208 skip = skipv[skipPos];
230 void Resize(
const std::vector<size_t> &newdims) {
238 template <
class indT =
size_t>
239 void Resize(std::initializer_list<indT> newdims) {
247 const T
at(
const std::vector<size_t> &indices)
const {
248 assert(indices.size() ==
dims.size());
254 assert(indices.size() ==
dims.size());
259 inline const T &
operator[](
const std::vector<size_t> &indices)
const {
264 assert(
values.size() > offset);
270 assert(
values.size() > offset);
277 template <
class indT =
size_t>
279 assert(args.size() ==
dims.size());
281 const std::vector<size_t> indices(args.begin(), args.end());
286 template <
class indT =
size_t>
287 const T &
operator[](std::initializer_list<indT> args)
const {
288 assert(args.size() ==
dims.size());
290 const std::vector<size_t> indices(args.begin(), args.end());
296 assert(indices.size() ==
dims.size());
301 inline const T &
operator()(
const std::vector<size_t> &indices)
const {
305 template <
class indT =
size_t>
307 assert(args.size() ==
dims.size());
309 const std::vector<size_t> indices(args.begin(), args.end());
314 template <
class indT =
size_t>
315 const T &
operator()(std::initializer_list<indT> args)
const {
316 assert(args.size() ==
dims.size());
318 const std::vector<size_t> indices(args.begin(), args.end());
332 for (
size_t i = 0; i <
GetSize(); ++i)
343 for (
size_t i = 0; i <
GetSize(); ++i)
354 for (
size_t i = 0; i <
GetSize(); ++i)
365 for (
size_t i = 0; i <
GetSize(); ++i)
374 for (
size_t i = 0; i <
GetSize(); ++i)
383 for (
size_t i = 0; i <
GetSize(); ++i)
392 for (
size_t i = 0; i <
GetSize(); ++i)
401 for (
size_t i = 0; i <
GetSize(); ++i)
502 std::vector<size_t> newdims(
dims.size() + other.
dims.size());
504 for (
size_t i = 0; i <
dims.size(); ++i) newdims[i] =
dims[i];
506 for (
size_t i = 0; i < other.
dims.size(); ++i)
507 newdims[
dims.size() + i] = other.
dims[i];
512 for (
size_t i = 0; i <
GetSize(); ++i)
513 for (
size_t j = 0; j < other.
GetSize(); ++j)
521 bool allowMultithreading =
true)
const {
522 assert(
dims[ind1] == other.
dims[ind2]);
524 std::vector<size_t> newdims;
526 const size_t newsize =
dims.size() + other.
dims.size() - 2;
528 newdims.resize(1, 1);
530 newdims.reserve(newsize);
532 for (
size_t i = 0; i <
dims.size(); ++i)
533 if (i != ind1) newdims.push_back(
dims[i]);
535 for (
size_t i = 0; i < other.
dims.size(); ++i)
536 if (i != ind2) newdims.push_back(other.
dims[i]);
545 std::vector<size_t> indices1(
dims.size());
546 std::vector<size_t> indices2(other.
dims.size());
548 for (
size_t offset = 0; offset <
sz; ++offset) {
552 for (
size_t i = 0; i <
dims.size(); ++i)
554 indices1[i] = indicesres[pos];
558 for (
size_t i = 0; i < other.
dims.size(); ++i)
560 indices2[i] = indicesres[pos];
566 for (
size_t i = 0; i <
dims[ind1]; ++i) {
567 indices1[ind1] = indices2[ind2] = i;
576#pragma omp parallel for num_threads(processor_count) \
577 schedule(static, OmpLimit / divSchedule)
578 for (
long long int offset = 0; offset < static_cast<long long int>(
sz);
581 std::vector<size_t> indices1(
dims.size());
582 std::vector<size_t> indices2(other.
dims.size());
585 for (
size_t i = 0; i <
dims.size(); ++i)
587 indices1[i] = indicesres[pos];
591 for (
size_t i = 0; i < other.
dims.size(); ++i)
593 indices2[i] = indicesres[pos];
599 for (
size_t i = 0; i <
dims[ind1]; ++i) {
600 indices1[ind1] = indices2[ind2] = i;
614 const std::vector<std::pair<size_t, size_t>> &indices,
615 bool allowMultithreading =
true)
const {
616 std::vector<size_t> newdims;
617 std::vector<size_t> contractDims;
619 std::unordered_set<size_t> indicesSet1;
620 std::unordered_set<size_t> indicesSet2;
622 for (
const auto &index : indices) {
623 assert(
dims[index.first] == other.
dims[index.second]);
625 indicesSet1.insert(index.first);
626 indicesSet2.insert(index.second);
627 contractDims.push_back(
dims[index.first]);
630 const size_t newsize =
dims.size() + other.
dims.size() - 2 * indices.size();
632 newdims.resize(1, 1);
634 newdims.reserve(newsize);
636 for (
size_t i = 0; i <
dims.size(); ++i)
637 if (indicesSet1.find(i) == indicesSet1.end())
638 newdims.push_back(
dims[i]);
640 for (
size_t i = 0; i < other.
dims.size(); ++i)
641 if (indicesSet2.find(i) == indicesSet2.end())
642 newdims.push_back(other.
dims[i]);
653 std::vector<size_t> dummyIndices(
654 contractDims.size(), 0);
656 std::vector<size_t> indices1(
dims.size());
657 std::vector<size_t> indices2(other.
dims.size());
659 for (
size_t offset = 0; offset <
sz; ++offset) {
663 for (
size_t i = 0; i <
dims.size(); ++i)
664 if (indicesSet1.find(i) == indicesSet1.end()) {
665 indices1[i] = indicesres[pos];
669 for (
size_t i = 0; i < other.
dims.size(); ++i)
670 if (indicesSet2.find(i) == indicesSet2.end()) {
671 indices2[i] = indicesres[pos];
678 for (
size_t i = 0; i < dummyIndices.size(); ++i)
679 indices1[indices[i].first] = indices2[indices[i].second] =
684 const T &val2 = other[indices2];
685 auto mulRes = val1 * val2;
686 result[offset] = result[offset] + std::move(mulRes);
692#pragma omp parallel for num_threads(processor_count) \
693 schedule(static, OmpLimit / divSchedule)
694 for (
long long int offset = 0; offset < static_cast<long long int>(
sz);
697 std::vector<size_t> indices1(
dims.size());
698 std::vector<size_t> indices2(other.
dims.size());
701 for (
size_t i = 0; i <
dims.size(); ++i)
702 if (indicesSet1.find(i) == indicesSet1.end()) {
703 indices1[i] = indicesres[pos];
707 for (
size_t i = 0; i < other.
dims.size(); ++i)
708 if (indicesSet2.find(i) == indicesSet2.end()) {
709 indices2[i] = indicesres[pos];
715 std::vector<size_t> dummyIndices(
721 for (
size_t i = 0; i < dummyIndices.size(); ++i)
722 indices1[indices[i].first] = indices2[indices[i].second] =
727 const T &val2 = other[indices2];
729 auto mulRes = val1 * val2;
730 result[offset] = result[offset] + std::move(mulRes);
740 assert(
dims.size() > 0);
744 if (
values.size() == 0)
return result;
746 size_t dimMin =
dims[0];
747 for (
size_t i = 1; i <
dims.size(); ++i)
748 if (
dims[i] < dimMin) dimMin =
dims[i];
750 for (
size_t i = 0; i < dimMin; ++i) {
751 const std::vector<size_t> indices(
dims.size(), i);
760 assert(ind1 <
dims.size() && ind2 <
dims.size());
762 std::vector<size_t> newdims;
764 size_t newsize =
dims.size() - 2;
768 newdims.push_back(1);
770 newdims.reserve(newsize);
772 for (
size_t i = 0; i <
dims.size(); ++i)
773 if (i != ind1 && i != ind2) newdims.push_back(
dims[i]);
776 size_t dimMin =
dims[ind1];
777 if (
dims[ind2] < dimMin) dimMin =
dims[ind2];
782 for (
size_t offset = 0; offset < result.
GetSize(); ++offset) {
783 std::vector<size_t> indices1(
dims.size());
787 for (
size_t i = 0; i <
dims.size(); ++i)
788 if (i != ind1 && i != ind2)
789 indices1[i] = indices[i - skip];
793 for (
size_t i = 0; i < dimMin; ++i) {
794 indices1[ind1] = indices1[ind2] = i;
805 std::vector<size_t> newdims;
807 size_t newsize =
dims.size() - tindices.size();
811 newdims.push_back(1);
813 newdims.reserve(newsize);
815 for (
size_t i = 0; i <
dims.size(); ++i) {
817 for (
size_t j = 0; j < tindices.size(); ++j)
818 if (i == tindices[j]) {
823 if (!skip) newdims.push_back(
dims[i]);
827 size_t dimMin =
dims[tindices[0]];
828 for (
size_t i = 1; i < tindices.size(); ++i)
829 if (
dims[tindices[i]] < dimMin) dimMin =
dims[tindices[i]];
834 for (
size_t offset = 0; offset < result.
GetSize(); ++offset) {
835 std::vector<size_t> indices1(
dims.size());
839 for (
size_t i = 0; i <
dims.size(); ++i) {
841 for (
size_t j = 0; j < tindices.size(); ++j)
842 if (i == tindices[j]) {
848 indices1[i] = indices[i - skipv];
853 for (
size_t i = 0; i < dimMin; ++i) {
854 for (
size_t j = 0; j < tindices.size(); ++j)
855 indices1[tindices[j]] = i;
865 template <
class indT =
size_t>
867 const std::vector<size_t> indices(args.begin(), args.end());
869 return Trace(indices);
873 assert(indices.size() ==
dims.size());
875 std::vector<size_t> newdims(
dims.size());
876 for (
size_t i = 0; i <
dims.size(); ++i) newdims[i] =
dims[indices[i]];
881 std::vector<size_t> indices2(
dims.size());
882 for (
size_t offset = 0; offset <
GetSize(); ++offset) {
885 for (
size_t i = 0; i <
dims.size(); ++i)
886 indices2[indices[i]] = indices1[i];
888 result[indices2] =
values[offset];
895 template <
class indT =
size_t>
897 const std::vector<size_t> indices(args.begin(), args.end());
908 for (
size_t offset = 0; offset <
GetSize(); ++offset)
909 result[offset] =
values[offset];
914 template <
class indT =
size_t>
916 const std::vector<size_t> newdims(args.begin(), args.end());
924 inline size_t GetOffset(
const std::vector<size_t> &indices)
const {
925 return GetFortranOffset(indices);
929 return IndexFromFortranOffset(offset);
935 inline size_t GetCPPOffset(
const std::vector<size_t> &indices)
const {
936 assert(indices.size() ==
dims.size());
940 for (
size_t i = 0; i <
dims.size(); ++i) {
941 assert(indices[i] <
dims[i]);
943 result = result *
dims[i] + indices[i];
949 inline std::vector<size_t> IndexFromCPPOffset(
size_t offset)
const {
950 std::vector<size_t> indices(
dims.size(), 0);
952 for (
int i =
static_cast<int>(
dims.size()) - 1; i >= 0; --i) {
953 indices[i] = offset %
dims[i];
955 if (0 == offset)
break;
966 inline size_t GetFortranOffset(
const std::vector<size_t> &indices)
const {
967 assert(indices.size() ==
dims.size());
971 for (
long long int i =
dims.size() - 1; i >= 0; --i) {
972 assert(indices[i] <
dims[i]);
974 result = result *
dims[i] + indices[i];
980 inline std::vector<size_t> IndexFromFortranOffset(
size_t offset)
const {
981 std::vector<size_t> indices(
dims.size(), 0);
983 for (
size_t i = 0; i <
dims.size(); ++i) {
984 indices[i] = offset %
dims[i];
986 if (0 == offset)
break;
Tensor< T, Storage > & operator+=(const Tensor< T, Storage > &other)
Tensor< T, Storage > operator+(const Tensor< T, Storage > &other) const
T & operator[](size_t offset)
Tensor(const std::vector< size_t > &dims, bool dummy=false)
T & operator[](std::initializer_list< indT > args)
T atOffset(size_t offset)
Tensor< T, Storage > & operator=(const Tensor< T, Storage > &other)
std::vector< size_t > dims
Tensor< T, Storage > & operator*=(const Tensor< T, Storage > &other)
bool IncrementIndexSkip(std::vector< size_t > &indices, const std::vector< int > &skipv) const
Tensor< T, Storage > operator/(const Tensor< T, Storage > &other) const
Tensor< T, Storage > operator/(const T &scalar) const
static int GetNumberOfThreads()
bool IncrementIndexSkip(std::vector< size_t > &indices, int skip) const
void Resize(const std::vector< size_t > &newdims)
Tensor< T, Storage > & operator-=(const Tensor< T, Storage > &other)
const Storage & GetValues() const
const T & operator[](const std::vector< size_t > &indices) const
bool IncrementIndex(std::vector< size_t > &indices) const
void serialize(Archive &ar, const unsigned int)
size_t GetDim(size_t index) const
static constexpr int divSchedule
Tensor(std::initializer_list< indT > dims, bool dummy=false)
Tensor< T, Storage > & operator/=(const double scalar)
Tensor< T, Storage > Trace(const std::vector< size_t > &tindices) const
static constexpr size_t OmpLimit
void Swap(Tensor< T, Storage > &other)
const T & operator[](std::initializer_list< indT > args) const
Tensor< T, Storage > & operator/=(const Tensor< T, Storage > &other)
Tensor(const std::vector< int > &dims, bool dummy=false)
Tensor< T, Storage > & operator-=(const double scalar)
Tensor< T, Storage > & operator+=(const T &scalar)
const T at(const std::vector< size_t > &indices) const
T & operator[](const std::vector< size_t > &indices)
const size_t GetRank() const
Tensor< T, Storage > & operator+=(const double scalar)
Tensor(const Tensor< T, Storage > &other)
Tensor< T, Storage > operator*(const Tensor< T, Storage > &other) const
virtual ~Tensor()=default
Tensor(Tensor< T, Storage > &&other) noexcept
const T & operator()(std::initializer_list< indT > args) const
const std::vector< size_t > & GetDims() const
T & operator()(std::initializer_list< indT > args)
Tensor< T, Storage > & operator*=(const T &scalar)
size_t GetOffset(const std::vector< size_t > &indices) const
Tensor< T, Storage > Reshape(const std::vector< size_t > &newdims) const
Tensor< T, Storage > & operator-=(const T &scalar)
Tensor< T, Storage > operator+(const T &scalar) const
const T & operator()(const std::vector< size_t > &indices) const
Tensor< T, Storage > & operator=(Tensor< T, Storage > &&other) noexcept
Tensor< T, Storage > operator*(const T &scalar) const
const T & operator[](size_t offset) const
Tensor< T, Storage > operator-(const T &scalar) const
Tensor< T, Storage > TensorProduct(const Tensor< T, Storage > &other) const
Tensor< T, Storage > operator-() const
Tensor< T, Storage > Reshape(std::initializer_list< indT > args) const
std::vector< size_t > IndexFromOffset(size_t offset) const
friend class boost::serialization::access
Tensor< T, Storage > & operator/=(const T &scalar)
Tensor< T, Storage > Shuffle(const std::vector< size_t > &indices) const
void Resize(std::initializer_list< indT > newdims)
Tensor< T, Storage > Trace(size_t ind1, size_t ind2) const
Tensor< T, Storage > Contract(const Tensor< T, Storage > &other, const std::vector< std::pair< size_t, size_t > > &indices, bool allowMultithreading=true) const
Tensor< T, Storage > Trace(std::initializer_list< indT > args) const
Tensor< T, Storage > Contract(const Tensor< T, Storage > &other, size_t ind1, size_t ind2, bool allowMultithreading=true) const
Tensor< T, Storage > & operator*=(const double scalar)
Tensor< T, Storage > Shuffle(std::initializer_list< indT > args) const
Tensor< T, Storage > operator-(const Tensor< T, Storage > &other) const
T & operator()(const std::vector< size_t > &indices)