43 template<
typename Derived>
45 const Eigen::MatrixBase<Derived>& A)
53 return rA.transpose();
63 template<
typename Derived>
65 const Eigen::MatrixBase<Derived>& A)
73 return rA.conjugate();
83 template<
typename Derived>
105 template<
typename Derived>
126 template<
typename Derived>
127 typename Derived::Scalar
trace(
const Eigen::MatrixBase<Derived>& A)
148 template<
typename Derived>
149 typename Derived::Scalar
det(
const Eigen::MatrixBase<Derived>& A)
160 return rA.determinant();
172 template<
typename Derived>
173 typename Derived::Scalar
logdet(
const Eigen::MatrixBase<Derived>& A)
189 Eigen::PartialPivLU<dyn_mat<typename Derived::Scalar>> lu(rA);
191 lu.matrixLU().template triangularView<Eigen::Upper>();
192 typename Derived::Scalar result = std::log(U(0, 0));
194 for (
idx i = 1; i < static_cast<idx>(rA.rows()); ++i)
195 result += std::log(U(i, i));
208 template<
typename Derived>
209 typename Derived::Scalar
sum(
const Eigen::MatrixBase<Derived>& A)
230 template<
typename Derived>
231 typename Derived::Scalar
prod(
const Eigen::MatrixBase<Derived>& A)
251 template<
typename Derived>
252 double norm(
const Eigen::MatrixBase<Derived>& A)
264 return (rA.template cast<cplx>()).norm();
275 template<
typename Derived>
276 std::pair<dyn_col_vect<cplx>,
cmat>
eig(
const Eigen::MatrixBase<Derived>& A)
291 Eigen::ComplexEigenSolver<cmat> es(rA.template cast<cplx>());
293 return std::make_pair(es.eigenvalues(), es.eigenvectors());
303 template<
typename Derived>
319 return eig(rA).first;
329 template<
typename Derived>
345 Eigen::ComplexEigenSolver<cmat> es(rA.template cast<cplx>());
347 return eig(rA).second;
358 template<
typename Derived>
359 std::pair<dyn_col_vect<double>,
cmat>
heig(
const Eigen::MatrixBase<Derived>& A)
374 Eigen::SelfAdjointEigenSolver<cmat> es(rA.template cast<cplx>());
376 return std::make_pair(es.eigenvalues(), es.eigenvectors());
386 template<
typename Derived>
402 return heig(rA).first;
412 template<
typename Derived>
429 return heig(rA).second;
441 template<
typename Derived>
442 std::tuple<cmat, dyn_col_vect<double>,
cmat>
443 svd(
const Eigen::MatrixBase<Derived>& A)
454 Eigen::JacobiSVD<dyn_mat<typename Derived::Scalar>>
456 Eigen::DecompositionOptions::ComputeFullU |
457 Eigen::DecompositionOptions::ComputeFullV);
459 return std::make_tuple(sv.matrixU(), sv.singularValues(), sv.matrixV());
469 template<
typename Derived>
481 Eigen::JacobiSVD<dyn_mat<typename Derived::Scalar>> sv(rA);
483 return sv.singularValues();
493 template<
typename Derived>
505 Eigen::JacobiSVD<dyn_mat<typename Derived::Scalar>>
506 sv(rA, Eigen::DecompositionOptions::ComputeFullU);
518 template<
typename Derived>
530 Eigen::JacobiSVD<dyn_mat<typename Derived::Scalar>>
531 sv(rA, Eigen::DecompositionOptions::ComputeFullV);
545 template<
typename Derived>
561 Eigen::ComplexEigenSolver<cmat> es(rA.template cast<cplx>());
564 for (
idx i = 0; i < static_cast<idx>(evals.rows()); ++i)
567 cmat evalsdiag = evals.asDiagonal();
569 return evects * evalsdiag * evects.inverse();
578 template<
typename Derived>
594 return funm(rA, &std::sqrt);
603 template<
typename Derived>
628 template<
typename Derived>
644 return funm(rA, &std::exp);
653 template<
typename Derived>
669 return funm(rA, &std::log);
678 template<
typename Derived>
694 return funm(rA, &std::sin);
703 template<
typename Derived>
719 return funm(rA, &std::cos);
733 template<
typename Derived>
751 if (real(z) == 0 && imag(z) == 0)
752 return cmat::Identity(rA.rows(), rA.rows());
754 Eigen::ComplexEigenSolver<cmat> es(rA.template cast<cplx>());
757 for (
idx i = 0; i < static_cast<idx>(evals.rows()); ++i)
760 cmat evalsdiag = evals.asDiagonal();
762 return evects * evalsdiag * evects.inverse();
777 template<
typename Derived>
807 for (; n > 0; n /= 2)
810 result = (result * cA).eval();
811 cA = (cA * cA).eval();
825 template<
typename Derived>
826 double schatten(
const Eigen::MatrixBase<Derived>& A,
double p)
855 template<
typename OutputScalar,
typename Derived>
858 const typename Derived::Scalar&))
872 #pragma omp parallel for collapse(2)
873 #endif // WITH_OPENMP_
875 for (
idx j = 0; j < static_cast<idx>(rA.cols()); ++j)
876 for (
idx i = 0; i < static_cast<idx>(rA.rows()); ++i)
877 result(i, j) = (*f)(rA(i, j));
910 template<
typename T,
typename ... Args>
925 template<
typename Derived>
939 for (
idx i = 1; i < As.size(); ++i)
941 result =
kron(result, As[i]);
959 template<
typename Derived>
961 const std::initializer_list<Derived>& As)
963 return kron(std::vector<Derived>(As));
975 template<
typename Derived>
992 std::vector<dyn_mat<typename Derived::Scalar>> As(n, rA);
1009 template<
typename T>
1025 template<
typename T,
typename ... Args>
1040 template<
typename Derived>
1048 for (
auto&& it : As)
1053 idx total_rows = 0, total_cols = 0;
1054 for (
idx i = 0; i < As.size(); ++i)
1056 total_rows +=
static_cast<idx>(As[i].rows());
1057 total_cols +=
static_cast<idx>(As[i].cols());
1062 idx cur_row = 0, cur_col = 0;
1063 for (
idx i = 0; i < As.size(); ++i)
1065 result.block(cur_row, cur_col, As[i].rows(), As[i].cols()) = As[i];
1066 cur_row +=
static_cast<idx>(As[i].rows());
1067 cur_col +=
static_cast<idx>(As[i].cols());
1085 template<
typename Derived>
1087 const std::initializer_list<Derived>& As)
1089 return dirsum(std::vector<Derived>(As));
1101 template<
typename Derived>
1103 const Eigen::MatrixBase<Derived>& A,
idx n)
1118 std::vector<dyn_mat<typename Derived::Scalar>> As(n, rA);
1134 template<
typename Derived>
1140 idx Arows =
static_cast<idx>(rA.rows());
1141 idx Acols =
static_cast<idx>(rA.cols());
1149 if (Arows * Acols != rows * cols)
1154 return Eigen::Map<dyn_mat<typename Derived::Scalar >>(
1155 const_cast<typename Derived::Scalar*
>(rA.data()), rows, cols);
1170 template<
typename Derived1,
typename Derived2>
1172 const Eigen::MatrixBase<Derived2>& B)
1180 if (!std::is_same<
typename Derived1::Scalar,
1181 typename Derived2::Scalar>::value)
1194 if (rA.rows() != rB.rows())
1198 return rA * rB - rB * rA;
1213 template<
typename Derived1,
typename Derived2>
1215 const Eigen::MatrixBase<Derived1>& A,
1216 const Eigen::MatrixBase<Derived2>& B)
1224 if (!std::is_same<
typename Derived1::Scalar,
1225 typename Derived2::Scalar>::value)
1239 if (rA.rows() != rB.rows())
1243 return rA * rB + rB * rA;
1256 template<
typename Derived>
1272 double normV =
norm(rV);
1274 return rV *
adjoint(rV) / (normV * normV);
1277 ::Zero(rV.rows(), rV.rows());
1288 template<
typename Derived>
1297 for (
auto&& it : Vs)
1306 for (
auto&& it : Vs)
1307 if (it.rows() != Vs[0].rows() || it.cols() != 1)
1318 std::vector<dyn_mat<typename Derived::Scalar>> outvecs;
1321 for (pos = 0; pos < Vs.size(); ++pos)
1325 outvecs.push_back(Vs[pos]);
1331 for (
idx i = pos + 1; i < Vs.size(); ++i)
1333 cut -=
prj(outvecs[i - 1 - pos]);
1335 outvecs.push_back(vi);
1341 for (
auto&& it : outvecs)
1343 double normV =
norm(it);
1346 result.col(cnt) = it / normV;
1351 return result.block(0, 0, Vs[0].rows(), cnt);
1362 template<
typename Derived>
1364 const std::initializer_list<Derived>& Vs)
1366 return grams(std::vector<Derived>(Vs));
1376 template<
typename Derived>
1387 std::vector<dyn_mat<typename Derived::Scalar>> input;
1389 for (
idx i = 0; i < static_cast<idx>(rA.cols()); ++i)
1390 input.push_back(rA.col(i));
1392 return grams<dyn_mat<typename Derived::Scalar>>
1413 if (n >= std::accumulate(std::begin(dims), std::end(dims),
1414 static_cast<idx>(1), std::multiplies<idx>()))
1422 return std::vector<idx>(result, result + dims.size());
1436 const std::vector<idx>& dims)
1443 for (
idx i = 0; i < dims.size(); ++i)
1444 if (midx[i] >= dims[i])
1465 const std::vector<idx>& dims)
1467 idx N = mask.size();
1469 idx D = std::accumulate(std::begin(dims), std::end(dims),
1470 static_cast<idx>(1), std::multiplies<idx>());
1481 if (mask.size() != dims.size())
1484 for (
idx i = 0; i < N; ++i)
1485 if (mask[i] >= dims[i])
1490 ket result = ket::Zero(D);
1511 idx N = mask.size();
1512 idx D =
static_cast<idx>(std::llround(std::pow(d, N)));
1525 for (
idx i = 0; i < N; ++i)
1531 ket result = ket::Zero(D);
1532 std::vector<idx> dims(N, d);
1554 const std::vector<idx>& dims)
1556 idx N = mask.size();
1558 idx D = std::accumulate(std::begin(dims), std::end(dims),
1559 static_cast<idx>(1), std::multiplies<idx>());
1570 if (mask.size() != dims.size())
1573 for (
idx i = 0; i < N; ++i)
1574 if (mask[i] >= dims[i])
1579 cmat result = cmat::Zero(D, D);
1581 result(pos, pos) = 1;
1602 idx N = mask.size();
1603 idx D =
static_cast<idx>(std::llround(std::pow(d, N)));
1616 for (
idx i = 0; i < N; ++i)
1622 cmat result = cmat::Zero(D, D);
1623 std::vector<idx> dims(N, d);
1625 result(pos, pos) = 1;
1638 template<
typename InputIterator>
1639 std::vector<double>
abssq(InputIterator first, InputIterator last)
1641 std::vector<double> weights(std::distance(first, last));
1642 std::transform(first, last, std::begin(weights),
1643 [](
cplx z) ->
double
1657 template<
typename Container>
1658 std::vector<double>
abssq(
const Container& c,
1659 typename std::enable_if<
1669 return abssq(std::begin(c), std::end(c));
1678 template<
typename Derived>
1679 std::vector<double>
abssq(
const Eigen::MatrixBase<Derived>& A)
1690 return abssq(rA.data(), rA.data() + rA.size());
1701 template<
typename InputIterator>
1702 typename std::iterator_traits<InputIterator>::value_type
1703 sum(InputIterator first, InputIterator last)
1706 typename std::iterator_traits<InputIterator>::value_type;
1708 return std::accumulate(first, last, static_cast<value_type>(0));
1718 template<
typename Container>
1719 typename Container::value_type
1723 return sum(std::begin(c), std::end(c));
1734 template<
typename InputIterator>
1735 typename std::iterator_traits<InputIterator>::value_type
1736 prod(InputIterator first, InputIterator last)
1739 typename std::iterator_traits<InputIterator>::value_type;
1741 return std::accumulate(first, last, static_cast<value_type>(1),
1742 std::multiplies<value_type>());
1752 template<
typename Container>
1753 typename Container::value_type
1757 return prod(std::begin(c), std::end(c));
1772 template<
typename Derived>
1774 const Eigen::MatrixBase<Derived>& A)
1793 for (
idx k = 0; k < static_cast<idx>(rA.rows()); ++k)
1795 if (std::abs(tmp_evals(k)) >
eps)
1797 result = tmp_evects.col(k);
1813 template<
typename T>
1818 if (N < subsys.size())
1822 std::vector<T> all(N);
1823 std::vector<T> subsys_bar(N - subsys.size());
1825 std::iota(std::begin(all), std::end(all), 0);
1826 std::sort(std::begin(subsys), std::end(subsys));
1827 std::set_difference(std::begin(all), std::end(all),
1828 std::begin(subsys), std::end(subsys),
1829 std::begin(subsys_bar));
1844 template<
typename Derived>
1846 const Eigen::MatrixBase<Derived>& A)
1857 std::vector<double> result(3);
1858 cmat X(2, 2), Y(2, 2), Z(2, 2);
1860 Y << 0, -1_i, 1_i, 0;
1862 result[0] = std::real(
trace(rA * X));
1863 result[1] = std::real(
trace(rA * Y));
1864 result[2] = std::real(
trace(rA * Z));
1883 throw Exception(
"qpp::bloch2rho",
"r is not a 3-dimensional vector!");
1886 cmat X(2, 2), Y(2, 2), Z(2, 2), Id2(2, 2);
1888 Y << 0, -1_i, 1_i, 0;
1892 return (Id2 + r[0] * X + r[1] * Y + r[2] * Z) / 2.;
bool check_nonzero_size(const T &x) noexcept
Definition: util.h:112
std::pair< dyn_col_vect< double >, cmat > heig(const Eigen::MatrixBase< Derived > &A)
Full eigen decomposition of Hermitian expression.
Definition: functions.h:359
Derived::Scalar logdet(const Eigen::MatrixBase< Derived > &A)
Logarithm of the determinant.
Definition: functions.h:173
constexpr idx maxn
Maximum number of allowed qu(d)its (subsystems)
Definition: constants.h:74
cmat absm(const Eigen::MatrixBase< Derived > &A)
Matrix absolute value.
Definition: functions.h:604
dyn_mat< typename Derived1::Scalar > kron2(const Eigen::MatrixBase< Derived1 > &A, const Eigen::MatrixBase< Derived2 > &B)
Definition: util.h:257
bool check_cvector(const Eigen::MatrixBase< Derived > &A)
Definition: util.h:105
cmat expm(const Eigen::MatrixBase< Derived > &A)
Matrix exponential.
Definition: functions.h:629
Derived::Scalar prod(const Eigen::MatrixBase< Derived > &A)
Element-wise product of A.
Definition: functions.h:231
constexpr double eps
Used to decide whether a number or expression in double precision is zero or not. ...
Definition: constants.h:67
dyn_mat< typename T::Scalar > dirsum(const T &head)
Direct sum.
Definition: functions.h:1010
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic > dyn_mat
Dynamic Eigen matrix over the field specified by Scalar.
Definition: types.h:78
std::vector< double > rho2bloch(const Eigen::MatrixBase< Derived > &A)
Computes the 3-dimensional real Bloch vector corresponding to the qubit density matrix A...
Definition: functions.h:1845
Eigen::VectorXcd ket
Complex (double precision) dynamic Eigen column vector.
Definition: types.h:51
dyn_mat< typename Derived1::Scalar > anticomm(const Eigen::MatrixBase< Derived1 > &A, const Eigen::MatrixBase< Derived2 > &B)
Anti-commutator.
Definition: functions.h:1214
Quantum++ main namespace.
Definition: codes.h:30
std::pair< dyn_col_vect< cplx >, cmat > eig(const Eigen::MatrixBase< Derived > &A)
Full eigen decomposition.
Definition: functions.h:276
double norm(const Eigen::MatrixBase< Derived > &A)
Frobenius norm.
Definition: functions.h:252
cmat spectralpowm(const Eigen::MatrixBase< Derived > &A, const cplx z)
Matrix power.
Definition: functions.h:734
dyn_col_vect< double > hevals(const Eigen::MatrixBase< Derived > &A)
Hermitian eigenvalues.
Definition: functions.h:387
dyn_mat< OutputScalar > cwise(const Eigen::MatrixBase< Derived > &A, OutputScalar(*f)(const typename Derived::Scalar &))
Functor.
Definition: functions.h:856
Checks whether T is compatible with an STL-like iterable container.
Definition: traits.h:52
dyn_col_vect< typename Derived::Scalar > rho2pure(const Eigen::MatrixBase< Derived > &A)
Finds the pure state representation of a matrix proportional to a projector onto a pure state...
Definition: functions.h:1773
dyn_mat< typename Derived::Scalar > dirsumpow(const Eigen::MatrixBase< Derived > &A, idx n)
Direct sum power.
Definition: functions.h:1102
std::vector< T > complement(std::vector< T > subsys, idx N)
Constructs the complement of a subsystem vector.
Definition: functions.h:1814
cmat evects(const Eigen::MatrixBase< Derived > &A)
Eigenvectors.
Definition: functions.h:330
cmat funm(const Eigen::MatrixBase< Derived > &A, cplx(*f)(const cplx &))
Functional calculus f(A)
Definition: functions.h:546
double schatten(const Eigen::MatrixBase< Derived > &A, double p)
Schatten matrix norm.
Definition: functions.h:826
cmat svdU(const Eigen::MatrixBase< Derived > &A)
Left singular vectors.
Definition: functions.h:494
cmat bloch2rho(const std::vector< double > &r)
Computes the density matrix corresponding to the 3-dimensional real Bloch vector r.
Definition: functions.h:1877
std::vector< idx > n2multiidx(idx n, const std::vector< idx > &dims)
Non-negative integer index to multi-index.
Definition: functions.h:1406
dyn_mat< typename Derived::Scalar > adjoint(const Eigen::MatrixBase< Derived > &A)
Adjoint.
Definition: functions.h:84
idx multiidx2n(const std::vector< idx > &midx, const std::vector< idx > &dims)
Multi-index to non-negative integer index.
Definition: functions.h:1435
dyn_col_vect< double > svals(const Eigen::MatrixBase< Derived > &A)
Singular values.
Definition: functions.h:470
dyn_mat< typename T::Scalar > kron(const T &head)
Kronecker product.
Definition: functions.h:895
std::complex< double > cplx
Complex number in double precision.
Definition: types.h:46
cmat mprj(const std::vector< idx > &mask, const std::vector< idx > &dims)
Projector onto multi-partite qudit ket.
Definition: functions.h:1553
Generates custom exceptions, used when validating function parameters.
Definition: exception.h:39
Derived::Scalar det(const Eigen::MatrixBase< Derived > &A)
Determinant.
Definition: functions.h:149
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 > kronpow(const Eigen::MatrixBase< Derived > &A, idx n)
Kronecker power.
Definition: functions.h:976
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
dyn_mat< typename Derived::Scalar > prj(const Eigen::MatrixBase< Derived > &V)
Projector.
Definition: functions.h:1257
dyn_col_vect< cplx > evals(const Eigen::MatrixBase< Derived > &A)
Eigenvalues.
Definition: functions.h:304
cmat hevects(const Eigen::MatrixBase< Derived > &A)
Hermitian eigenvectors.
Definition: functions.h:413
dyn_mat< typename Derived::Scalar > conjugate(const Eigen::MatrixBase< Derived > &A)
Complex conjugate.
Definition: functions.h:64
constexpr double infty
Used to denote infinity in double precision.
Definition: constants.h:88
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 svdV(const Eigen::MatrixBase< Derived > &A)
Right singular vectors.
Definition: functions.h:519
dyn_mat< typename Derived1::Scalar > dirsum2(const Eigen::MatrixBase< Derived1 > &A, const Eigen::MatrixBase< Derived2 > &B)
Definition: util.h:300
cmat cosm(const Eigen::MatrixBase< Derived > &A)
Matrix cos.
Definition: functions.h:704
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > dyn_col_vect
Dynamic Eigen column vector over the field specified by Scalar.
Definition: types.h:90
cmat sinm(const Eigen::MatrixBase< Derived > &A)
Matrix sin.
Definition: functions.h:679
std::vector< double > abssq(InputIterator first, InputIterator last)
Computes the absolute values squared of an STL-like range of complex numbers.
Definition: functions.h:1639
dyn_mat< typename Derived1::Scalar > comm(const Eigen::MatrixBase< Derived1 > &A, const Eigen::MatrixBase< Derived2 > &B)
Commutator.
Definition: functions.h:1171
Eigen::MatrixXcd cmat
Complex (double precision) dynamic Eigen matrix.
Definition: types.h:61
ket mket(const std::vector< idx > &mask, const std::vector< idx > &dims)
Multi-partite qudit ket.
Definition: functions.h:1464
dyn_mat< typename Derived::Scalar > inverse(const Eigen::MatrixBase< Derived > &A)
Inverse.
Definition: functions.h:106
cmat logm(const Eigen::MatrixBase< Derived > &A)
Matrix logarithm.
Definition: functions.h:654
bool check_dims(const std::vector< idx > &dims)
Definition: util.h:125
cmat sqrtm(const Eigen::MatrixBase< Derived > &A)
Matrix square root.
Definition: functions.h:579
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
std::tuple< cmat, dyn_col_vect< double >, cmat > svd(const Eigen::MatrixBase< Derived > &A)
Full singular value decomposition.
Definition: functions.h:443
dyn_mat< typename Derived::Scalar > grams(const std::vector< Derived > &Vs)
Gram-Schmidt orthogonalization.
Definition: functions.h:1289
bool check_qubit_matrix(const Eigen::MatrixBase< Derived > &A) noexcept
Definition: util.h:213
idx multiidx2n(const idx *midx, idx numdims, const idx *dims) noexcept
Definition: util.h:61