31 #if (__GNUC__ && !__clang__)
32 #pragma GCC diagnostic ignored "-Wunused-but-set-variable"
53 template<
typename Derived1,
typename Derived2>
55 const Eigen::MatrixBase<Derived1>& state,
56 const Eigen::MatrixBase<Derived2>& A,
57 const std::vector<idx>& ctrl,
58 const std::vector<idx>& subsys,
59 const std::vector<idx>& dims)
61 const typename Eigen::MatrixBase<Derived1>::EvalReturnType& rstate
68 if (!std::is_same<
typename Derived1::Scalar,
69 typename Derived2::Scalar>::value)
86 idx d = ctrl.size() > 0 ? dims[ctrl[0]] : 1;
87 for (
idx i = 1; i < ctrl.size(); ++i)
88 if (dims[ctrl[i]] != d)
102 std::vector<idx> subsys_dims(subsys.size());
103 for (
idx i = 0; i < subsys.size(); ++i)
104 subsys_dims[i] = dims[subsys[i]];
109 std::vector<idx> ctrlgate = ctrl;
110 ctrlgate.insert(std::end(ctrlgate), std::begin(subsys), std::end(subsys));
111 std::sort(std::begin(ctrlgate), std::end(ctrlgate));
121 std::vector<dyn_mat<typename Derived1::Scalar>> Ai;
122 std::vector<dyn_mat<typename Derived1::Scalar>> Aidagger;
123 for (
idx i = 0; i < std::max(d, static_cast<idx>(2)); ++i)
125 Ai.push_back(
powm(rA, i));
129 idx D =
static_cast<idx>(rstate.rows());
131 idx ctrlsize = ctrl.size();
132 idx ctrlgatesize = ctrlgate.size();
133 idx subsyssize = subsys.size();
135 idx Dctrl =
static_cast<idx>(std::llround(std::pow(d, ctrlsize)));
136 idx DA =
static_cast<idx>(rA.rows());
144 std::vector<idx> ctrlgate_bar =
complement(ctrlgate, N);
146 idx ctrlgate_barsize = ctrlgate_bar.size();
149 for (
idx i = 0; i < ctrlgate_barsize; ++i)
150 DCTRLA_bar *= dims[ctrlgate_bar[i]];
152 for (
idx k = 0; k < N; ++k)
154 for (
idx k = 0; k < subsyssize; ++k)
155 CdimsA[k] = dims[subsys[k]];
156 for (
idx k = 0; k < ctrlsize; ++k)
158 for (
idx k = 0; k < ctrlgate_barsize; ++k)
159 CdimsCTRLA_bar[k] = dims[ctrlgate_bar[k]];
163 auto coeff_idx_ket = [&](
idx i_,
idx m_,
idx r_) noexcept
164 -> std::pair<typename Derived1::Scalar, idx>
167 typename Derived1::Scalar coeff = 0;
176 for (
idx k = 0; k < ctrlsize; ++k)
183 CdimsCTRLA_bar, CmidxCTRLA_bar);
184 for (
idx k = 0; k < N - ctrlgatesize; ++k)
186 Cmidx[ctrlgate_bar[k]] = CmidxCTRLA_bar[k];
191 for (
idx k = 0; k < subsyssize; ++k)
193 Cmidx[subsys[k]] = CmidxA[k];
200 for (
idx n_ = 0; n_ < DA; ++n_)
203 for (
idx k = 0; k < subsyssize; ++k)
205 Cmidx[subsys[k]] = CmidxA[k];
207 coeff += Ai[i_](m_, n_) *
211 return std::make_pair(coeff, indx);
217 auto coeff_idx_rho = [&](
idx i1_,
idx m1_,
219 -> std::tuple<typename Derived1::Scalar, idx, idx>
223 typename Derived1::Scalar coeff = 0, lhs = 1, rhs = 1;
238 CdimsCTRL, CmidxCTRLrow);
240 CdimsCTRL, CmidxCTRLcol);
242 for (
idx k = 0; k < ctrlsize; ++k)
244 Cmidxrow[ctrl[k]] = CmidxCTRLrow[k];
245 Cmidxcol[ctrl[k]] = CmidxCTRLcol[k];
250 CdimsCTRLA_bar, CmidxCTRLA_barrow);
252 CdimsCTRLA_bar, CmidxCTRLA_barcol);
253 for (
idx k = 0; k < N - ctrlgatesize; ++k)
255 Cmidxrow[ctrlgate_bar[k]] = CmidxCTRLA_barrow[k];
256 Cmidxcol[ctrlgate_bar[k]] = CmidxCTRLA_barcol[k];
262 for (
idx k = 0; k < subsys.size(); ++k)
264 Cmidxrow[subsys[k]] = CmidxArow[k];
265 Cmidxcol[subsys[k]] = CmidxAcol[k];
273 bool all_ctrl_rows_equal =
true;
274 bool all_ctrl_cols_equal =
true;
276 idx first_ctrl_row, first_ctrl_col;
279 first_ctrl_row = CmidxCTRLrow[0];
280 first_ctrl_col = CmidxCTRLcol[0];
283 first_ctrl_row = first_ctrl_col = 1;
286 for (
idx k = 1; k < ctrlsize; ++k)
288 if (CmidxCTRLrow[k] != first_ctrl_row)
290 all_ctrl_rows_equal =
false;
294 for (
idx k = 1; k < ctrlsize; ++k)
296 if (CmidxCTRLcol[k] != first_ctrl_col)
298 all_ctrl_cols_equal =
false;
304 for (
idx n1_ = 0; n1_ < DA; ++n1_)
307 for (
idx k = 0; k < subsyssize; ++k)
309 Cmidxrow[subsys[k]] = CmidxArow[k];
313 if (all_ctrl_rows_equal)
315 lhs = Ai[first_ctrl_row](m1_, n1_);
318 lhs = (m1_ == n1_) ? 1 : 0;
321 for (
idx n2_ = 0; n2_ < DA; ++n2_)
324 for (
idx k = 0; k < subsyssize; ++k)
326 Cmidxcol[subsys[k]] = CmidxAcol[k];
329 if (all_ctrl_cols_equal)
331 rhs = Aidagger[first_ctrl_col](n2_, m2_);
334 rhs = (n2_ == m2_) ? 1 : 0;
339 coeff += lhs * rstate(idxrowtmp, idxcoltmp) * rhs;
343 return std::make_tuple(coeff, idxrow, idxcol);
359 #pragma omp parallel for collapse(2)
360 #endif // WITH_OPENMP_
361 for (
idx m = 0; m < DA; ++m)
362 for (
idx r = 0; r < DCTRLA_bar; ++r)
366 result(coeff_idx_ket(1, m, r).second) =
367 coeff_idx_ket(1, m, r).first;
369 for (
idx i = 0; i < d; ++i)
371 result(coeff_idx_ket(i, m, r).second) =
372 coeff_idx_ket(i, m, r).first;
392 #pragma omp parallel for collapse(4)
393 #endif // WITH_OPENMP_
394 for (
idx m1 = 0; m1 < DA; ++m1)
395 for (
idx r1 = 0; r1 < DCTRLA_bar; ++r1)
396 for (
idx m2 = 0; m2 < DA; ++m2)
397 for (
idx r2 = 0; r2 < DCTRLA_bar; ++r2)
400 auto coeff_idxes = coeff_idx_rho(1, m1, r1,
402 result(std::get<1>(coeff_idxes),
403 std::get<2>(coeff_idxes)) =
404 std::get<0>(coeff_idxes);
407 for (
idx i1 = 0; i1 < Dctrl; ++i1)
408 for (
idx i2 = 0; i2 < Dctrl; ++i2)
410 auto coeff_idxes = coeff_idx_rho(
413 result(std::get<1>(coeff_idxes),
414 std::get<2>(coeff_idxes)) =
415 std::get<0>(coeff_idxes);
442 template<
typename Derived1,
typename Derived2>
444 const Eigen::MatrixBase<Derived1>& state,
445 const Eigen::MatrixBase<Derived2>& A,
446 const std::vector<idx>& ctrl,
447 const std::vector<idx>& subsys,
450 const typename Eigen::MatrixBase<Derived1>::EvalReturnType& rstate
466 std::vector<idx> dims(N, d);
468 return applyCTRL(rstate, rA, ctrl, subsys, dims);
484 template<
typename Derived1,
typename Derived2>
486 const Eigen::MatrixBase<Derived1>& state,
487 const Eigen::MatrixBase<Derived2>& A,
488 const std::vector<idx>& subsys,
489 const std::vector<idx>& dims)
491 const typename Eigen::MatrixBase<Derived1>::EvalReturnType& rstate
498 if (!std::is_same<
typename Derived1::Scalar,
499 typename Derived2::Scalar>::value)
523 std::vector<idx> subsys_dims(subsys.size());
524 for (
idx i = 0; i < subsys.size(); ++i)
525 subsys_dims[i] = dims[subsys[i]];
539 return applyCTRL(rstate, rA, {}, subsys, dims);
550 return applyCTRL(rstate, rA, {}, subsys, dims);
571 template<
typename Derived1,
typename Derived2>
573 const Eigen::MatrixBase<Derived1>& state,
574 const Eigen::MatrixBase<Derived2>& A,
575 const std::vector<idx>& subsys,
578 const typename Eigen::MatrixBase<Derived1>::EvalReturnType& rstate
594 std::vector<idx> dims(N, d);
596 return apply(rstate, rA, subsys, dims);
607 template<
typename Derived>
609 const std::vector<cmat>& Ks)
611 const cmat& rrho = rho.derived();
623 if (Ks[0].rows() != rrho.rows())
627 if (it.rows() != Ks[0].rows() || it.cols() != Ks[0].rows())
631 cmat result = cmat::Zero(rrho.rows(), rrho.rows());
634 #pragma omp parallel for
635 #endif // WITH_OPENMP_
636 for (
idx i = 0; i < Ks.size(); ++i)
640 #endif // WITH_OPENMP_
642 result += Ks[i] * rrho *
adjoint(Ks[i]);
659 template<
typename Derived>
661 const std::vector<cmat>& Ks,
662 const std::vector<idx>& subsys,
663 const std::vector<idx>& dims)
665 const cmat& rrho = rho.derived();
691 std::vector<idx> subsys_dims(subsys.size());
692 for (
idx i = 0; i < subsys.size(); ++i)
693 subsys_dims[i] = dims[subsys[i]];
704 if (it.rows() != Ks[0].rows() || it.cols() != Ks[0].rows())
708 cmat result = cmat::Zero(rrho.rows(), rrho.rows());
710 for (
idx i = 0; i < Ks.size(); ++i)
711 result +=
apply(rrho, Ks[i], subsys, dims);
726 template<
typename Derived>
728 const std::vector<cmat>& Ks,
729 const std::vector<idx>& subsys,
732 const cmat& rrho = rho.derived();
746 std::vector<idx> dims(N, d);
748 return apply(rrho, Ks, subsys, dims);
774 if (it.rows() != Ks[0].rows() || it.cols() != Ks[0].rows())
779 idx D =
static_cast<idx>(Ks[0].rows());
781 cmat result(D * D, D * D);
782 cmat MN = cmat::Zero(D, D);
783 bra A = bra::Zero(D);
784 ket B = ket::Zero(D);
785 cmat EMN = cmat::Zero(D, D);
788 #pragma omp parallel for collapse(2)
789 #endif // WITH_OPENMP_
790 for (
idx m = 0; m < D; ++m)
792 for (
idx n = 0; n < D; ++n)
796 #endif // WITH_OPENMP_
800 for (
idx i = 0; i < Ks.size(); ++i)
801 EMN += Ks[i] * MN *
adjoint(Ks[i]);
804 for (
idx a = 0; a < D; ++a)
807 for (
idx b = 0; b < D; ++b)
811 result(a * D + b, m * D + n) =
812 static_cast<cmat>(A * EMN * B).value();
817 EMN = cmat::Zero(D, D);
852 if (it.rows() != Ks[0].rows() || it.cols() != Ks[0].rows())
857 idx D =
static_cast<idx>(Ks[0].rows());
861 cmat MES = cmat::Zero(D * D, 1);
862 for (
idx a = 0; a < D; ++a)
867 cmat result = cmat::Zero(D * D, D * D);
870 #pragma omp parallel for
871 #endif // WITH_OPENMP_
872 for (
idx i = 0; i < Ks.size(); ++i)
876 #endif // WITH_OPENMP_
878 result +=
kron(cmat::Identity(D, D), Ks[i]) * Omega
910 if (D * D != static_cast<idx>(A.rows()))
916 std::vector<cmat> result;
918 for (
idx i = 0; i < D * D; ++i)
920 if (std::abs(ev(i)) >
eps)
922 std::sqrt(std::abs(ev(i))) *
reshape(evec.col(i), D, D));
946 if (D * D != static_cast<idx>(A.rows()))
950 cmat result(D * D, D * D);
953 #pragma omp parallel for collapse(4)
954 #endif // WITH_OPENMP_
955 for (
idx a = 0; a < D; ++a)
956 for (
idx b = 0; b < D; ++b)
957 for (
idx m = 0; m < D; ++m)
958 for (
idx n = 0; n < D; ++n)
959 result(a * D + b, m * D + n) = A(m * D + a, n * D + b);
982 if (D * D != static_cast<idx>(A.rows()))
986 cmat result(D * D, D * D);
989 #pragma omp parallel for collapse(4)
990 #endif // WITH_OPENMP_
991 for (
idx a = 0; a < D; ++a)
992 for (
idx b = 0; b < D; ++b)
993 for (
idx m = 0; m < D; ++m)
994 for (
idx n = 0; n < D; ++n)
995 result(m * D + a, n * D + b) = A(a * D + b, m * D + n);
1014 template<
typename Derived>
1016 const std::vector<idx>& dims)
1018 const typename Eigen::MatrixBase<Derived>::EvalReturnType& rA
1032 if (dims.size() != 2)
1050 auto worker = [=](
idx i,
idx j) noexcept
1051 ->
typename Derived::Scalar
1053 typename Derived::Scalar
sum = 0;
1054 for (
idx m = 0; m < DA; ++m)
1055 sum += rA(m * DB + i) * std::conj(rA(m * DB + j));
1061 #pragma omp parallel for collapse(2)
1062 #endif // WITH_OPENMP_
1063 for (
idx j = 0; j < DB; ++j)
1064 for (
idx i = 0; i < DB; ++i)
1065 result(i, j) = worker(i, j);
1077 auto worker = [=](
idx i,
idx j) noexcept
1078 ->
typename Derived::Scalar
1080 typename Derived::Scalar
sum = 0;
1081 for (
idx m = 0; m < DA; ++m)
1082 sum += rA(m * DB + i, m * DB + j);
1088 #pragma omp parallel for collapse(2)
1089 #endif // WITH_OPENMP_
1090 for (
idx j = 0; j < DB; ++j)
1091 for (
idx i = 0; i < DB; ++i)
1092 result(i, j) = worker(i, j);
1116 template<
typename Derived>
1118 const std::vector<idx>& dims)
1120 const typename Eigen::MatrixBase<Derived>::EvalReturnType& rA
1134 if (dims.size() != 2)
1152 auto worker = [=](
idx i,
idx j) noexcept
1153 ->
typename Derived::Scalar
1155 typename Derived::Scalar
sum = 0;
1156 for (
idx m = 0; m < DB; ++m)
1157 sum += rA(i * DB + m) * std::conj(rA(j * DB + m));
1163 #pragma omp parallel for collapse(2)
1164 #endif // WITH_OPENMP_
1165 for (
idx j = 0; j < DA; ++j)
1166 for (
idx i = 0; i < DA; ++i)
1167 result(i, j) = worker(i, j);
1180 #pragma omp parallel for collapse(2)
1181 #endif // WITH_OPENMP_
1182 for (
idx j = 0; j < DA; ++j)
1183 for (
idx i = 0; i < DA; ++i)
1184 result(i, j) =
trace(rA.block(i * DB, j * DB, DB, DB));
1208 template<
typename Derived>
1210 const std::vector<idx>& subsys,
1211 const std::vector<idx>& dims)
1213 const typename Eigen::MatrixBase<Derived>::EvalReturnType& rA
1232 idx D =
static_cast<idx>(rA.rows());
1233 idx N = dims.size();
1234 idx Nsubsys = subsys.size();
1235 idx Nsubsys_bar = N - Nsubsys;
1237 for (
idx i = 0; i < Nsubsys; ++i)
1238 Dsubsys *= dims[subsys[i]];
1239 idx Dsubsys_bar = D / Dsubsys;
1249 std::vector<idx> subsys_bar =
complement(subsys, N);
1250 std::copy(std::begin(subsys_bar), std::end(subsys_bar),
1251 std::begin(Csubsys_bar));
1253 for (
idx i = 0; i < N; ++i)
1257 for (
idx i = 0; i < Nsubsys; ++i)
1259 Csubsys[i] = subsys[i];
1260 Cdimssubsys[i] = dims[subsys[i]];
1262 for (
idx i = 0; i < Nsubsys_bar; ++i)
1264 Cdimssubsys_bar[i] = dims[subsys_bar[i]];
1278 if (subsys.size() == dims.size())
1280 result(0, 0) = (
adjoint(rA) * rA).value();
1284 if (subsys.size() == 0)
1287 auto worker = [=, &Cmidxcolsubsys_bar](
idx i) noexcept
1288 ->
typename Derived::Scalar
1299 Cdimssubsys_bar, Cmidxrowsubsys_bar);
1301 for (
idx k = 0; k < Nsubsys_bar; ++k)
1303 Cmidxrow[Csubsys_bar[k]] = Cmidxrowsubsys_bar[k];
1304 Cmidxcol[Csubsys_bar[k]] = Cmidxcolsubsys_bar[k];
1306 typename Derived::Scalar sm = 0;
1307 for (
idx a = 0; a < Dsubsys; ++a)
1312 for (
idx k = 0; k < Nsubsys; ++k)
1313 Cmidxrow[Csubsys[k]] = Cmidxcol[Csubsys[k]]
1325 for (
idx j = 0; j < Dsubsys_bar; ++j)
1329 Cdimssubsys_bar, Cmidxcolsubsys_bar);
1331 #pragma omp parallel for
1332 #endif // WITH_OPENMP_
1333 for (
idx i = 0; i < Dsubsys_bar; ++i)
1335 result(i, j) = worker(i);
1349 if (subsys.size() == dims.size())
1351 result(0, 0) = rA.trace();
1355 if (subsys.size() == 0)
1358 auto worker = [=, &Cmidxcolsubsys_bar](
idx i) noexcept
1359 ->
typename Derived::Scalar
1370 Cdimssubsys_bar, Cmidxrowsubsys_bar);
1372 for (
idx k = 0; k < Nsubsys_bar; ++k)
1374 Cmidxrow[Csubsys_bar[k]] = Cmidxrowsubsys_bar[k];
1375 Cmidxcol[Csubsys_bar[k]] = Cmidxcolsubsys_bar[k];
1377 typename Derived::Scalar sm = 0;
1378 for (
idx a = 0; a < Dsubsys; ++a)
1383 for (
idx k = 0; k < Nsubsys; ++k)
1384 Cmidxrow[Csubsys[k]] = Cmidxcol[Csubsys[k]]
1395 for (
idx j = 0; j < Dsubsys_bar; ++j)
1399 Cdimssubsys_bar, Cmidxcolsubsys_bar);
1401 #pragma omp parallel for
1402 #endif // WITH_OPENMP_
1403 for (
idx i = 0; i < Dsubsys_bar; ++i)
1405 result(i, j) = worker(i);
1431 template<
typename Derived>
1433 const std::vector<idx>& subsys,
1436 const typename Eigen::MatrixBase<Derived>::EvalReturnType& rA = A.derived();
1450 std::vector<idx> dims(N, d);
1452 return ptrace(rA, subsys, dims);
1468 template<
typename Derived>
1470 const Eigen::MatrixBase<Derived>& A,
1471 const std::vector<idx>& subsys,
1472 const std::vector<idx>& dims)
1474 const typename Eigen::MatrixBase<Derived>::EvalReturnType& rA = A.derived();
1492 idx D =
static_cast<idx>(rA.rows());
1493 idx N = dims.size();
1494 idx Nsubsys = subsys.size();
1500 for (
idx i = 0; i < N; ++i)
1502 for (
idx i = 0; i < Nsubsys; ++i)
1503 Csubsys[i] = subsys[i];
1515 if (subsys.size() == dims.size())
1518 if (subsys.size() == 0)
1521 auto worker = [=, &Cmidxcol](
idx i) noexcept
1522 ->
typename Derived::Scalar
1528 for (
idx k = 0; k < N; ++k)
1529 midxcoltmp[k] = Cmidxcol[k];
1534 for (
idx k = 0; k < Nsubsys; ++k)
1535 std::swap(midxcoltmp[Csubsys[k]], midxrow[Csubsys[k]]);
1543 for (
idx j = 0; j < D; ++j)
1549 #pragma omp parallel for
1550 #endif // WITH_OPENMP_
1551 for (
idx i = 0; i < D; ++i)
1552 result(i, j) = worker(i);
1565 if (subsys.size() == dims.size())
1566 return rA.transpose();
1568 if (subsys.size() == 0)
1571 auto worker = [=, &Cmidxcol](
idx i) noexcept
1572 ->
typename Derived::Scalar
1578 for (
idx k = 0; k < N; ++k)
1579 midxcoltmp[k] = Cmidxcol[k];
1584 for (
idx k = 0; k < Nsubsys; ++k)
1585 std::swap(midxcoltmp[Csubsys[k]], midxrow[Csubsys[k]]);
1592 for (
idx j = 0; j < D; ++j)
1598 #pragma omp parallel for
1599 #endif // WITH_OPENMP_
1600 for (
idx i = 0; i < D; ++i)
1601 result(i, j) = worker(i);
1625 template<
typename Derived>
1627 const Eigen::MatrixBase<Derived>& A,
1628 const std::vector<idx>& subsys,
1631 const typename Eigen::MatrixBase<Derived>::EvalReturnType& rA = A.derived();
1645 std::vector<idx> dims(N, d);
1662 template<
typename Derived>
1664 const Eigen::MatrixBase<Derived>& A,
1665 const std::vector<idx>& perm,
1666 const std::vector<idx>& dims)
1668 const typename Eigen::MatrixBase<Derived>::EvalReturnType& rA = A.derived();
1685 if (perm.size() != dims.size())
1690 idx D =
static_cast<idx>(rA.rows());
1691 idx N = dims.size();
1707 for (
idx i = 0; i < N; ++i)
1712 result.resize(D, 1);
1714 auto worker = [&Cdims, &Cperm, N](
idx i) noexcept
1726 for (
idx k = 0; k < N; ++k)
1728 permdims[k] = Cdims[Cperm[k]];
1729 midxtmp[k] = midx[Cperm[k]];
1736 #pragma omp parallel for
1737 #endif // WITH_OPENMP_
1738 for (
idx i = 0; i < D; ++i)
1739 result(worker(i)) = rA(i);
1755 for (
idx i = 0; i < N; ++i)
1758 Cdims[i + N] = dims[i];
1760 Cperm[i + N] = perm[i] + N;
1762 result.resize(D * D, 1);
1765 Eigen::Map<dyn_mat<typename Derived::Scalar >>(
1766 const_cast<typename Derived::Scalar*
>(rA.data()), D *
1770 auto worker = [&Cdims, &Cperm, N](
idx i) noexcept
1782 for (
idx k = 0; k < 2 * N; ++k)
1784 permdims[k] = Cdims[Cperm[k]];
1785 midxtmp[k] = midx[Cperm[k]];
1792 #pragma omp parallel for
1793 #endif // WITH_OPENMP_
1794 for (
idx i = 0; i < D * D; ++i)
1795 result(worker(i)) = rA(i);
1817 template<
typename Derived>
1819 const Eigen::MatrixBase<Derived>& A,
1820 const std::vector<idx>& perm,
1823 const typename Eigen::MatrixBase<Derived>::EvalReturnType& rA = A.derived();
1837 std::vector<idx> dims(N, d);
bool check_nonzero_size(const T &x) noexcept
Definition: util.h:112
constexpr idx maxn
Maximum number of allowed qu(d)its (subsystems)
Definition: constants.h:74
bool check_subsys_match_dims(const std::vector< idx > &subsys, const std::vector< idx > &dims)
Definition: util.h:183
bool check_cvector(const Eigen::MatrixBase< Derived > &A)
Definition: util.h:105
std::vector< cmat > choi2kraus(const cmat &A)
Orthogonal Kraus operators from Choi matrix.
Definition: operations.h:899
Eigen::MatrixXd dmat
Real (double precision) dynamic Eigen matrix.
Definition: types.h:66
constexpr double eps
Used to decide whether a number or expression in double precision is zero or not. ...
Definition: constants.h:67
Eigen::RowVectorXcd bra
Complex (double precision) dynamic Eigen row vector.
Definition: types.h:56
dyn_mat< typename Derived::Scalar > syspermute(const Eigen::MatrixBase< Derived > &A, const std::vector< idx > &perm, const std::vector< idx > &dims)
Subsystem permutation.
Definition: operations.h:1663
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic > dyn_mat
Dynamic Eigen matrix over the field specified by Scalar.
Definition: types.h:78
bool check_dims_match_mat(const std::vector< idx > &dims, const Eigen::MatrixBase< Derived > &A)
Definition: util.h:141
Eigen::VectorXcd ket
Complex (double precision) dynamic Eigen column vector.
Definition: types.h:51
dyn_mat< typename Derived1::Scalar > apply(const Eigen::MatrixBase< Derived1 > &state, const Eigen::MatrixBase< Derived2 > &A, const std::vector< idx > &subsys, const std::vector< idx > &dims)
Applies the gate A to the part subsys of the multi-partite state vector or density matrix state...
Definition: operations.h:485
Quantum++ main namespace.
Definition: codes.h:30
bool check_perm(const std::vector< idx > &perm)
Definition: util.h:242
idx get_num_subsys(idx sz, idx d)
Definition: util.h:355
dyn_col_vect< double > hevals(const Eigen::MatrixBase< Derived > &A)
Hermitian eigenvalues.
Definition: functions.h:387
bool check_dims_match_cvect(const std::vector< idx > &dims, const Eigen::MatrixBase< Derived > &V)
Definition: util.h:152
cmat choi2super(const cmat &A)
Converts Choi matrix to superoperator matrix.
Definition: operations.h:935
cmat super2choi(const cmat &A)
Converts superoperator matrix to Choi matrix.
Definition: operations.h:971
std::vector< T > complement(std::vector< T > subsys, idx N)
Constructs the complement of a subsystem vector.
Definition: functions.h:1814
dyn_mat< typename Derived::Scalar > ptrace1(const Eigen::MatrixBase< Derived > &A, const std::vector< idx > &dims)
Partial trace.
Definition: operations.h:1015
dyn_mat< typename Derived::Scalar > adjoint(const Eigen::MatrixBase< Derived > &A)
Adjoint.
Definition: functions.h:84
dyn_mat< typename T::Scalar > kron(const T &head)
Kronecker product.
Definition: functions.h:895
dyn_mat< typename Derived::Scalar > ptrace2(const Eigen::MatrixBase< Derived > &A, const std::vector< idx > &dims)
Partial trace.
Definition: operations.h:1117
Generates custom exceptions, used when validating function parameters.
Definition: exception.h:39
dyn_mat< typename Derived::Scalar > transpose(const Eigen::MatrixBase< Derived > &A)
Transpose.
Definition: functions.h:44
void n2multiidx(idx n, idx numdims, const idx *dims, idx *result) noexcept
Definition: util.h:48
dyn_mat< typename Derived::Scalar > powm(const Eigen::MatrixBase< Derived > &A, idx n)
Fast matrix power based on the SQUARE-AND-MULTIPLY algorithm.
Definition: functions.h:778
idx get_dim_subsys(idx sz, idx N)
Definition: util.h:363
cmat hevects(const Eigen::MatrixBase< Derived > &A)
Hermitian eigenvectors.
Definition: functions.h:413
cmat kraus2choi(const std::vector< cmat > &Ks)
Choi matrix.
Definition: operations.h:840
dyn_mat< typename Derived::Scalar > ptrace(const Eigen::MatrixBase< Derived > &A, const std::vector< idx > &subsys, const std::vector< idx > &dims)
Partial trace.
Definition: operations.h:1209
bool check_square_mat(const Eigen::MatrixBase< Derived > &A)
Definition: util.h:84
Derived::Scalar trace(const Eigen::MatrixBase< Derived > &A)
Trace.
Definition: functions.h:127
cmat kraus2super(const std::vector< cmat > &Ks)
Superoperator matrix.
Definition: operations.h:762
dyn_mat< typename Derived1::Scalar > applyCTRL(const Eigen::MatrixBase< Derived1 > &state, const Eigen::MatrixBase< Derived2 > &A, const std::vector< idx > &ctrl, const std::vector< idx > &subsys, const std::vector< idx > &dims)
Applies the controlled-gate A to the part subsys of the multi-partite state vector or density matrix ...
Definition: operations.h:54
dyn_mat< typename Derived::Scalar > ptranspose(const Eigen::MatrixBase< Derived > &A, const std::vector< idx > &subsys, const std::vector< idx > &dims)
Partial transpose.
Definition: operations.h:1469
Eigen::MatrixXcd cmat
Complex (double precision) dynamic Eigen matrix.
Definition: types.h:61
bool check_dims(const std::vector< idx > &dims)
Definition: util.h:125
std::size_t idx
Non-negative integer index.
Definition: types.h:36
Derived::Scalar sum(const Eigen::MatrixBase< Derived > &A)
Element-wise sum of A.
Definition: functions.h:209
dyn_mat< typename Derived::Scalar > reshape(const Eigen::MatrixBase< Derived > &A, idx rows, idx cols)
Reshape.
Definition: functions.h:1135
idx multiidx2n(const idx *midx, idx numdims, const idx *dims) noexcept
Definition: util.h:61