42 template<
typename Derived>
43 double entropy(
const Eigen::MatrixBase<Derived>& A)
60 for (
idx i = 0; i < static_cast<idx>(ev.rows()); ++i)
62 result -= ev(i) * std::log2(ev(i));
73 inline double entropy(
const std::vector<double>& prob)
83 for (
idx i = 0; i < prob.size(); ++i)
84 if (std::abs(prob[i]) != 0)
85 result -= std::abs(prob[i]) * std::log2(std::abs(prob[i]));
102 template<
typename Derived>
103 double renyi(
const Eigen::MatrixBase<Derived>& A,
double alpha)
122 return std::log2(rA.rows());
128 return -std::log2(
svals(rA)[0]);
132 for (
idx i = 0; i < static_cast<idx>(sv.rows()); ++i)
134 result += std::pow(sv(i), alpha);
136 return std::log2(result) / (1 - alpha);
151 inline double renyi(
const std::vector<double>& prob,
double alpha)
163 return std::log2(prob.size());
172 for (
idx i = 0; i < prob.size(); ++i)
173 if (std::abs(prob[i]) > max)
174 max = std::abs(prob[i]);
176 return -std::log2(max);
180 for (
idx i = 0; i < prob.size(); ++i)
181 if (std::abs(prob[i]) != 0)
182 result += std::pow(std::abs(prob[i]), alpha);
184 return std::log2(result) / (1 - alpha);
199 template<
typename Derived>
200 double tsallis(
const Eigen::MatrixBase<Derived>& A,
double q)
219 return entropy(rA) * std::log(2.);
223 for (
idx i = 0; i < static_cast<idx>(ev.rows()); ++i)
225 result += std::pow(ev(i), q);
227 return (result - 1) / (1 - q);
242 inline double tsallis(
const std::vector<double>& prob,
double q)
255 return entropy(prob) * std::log(2.);
258 for (
idx i = 0; i < prob.size(); ++i)
259 if (std::abs(prob[i]) != 0)
260 result += std::pow(std::abs(prob[i]), q);
262 return (result - 1) / (1 - q);
274 template<
typename Derived>
276 const std::vector<idx>& subsysA,
277 const std::vector<idx>& subsysB,
278 const std::vector<idx>& dims)
310 std::vector<idx> full_system(dims.size());
311 std::iota(std::begin(full_system), std::end(full_system), 0);
314 std::vector<idx> subsysAsorted{subsysA};
315 std::vector<idx> subsysBsorted{subsysB};
318 std::sort(std::begin(subsysAsorted), std::end(subsysAsorted));
319 std::sort(std::begin(subsysBsorted), std::end(subsysBsorted));
322 std::vector<idx> subsysAB;
323 std::set_union(std::begin(subsysAsorted), std::end(subsysAsorted),
324 std::begin(subsysBsorted), std::end(subsysBsorted),
325 std::back_inserter(subsysAB));
328 std::vector<idx> subsysA_bar =
complement(subsysA, dims.size());
329 std::vector<idx> subsysB_bar =
complement(subsysB, dims.size());;
330 std::vector<idx> subsysAB_bar =
complement(subsysAB, dims.size());
348 template<
typename Derived>
350 const std::vector<idx>& subsysA,
351 const std::vector<idx>& subsysB,
368 std::vector<idx> dims(N, d);
bool check_nonzero_size(const T &x) noexcept
Definition: util.h:127
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:275
bool check_subsys_match_dims(const std::vector< idx > &subsys, const std::vector< idx > &dims)
Definition: util.h:198
double tsallis(const Eigen::MatrixBase< Derived > &A, double q)
Tsallis- entropy of the density matrix A, for .
Definition: entropies.h:200
Eigen::MatrixXd dmat
Real (double precision) dynamic Eigen matrix.
Definition: types.h:66
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:156
Quantum++ main namespace.
Definition: codes.h:30
idx get_num_subsys(idx sz, idx d)
Definition: util.h:370
std::vector< T > complement(std::vector< T > subsys, idx N)
Constructs the complement of a subsystem vector.
Definition: functions.h:1814
dyn_col_vect< double > svals(const Eigen::MatrixBase< Derived > &A)
Singular values.
Definition: functions.h:470
double entropy(const Eigen::MatrixBase< Derived > &A)
von-Neumann entropy of the density matrix A
Definition: entropies.h:43
Generates custom exceptions, used when validating function parameters.
Definition: exception.h:39
double renyi(const Eigen::MatrixBase< Derived > &A, double alpha)
Renyi- entropy of the density matrix A, for .
Definition: entropies.h:103
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:1280
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:99
Eigen::MatrixXcd cmat
Complex (double precision) dynamic Eigen matrix.
Definition: types.h:61
bool check_dims(const std::vector< idx > &dims)
Definition: util.h:140
std::size_t idx
Non-negative integer index.
Definition: types.h:36