10#include <Eigen/SparseCore>
14#include <oneapi/tbb.h>
18template <
typename Scalar>
19Eigen::SparseMatrix<Scalar, Eigen::RowMajor>
22 const Eigen::SparseMatrix<Scalar, Eigen::RowMajor> &matrix1,
23 const Eigen::SparseMatrix<Scalar, Eigen::RowMajor> &matrix2) {
24 oneapi::tbb::concurrent_vector<Eigen::Triplet<Scalar>> triplets;
27 oneapi::tbb::parallel_for(
28 oneapi::tbb::blocked_range<Eigen::Index>(0, matrix1.outerSize()), [&](
const auto &range) {
29 for (Eigen::Index row1 = range.begin(); row1 != range.end(); ++row1) {
31 const auto &range_row2 = basis_final->get_index_range(row1);
34 for (auto row2 = static_cast<Eigen::Index>(range_row2.min());
35 row2 < static_cast<Eigen::Index>(range_row2.max()); ++row2) {
37 Eigen::Index row = basis_final->get_ket_index_from_tuple(row1, row2);
43 for (typename Eigen::SparseMatrix<Scalar, Eigen::RowMajor>::InnerIterator it1(
47 Eigen::Index col1 = it1.col();
48 Scalar value1 = it1.value();
52 Eigen::Index begin_idxptr2 = matrix2.outerIndexPtr()[row2];
53 Eigen::Index end_idxptr2 = matrix2.outerIndexPtr()[row2 + 1];
57 const auto &range_col2 = basis_initial->get_index_range(it1.index());
59 std::distance(matrix2.innerIndexPtr() + begin_idxptr2,
60 std::lower_bound(matrix2.innerIndexPtr() + begin_idxptr2,
61 matrix2.innerIndexPtr() + end_idxptr2,
67 for (Eigen::Index idxptr2 = begin_idxptr2; idxptr2 < end_idxptr2;
70 Eigen::Index col2 = matrix2.innerIndexPtr()[idxptr2];
71 if (col2 >= static_cast<Eigen::Index>(range_col2.max())) {
75 Eigen::Index col = basis_initial->get_ket_index_from_tuple(col1, col2);
80 Scalar value2 = matrix2.valuePtr()[idxptr2];
83 Scalar value = value1 * value2;
84 triplets.emplace_back(row, col, value);
92 Eigen::SparseMatrix<Scalar, Eigen::RowMajor> matrix(basis_final->get_number_of_kets(),
93 basis_initial->get_number_of_kets());
94 matrix.setFromTriplets(triplets.begin(), triplets.end());
95 matrix.makeCompressed();
101template Eigen::SparseMatrix<double, Eigen::RowMajor>
103 const std::shared_ptr<
const BasisPair<double>> &,
104 const Eigen::SparseMatrix<double, Eigen::RowMajor> &,
105 const Eigen::SparseMatrix<double, Eigen::RowMajor> &);
106template Eigen::SparseMatrix<std::complex<double>, Eigen::RowMajor>
108 const std::shared_ptr<
const BasisPair<std::complex<double>>> &,
109 const Eigen::SparseMatrix<std::complex<double>, Eigen::RowMajor> &,
110 const Eigen::SparseMatrix<std::complex<double>, Eigen::RowMajor> &);
Eigen::SparseMatrix< Scalar, Eigen::RowMajor > calculate_tensor_product(const std::shared_ptr< const BasisPair< Scalar > > &basis_initial, const std::shared_ptr< const BasisPair< Scalar > > &basis_final, const Eigen::SparseMatrix< Scalar, Eigen::RowMajor > &matrix1, const Eigen::SparseMatrix< Scalar, Eigen::RowMajor > &matrix2)