Skip to content

Commit

Permalink
Added TraceDistanceCoherence and RelEntCoherence
Browse files Browse the repository at this point in the history
  • Loading branch information
nathanieljohnston committed Jan 12, 2016
1 parent 8f4e3c7 commit 60d0d14
Show file tree
Hide file tree
Showing 4 changed files with 150 additions and 3 deletions.
5 changes: 4 additions & 1 deletion CHANGELOG.md
Original file line number Diff line number Diff line change
@@ -1,7 +1,7 @@
# Change Log
All notable changes to QETLAB will be documented in this file.

## Changes since 2015-04-13
## [0.9] - 2016-01-12
### Added
- BCSGameLB: Computes a lower bound on the quantum value of a binary constraint system (BCS) game.
- BCSGameValue: Computes the maximum value of a binary constraint system (BCS) game. In classical and no-signalling settings, the value computed is exact, but the quantum value is just an upper bound.
Expand All @@ -10,11 +10,14 @@ All notable changes to QETLAB will be documented in this file.
- L1NormCoherence: Computes the l1-norm of coherence of a quantum state.
- NonlocalGameLB: Computes a lower bound on the quantum value of a two-player non-local game.
- RandomPOVM: Computes a random POVM of a specified size and with a specified number of outcomes.
- RelEntCoherence: Computes the relative entropy of coherence of a quantum state.
- RobustnessCoherence: Computes the robustness of coherence of a quantum state.
- TraceDistanceCoherence: Computes the trace distance of coherence of a quantum state.
- helpers/bcs_to_nonlocal: Converts a description of a binary constraint system (BCS) game into a form that can be presented as a general non-local game.
- helpers/pure_to_mixed: Converts a state vector or density matrix representation of a state to a density matrix.

### Changed
- Entropy: Fixed a bug that would cause NaN output for some low-rank input states.
- kpNorm: Can now be used as the objective function or as a constraint in a CVX optimization problem, regardless of k and p (only certain special values of k and p were supported previously).
- kpNormDual: Can now be used as the objective function or as a constraint in a CVX optimization problem, regardless of k and p (only certain special values of k and p were supported previously).
- NonlocalGameValue: Added the REPT optional input argument, which lets the user specify the number of times that the non-local game will be repeated in parallel.
Expand Down
5 changes: 3 additions & 2 deletions Entropy.m
Original file line number Diff line number Diff line change
Expand Up @@ -17,7 +17,7 @@

% requires: nothing
% author: Nathaniel Johnston ([email protected])
% last updated: November 27, 2014
% last updated: January 12, 2016

function ent = Entropy(rho,varargin)

Expand All @@ -28,13 +28,14 @@

% If alpha == 1, compute the von Neumann entropy
if(abs(alpha - 1) <= eps^(3/4))
lam(lam==0) = 1; % handle zero entries better: we want 0*log(0) = 0, not NaN
if(base == 2)
ent = -sum(real(lam.*log2(lam)));
else
ent = -sum(real(lam.*log(lam)))/log(base);
end
elseif(alpha >= 0)
if(alpha < Inf) % Renty-alpha entropy with ALPHA < Inf
if(alpha < Inf) % Renyi-alpha entropy with ALPHA < Inf
ent = log(sum(lam.^alpha))/(log(base)*(1-alpha));

% Check whether or not we ran into numerical problems due to ALPHA
Expand Down
19 changes: 19 additions & 0 deletions RelEntCoherence.m
Original file line number Diff line number Diff line change
@@ -0,0 +1,19 @@
%% RelEntCoherence Computes the relative entropy of coherence of a quantum state
% This function has one required argument:
% RHO: a density matrix or pure state vector
%
% REC = RelEntCoherence(RHO) is the relative entropy of coherence, which
% equals S(diag(RHO)) - S(RHO), where S(.) is the von Neumann entropy.
%
% URL: http://www.qetlab.com/RelEntCoherence

% requires: Entropy.m, pure_to_mixed.m
% author: Nathaniel Johnston ([email protected])
% package: QETLAB
% last updated: January 12, 2016

function REC = RelEntCoherence(rho)

rho = pure_to_mixed(rho); % Let the user specify rho as either a density matrix or a pure state vector

REC = Entropy(diag(diag(rho))) - Entropy(rho);
124 changes: 124 additions & 0 deletions TraceDistanceCoherence.m
Original file line number Diff line number Diff line change
@@ -0,0 +1,124 @@
%% TraceDistanceCoherence Computes the trace distance of coherence of a quantum state
% This function has one required argument:
% RHO: a pure state vector or a density matrix
%
% [TDC,D] = TraceDistanceCoherence(RHO) is the trace distance of
% coherence of the quantum state (density matrix) RHO (i.e., it is the
% minimal distance in the trace norm between RHO and a diagonal density
% matrix). The optional output argument D is a vector containing the
% entries of the diagonal state that is closest to RHO
% (in other words, TraceNorm(RHO - diag(D)) == TDC).
%
% URL: http://www.qetlab.com/TraceDistanceCoherence

% requires: cvx (http://cvxr.com/cvx/)
% author: Nathaniel Johnston ([email protected])
% package: QETLAB
% last updated: January 12, 2016

function [TDC,D] = TraceDistanceCoherence(rho)

% Compute the size of RHO. If it's a pure state (vector), we can
% compute it quickly. If not, we have to use a semidefinite program.
[m,n] = size(rho);
ispure = (min(m,n) == 1);

% If the state is single-qubit, we can compute it faster.
if(max(m,n) == 1)
TDC = 0;
D = 1;
return
elseif(m == 2 && n == 2)
TDC = 2*abs(rho(1,2));
D = diag(diag(rho));
return
else
if(m == n && rank(rho) == 1)
[~,~,v] = svd(rho); % v(:,1) is the pure state form of rho
rho = v(:,1);
ispure = 1;
end

% Fore pure states, we can compute it much faster than via SDP.
% The algorithm is somewhat complicated though, so we contain it in
% the auxiliary function TraceCoherencePure.
if(ispure)
[TDC,D] = TraceCoherencePure(rho);
return
end
end

% In general, trace distance of coherence is computed by semidefinite
% programming.
cvx_begin sdp quiet
cvx_precision best
variable d(n)

minimize norm_nuc(rho - diag(d))

subject to
d >= 0;
sum(d) == 1;
cvx_end

D = d/sum(d); % divide by sum(d) for numerical reasons
TDC = real(cvx_optval);
end

% The following function computes the trace distance coherence of a pure
% state using a method that is much faster than the naive SDP.
function [val,D] = TraceCoherencePure(x)
x = abs(x(:));
[xs,ind] = sort(x,'descend'); % sort the entries of x in decending order

n = length(xs);

% Take care of a degenerate case.
if(n == 1)
val = 0;
D = 1;
return
end

% Now loop through to try to find k. Use binary search to make this as
% fast as possible.
klow = 1;
khigh = n;
while 1==1 % loop until we explicitly break out
k = ceil((khigh+klow)/2);
s = sum(xs(1:k));
m = sum(xs((k+1):n).^2);
r = (s^2-1-k*m);
q = (r + sqrt(r^2 + 4*k*m*s^2))/(2*k*s);

if(q < xs(k))
% We have found a potential k!
% Now compute D.
D = zeros(n,1);
for j = 1:k
D(j) = (xs(j) - q)/(s - k*q);
end

D = D(perm_inv(ind));

% We could compute val = sum(svd(x*x' - D)), but the following
% method is much quicker.
val = 2*q/(s - k*q);

% We want to increase klow (i.e., look at the top half of
% our range).
klow = k;
end

if(k == khigh)
% We're done!
return
end

if(q >= xs(k))
% We want to decrease khigh (i.e., look at the bottom half of
% our range).
khigh = k;
end
end
end

0 comments on commit 60d0d14

Please sign in to comment.