42 template <
typename Derived>
43 double entropy(
const Eigen::MatrixBase<Derived>& A) {
59 for (
idx i = 0; i < static_cast<idx>(ev.rows()); ++i)
61 result -= ev(i) * std::log2(ev(i));
72 inline double entropy(
const std::vector<double>& prob) {
81 for (
idx i = 0; i < prob.size(); ++i)
82 if (std::abs(prob[i]) != 0)
83 result -= std::abs(prob[i]) * std::log2(std::abs(prob[i]));
100 template <
typename Derived>
101 double renyi(
const Eigen::MatrixBase<Derived>& A,
double alpha) {
119 return std::log2(rA.rows());
125 return -std::log2(svals(rA)[0]);
129 for (
idx i = 0; i < static_cast<idx>(sv.rows()); ++i)
131 result += std::pow(sv(i), alpha);
133 return std::log2(result) / (1 - alpha);
148 inline double renyi(
const std::vector<double>& prob,
double alpha) {
159 return std::log2(prob.size());
168 for (
idx i = 0; i < prob.size(); ++i)
169 if (std::abs(prob[i]) > max)
170 max = std::abs(prob[i]);
172 return -std::log2(max);
176 for (
idx i = 0; i < prob.size(); ++i)
177 if (std::abs(prob[i]) != 0)
178 result += std::pow(std::abs(prob[i]), alpha);
180 return std::log2(result) / (1 - alpha);
195 template <
typename Derived>
196 double tsallis(
const Eigen::MatrixBase<Derived>& A,
double q) {
214 return entropy(rA) * std::log(2);
218 for (
idx i = 0; i < static_cast<idx>(ev.rows()); ++i)
220 result += std::pow(ev(i), q);
222 return (result - 1) / (1 - q);
237 inline double tsallis(
const std::vector<double>& prob,
double q) {
249 return entropy(prob) * std::log(2);
252 for (
idx i = 0; i < prob.size(); ++i)
253 if (std::abs(prob[i]) != 0)
254 result += std::pow(std::abs(prob[i]), q);
256 return (result - 1) / (1 - q);
268 template <
typename Derived>
270 const std::vector<idx>& subsysA,
271 const std::vector<idx>& subsysB,
272 const std::vector<idx>& dims) {
300 std::vector<idx> full_system(dims.size());
301 std::iota(std::begin(full_system), std::end(full_system), 0);
304 std::vector<idx> subsysAsorted{subsysA};
305 std::vector<idx> subsysBsorted{subsysB};
308 std::sort(std::begin(subsysAsorted), std::end(subsysAsorted));
309 std::sort(std::begin(subsysBsorted), std::end(subsysBsorted));
312 std::vector<idx> subsysAB;
313 std::set_union(std::begin(subsysAsorted), std::end(subsysAsorted),
314 std::begin(subsysBsorted), std::end(subsysBsorted),
315 std::back_inserter(subsysAB));
318 std::vector<idx> subsysA_bar = complement(subsysA, dims.size());
319 std::vector<idx> subsysB_bar = complement(subsysB, dims.size());
321 std::vector<idx> subsysAB_bar = complement(subsysAB, dims.size());
339 template <
typename Derived>
341 const std::vector<idx>& subsysA,
342 const std::vector<idx>& subsysB,
idx d = 2) {
357 std::vector<idx> dims(N, d);
bool check_nonzero_size(const T &x) noexcept
Definition: util.h:123
Dimension(s) mismatch matrix size exception.
Definition: exception.h:298
double qmutualinfo(const Eigen::MatrixBase< Derived > &A, const std::vector< idx > &subsysA, const std::vector< idx > &subsysB, const std::vector< idx > &dims)
Quantum mutual information between 2 subsystems of a composite system.
Definition: entropies.h:269
bool check_subsys_match_dims(const std::vector< idx > &subsys, const std::vector< idx > &dims)
Definition: util.h:210
double tsallis(const Eigen::MatrixBase< Derived > &A, double q)
Tsallis- entropy of the density matrix A, for .
Definition: entropies.h:196
Eigen::MatrixXd dmat
Real (double precision) dynamic Eigen matrix.
Definition: types.h:69
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic > dyn_mat
Dynamic Eigen matrix over the field specified by Scalar.
Definition: types.h:81
bool check_dims_match_mat(const std::vector< idx > &dims, const Eigen::MatrixBase< Derived > &A)
Definition: util.h:150
Subsystems mismatch dimensions exception.
Definition: exception.h:363
Quantum++ main namespace.
Definition: codes.h:35
idx get_num_subsys(idx sz, idx d)
Definition: util.h:366
Invalid dimension(s) exception.
Definition: exception.h:267
double entropy(const Eigen::MatrixBase< Derived > &A)
von-Neumann entropy of the density matrix A
Definition: entropies.h:43
double renyi(const Eigen::MatrixBase< Derived > &A, double alpha)
Renyi- entropy of the density matrix A, for .
Definition: entropies.h:101
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:1176
constexpr double infty
Used to denote infinity in double precision.
Definition: constants.h:87
bool check_square_mat(const Eigen::MatrixBase< Derived > &A)
Definition: util.h:99
Parameter out of range exception.
Definition: exception.h:513
Eigen::MatrixXcd cmat
Complex (double precision) dynamic Eigen matrix.
Definition: types.h:64
bool check_dims(const std::vector< idx > &dims)
Definition: util.h:134
std::size_t idx
Non-negative integer index.
Definition: types.h:39
Matrix is not square exception.
Definition: exception.h:147
Object has zero size exception.
Definition: exception.h:132