List of functions
Basic
Ket.ket
— Functionket([T=ComplexF64,] i::Integer, d::Integer = 2)
Produces a ket of dimension d
with nonzero element i
.
Ket.ketbra
— Functionketbra(v::AbstractVector)
Produces a ketbra of vector v
.
Ket.proj
— Functionproj([T=ComplexF64,] i::Integer, d::Integer = 2)
Produces a projector onto the basis state i
in dimension d
.
Ket.shift
— Functionshift([T=ComplexF64,] d::Integer, p::Integer = 1)
Constructs the shift operator X of dimension d
to the power p
.
Reference: Generalized Clifford algebra
Ket.clock
— Functionclock([T=ComplexF64,] d::Integer, p::Integer = 1)
Constructs the clock operator Z of dimension d
to the power p
.
Reference: Generalized Clifford algebra
Ket.pauli
— Functionpauli([T=ComplexF64,], ind::Vector{<:Integer})
Constructs the Pauli matrices: 0 or "I" for the identity, 1 or "X" for the Pauli X operation, 2 or "Y" for the Pauli Y operator, and 3 or "Z" for the Pauli Z operator. Vectors of integers between 0 and 3 or strings of I, X, Y, Z automatically generate Kronecker products of the corresponding operators.
Ket.gell_mann
— Functiongell_mann([T=ComplexF64,], d::Integer = 3)
Constructs the set G
of generalized Gell-Mann matrices in dimension d
such that G[1] = I
and G[i]*G[j] = 2 δ_ij
.
Reference: Generalizations of Pauli matrices
gell_mann([T=ComplexF64,], i::Integer, j::Integer, d::Integer = 3)
Constructs the set i
,j
th Gell-Mann matrix of dimension d
.
Reference: Generalizations of Pauli matrices
Ket.gell_mann!
— Functiongell_mann!(res::AbstractMatrix{T}, i::Integer, j::Integer, d::Integer = 3)
In-place version of gell_mann
.
Ket.partial_trace
— Functionpartial_trace(X::AbstractMatrix, remove::AbstractVector, dims::Vector)
Takes the partial trace of matrix X
with subsystem dimensions dims
over the subsystems in remove
.
partial_trace(X::AbstractMatrix, remove::Integer, dims::Vector)
Takes the partial trace of matrix X
with subsystem dimensions dims
over the subsystem remove
.
partial_trace(X::AbstractMatrix, remove::Integer)
Takes the partial trace of matrix X
over the subsystems remove
assuming two equally-sized subsystems.
Ket.partial_transpose
— Functionpartial_transpose(X::AbstractMatrix, transp::AbstractVector, dims::Vector)
Takes the partial transpose of matrix X
with subsystem dimensions dims
on the subsystems in transp
.
partial_transpose(X::AbstractMatrix, transp::Integer, dims::Vector)
Takes the partial transpose of matrix X
with subsystem dimensions dims
on the subsystem transp
.
partial_transpose(X::AbstractMatrix, transp)
Takes the partial transpose of matrix X
on the subsystems in transp
assuming two equally-sized subsystems.
Ket.permute_systems!
— Functionpermute_systems!(X::AbstractVector, perm::Vector, dims::Vector)
Permutes the order of the subsystems of vector X
with subsystem dimensions dims
in-place according to the permutation perm
.
permute_systems!(X::AbstractVector, perm::Vector)
Permutes the order of the subsystems of vector X
, which is composed by two subsystems of equal dimensions, in-place according to the permutation perm
.
Ket.permute_systems
— Functionpermute_systems(X::AbstractMatrix, perm::Vector, dims::Vector)
Permutes the order of the subsystems of the square matrix X
, which is composed by square subsystems of dimensions dims
, according to the permutation perm
.
permute_systems(X::AbstractMatrix, perm::Vector)
Permutes the order of the subsystems of the square matrix X
, which is composed by two square subsystems of equal dimensions, according to the permutation perm
.
permute_systems(X::AbstractMatrix, perm::Vector, dims::Matrix)
Permutes the order of the subsystems of the matrix X
, which is composed by subsystems of dimensions dims
, according to the permutation perm
. dims
should be a n x 2 matrix where dims[i, 1]
is the number of rows of subsystem i, and dims[i,2]
is its number of columns.
Ket.cleanup!
— Functioncleanup!(M::AbstractArray{T}; tol = Base.rtoldefault(real(T)))
Zeroes out real or imaginary parts of M
that are smaller than tol
.
Ket.symmetric_projection
— Functionsymmetric_projection(dim::Integer, n::Integer; partial::Bool=true)
Orthogonal projection onto the symmetric subspace of n
copies of a dim
-dimensional space. By default (partial=true
) it returns an isometry (say, V
) encoding the symmetric subspace. If partial=false
, then it returns the actual projection V * V'
.
Reference: Watrous' book, Sec. 7.1.1
Ket.orthonormal_range
— Functionorthonormal_range(A::AbstractMatrix{T}; mode::Integer=nothing, tol::T=nothing, sp::Bool=true) where {T<:Number}
Orthonormal basis for the range of A
. When A
is sparse (or mode = 0
), uses a QR factorization and returns a sparse result, otherwise uses an SVD and returns a dense matrix (mode = 1
). Input A
will be overwritten during the factorization. Tolerance tol
is used to compute the rank and is automatically set if not provided.
Ket.permutation_matrix
— Functionpermutation_matrix(dims::Union{Integer,Vector{<:Integer}}, perm::Vector{<:Integer})
Unitary that permutes subsystems of dimension dims
according to the permutation perm
. If dims
is an Integer, assumes there are length(perm)
subsystems of equal dimensions dims
.
Entropy
Ket.entropy
— Functionentropy([base=2,] ρ::AbstractMatrix)
Computes the von Neumann entropy -tr(ρ log ρ) of a positive semidefinite operator ρ
using a base base
logarithm.
Reference: von Neumann entropy.
entropy([base=2,] p::AbstractVector)
Computes the Shannon entropy -Σᵢpᵢlog(pᵢ) of a non-negative vector p
using a base base
logarithm.
Reference: Entropy (information theory).
Ket.binary_entropy
— Functionbinary_entropy([base=2,] p::Real)
Computes the Shannon entropy -p log(p) - (1-p)log(1-p) of a probability p
using a base base
logarithm.
Reference: Entropy (information theory).
Ket.relative_entropy
— Functionrelative_entropy([base=2,] ρ::AbstractMatrix, σ::AbstractMatrix)
Computes the (quantum) relative entropy tr(ρ
(log ρ
- log σ
)) between positive semidefinite matrices ρ
and σ
using a base base
logarithm. Note that the support of ρ
must be contained in the support of σ
but for efficiency this is not checked.
Reference: Quantum relative entropy.
relative_entropy([base=2,] p::AbstractVector, q::AbstractVector)
Computes the relative entropy D(p
||q
) = Σᵢpᵢlog(pᵢ/qᵢ) between two non-negative vectors p
and q
using a base base
logarithm. Note that the support of p
must be contained in the support of q
but for efficiency this is not checked.
Reference: Relative entropy.
Ket.binary_relative_entropy
— Functionbinary_relative_entropy([base=2,] p::Real, q::Real)
Computes the binary relative entropy D(p
||q
) = p log(p/q) + (1-p) log((1-p)/(1-q)) between two probabilities p
and q
using a base base
logarithm.
Reference: Relative entropy.
Ket.conditional_entropy
— Functionconditional_entropy([base=2,] pAB::AbstractMatrix)
Computes the conditional Shannon entropy H(A|B) of the joint probability distribution pAB
using a base base
logarithm.
Reference: Conditional entropy.
conditional_entropy([base=2,], rho::AbstractMatrix, csys::AbstractVector, dims::AbstractVector)
Computes the conditional von Neumann entropy of rho
with subsystem dimensions dims
and conditioning systems csys
, using a base base
logarithm.
Reference: Conditional quantum entropy.
Entanglement
Ket.schmidt_decomposition
— Functionschmidt_decomposition(ψ::AbstractVector, dims::AbstractVector{<:Integer})
Produces the Schmidt decomposition of ψ
with subsystem dimensions dims
. Returns the (sorted) Schmidt coefficients λ and isometries U, V such that kron(U', V')*ψ
is of Schmidt form.
Reference: Schmidt decomposition.
schmidt_decomposition(ψ::AbstractVector)
Produces the Schmidt decomposition of ψ
assuming equally-sized subsystems. Returns the (sorted) Schmidt coefficients λ and isometries U, V such that kron(U', V')*ψ
is of Schmidt form.
Reference: Schmidt decomposition.
Ket.entanglement_entropy
— Functionentanglement_entropy(ψ::AbstractVector, dims::AbstractVector{<:Integer})
Computes the relative entropy of entanglement of a bipartite pure state ψ
with subsystem dimensions dims
.
entanglement_entropy(ψ::AbstractVector)
Computes the relative entropy of entanglement of a bipartite pure state ψ
assuming equally-sized subsystems.
entanglement_entropy(ρ::AbstractMatrix, dims::AbstractVector)
Lower bounds the relative entropy of entanglement of a bipartite state ρ
with subsystem dimensions dims
.
entanglement_entropy(ρ::AbstractMatrix)
Lower bounds the relative entropy of entanglement of a bipartite state ρ
assuming equally-sized subsystems.
Ket.entanglement_dps
— Functionentanglement_dps(ρ::AbstractMatrix{T}, n::Integer = 2, d1::Union{Integer,Nothing} = nothing; sn::Integer = 1, ppt::Bool = true, verbose::Bool = false,
-solver = Hypatia.Optimizer{_solver_type(T)}, witness::Bool = true, noise::Union{AbstractMatrix,Nothing} = nothing) where {T<:Number}
Test whether ρ
has a symmetric extension with n
copies of the second subsystem. If yes, returns true
and the extension found. Otherwise, returns false
and a witness W
such that tr(W * ρ) = -1
and tr(W * σ)
is nonnegative in any symmetrically-extendable σ
.
When witness = false
, it returns the maximum visibility of ρ
such that it has a symmetric extension when mixed with noise
(default is white noise). If this visibility is < 1, then the state is entangled, otherwise the test is inconclusive. No witness is returned in this case.
For sn > 1
, it tests whether the state has a Schmidt number greater than sn
. For example, if sn = 1
(default) and the state does not have an extension, then it has Schmidt number > 1 (is entangled). Likewise, if sn = 2
and it does not have an extension, then ρ
has Schidt number at least 3. The witness
parameter works as above. No witness is returned in this case.
References: Doherty, Parrilo, Spedalieri arXiv:quant-ph/0308032 Hulpke, Bruss, Lewenstein, Sanpera arXiv:quant-ph/0401118 Weilenmann, Dive, Trillo, Aguilar, Navascués arXiv:1912.10056
Measurements
Ket.sic_povm
— Functionsic_povm([T=ComplexF64,] d::Integer)
Constructs a vector of d²
vectors |vᵢ⟩ such that |vᵢ⟩⟨vᵢ| forms a SIC-POVM of dimension d
. This construction is based on the Weyl-Heisenberg fiducial.
Reference: Appleby, Yadsan-Appleby, Zauner, arXiv:1209.1813
Ket.test_sic
— Functiontest_sic(vecs)
Tests whether vecs
is a vector of d²
vectors |vᵢ⟩ such that |vᵢ⟩⟨vᵢ| forms a SIC-POVM of dimension d
.
Ket.test_povm
— Functiontest_povm(A::Vector{<:AbstractMatrix{T}})
Checks if the measurement defined by A is valid (hermitian, semi-definite positive, and normalized).
Ket.dilate_povm
— Functiondilate_povm(vecs::Vector{Vector{T}})
Does the Naimark dilation of a rank-1 POVM given as a vector of vectors. This is the minimal dilation.
dilate_povm(E::Vector{<:AbstractMatrix})
Does the Naimark dilation of a POVM given as a vector of matrices. This always works, but is wasteful if the POVM elements are not full rank.
Ket.povm
— Functionpovm(B::Vector{<:AbstractMatrix{T}})
Creates a set of (projective) measurements from a set of bases given as unitary matrices.
povm(A::Array{T, 4}, n::Vector{Int64})
Converts a set of measurements in the common tensor format into a matrix of matrices. The second argument is fixed by the size of A
but can also contain custom number of outcomes.
Ket.mub
— Functionmub([T=ComplexF64,] d::Integer)
Construction of the standard complete set of MUBs. The output contains 1+minᵢ pᵢ^rᵢ bases, where d
= ∏ᵢ pᵢ^rᵢ.
Reference: Durt, Englert, Bengtsson, Życzkowski, arXiv:1004.3348.
Ket.test_mub
— FunctionCheck whether the input is indeed mutually unbiased
Nonlocality
Ket.chsh
— Functionchsh([T=Float64,] d::Integer = 2)
CHSH-d nonlocal game in full probability notation. If T
is an integer type the game is unnormalized.
Reference: Buhrman and Massar, arXiv:quant-ph/0409066.
Ket.cglmp
— Functioncglmp([T=Float64,] d::Integer)
CGLMP nonlocal game in full probability notation. If T
is an integer type the game is unnormalized.
References: arXiv:quant-ph/0106024 for the original game, and arXiv:2005.13418 for the form presented here.
Ket.local_bound
— Functionlocal_bound(G::Array{T,4})
Computes the local bound of a bipartite Bell functional G
, written in full probability notation as a 4-dimensional array.
Reference: Araújo, Hirsch, and Quintino, arXiv:2005.13418.
Ket.tsirelson_bound
— Functiontsirelson_bound(CG::Matrix, scenario::Vector, level::Integer)
Upper bounds the Tsirelson bound of a bipartite Bell funcional game CG
, written in Collins-Gisin notation. scenario
is vector detailing the number of inputs and outputs, in the order [oa, ob, ia, ib]. level
is an integer determining the level of the NPA hierarchy.
This function requires Moment. It is only available if you first do "import MATLAB" or "using MATLAB".
Ket.correlation_tensor
— Functioncorrelation_tensor(p::AbstractArray{T, N2}; marg::Bool = true)
Applies N sets of measurements onto a state rho
to form a probability array. Convert a 2x...x2xmx...xm probability array into
- a mx...xm correlation array (no marginals)
- a (m+1)x...x(m+1) correlation array (marginals).
Ket.probability_tensor
— Functionprobability_tensor(rho::LA.Hermitian, all_Aax::Vector{Measurement}...)
Applies N sets of measurements onto a state rho
to form a probability array.
Ket.fp2cg
— Functionfp2cg(V::Array{T,4}) where {T <: Real}
Takes a bipartite Bell functional V
in full probability notation and transforms it to Collins-Gisin notation.
Norms
Ket.trace_norm
— Functiontrace_norm(X::AbstractMatrix)
Computes trace norm of matrix X
.
Ket.kyfan_norm
— Functionkyfan_norm(X::AbstractMatrix, k::Integer, p::Real = 2)
Computes Ky-Fan (k
,p
) norm of matrix X
.
Ket.schatten_norm
— Functionschatten_norm(X::AbstractMatrix, p::Real)
Computes Schatten p
-norm of matrix X
.
Ket.diamond_norm
— Functiondiamond_norm(J::AbstractMatrix, dims::AbstractVector)
Computes the diamond norm of the supermap J
given in the Choi-Jamiołkowski representation, with subsystem dimensions dims
.
Reference: Diamond norm
diamond_norm(K::Vector{<:AbstractMatrix})
Computes the diamond norm of the CP map given by the Kraus operators K
.
Random
Ket.random_state
— Functionrandom_state([T=ComplexF64,] d::Integer, k::Integer = d)
Produces a uniformly distributed random quantum state in dimension d
with rank k
.
Reference: Życzkowski and Sommers, arXiv:quant-ph/0012101.
Ket.random_state_ket
— Functionrandom_state_ket([T=ComplexF64,] d::Integer)
Produces a Haar-random quantum state vector in dimension d
.
Reference: Życzkowski and Sommers, arXiv:quant-ph/0012101.
Ket.random_unitary
— Functionrandom_unitary([T=ComplexF64,] d::Integer)
Produces a Haar-random unitary matrix in dimension d
. If T
is a real type the output is instead a Haar-random (real) orthogonal matrix.
Reference: Stewart, doi:10.1137/0717034.
Ket.random_povm
— Functionrandom_povm([T=ComplexF64,] d::Integer, n::Integer, r::Integer)
Produces a random POVM of dimension d
with n
outcomes and rank min(k, d)
.
Reference: Heinosaari et al., arXiv:1902.04751.
Ket.random_probability
— Functionrandom_probability([T=Float64,] d::Integer)
Produces a random probability vector of dimension d
uniformly distributed on the simplex.
Reference: Dirichlet distribution
States
Ket.state_phiplus_ket
— Functionstate_phiplus_ket([T=ComplexF64,] d::Integer = 2)
Produces the vector of the maximally entangled state Φ⁺ of local dimension d
.
Ket.state_phiplus
— Functionstate_phiplus([T=ComplexF64,] d::Integer = 2; v::Real = 1)
Produces the maximally entangled state Φ⁺ of local dimension d
with visibility v
.
Ket.isotropic
— Functionisotropic(v::Real, d::Integer = 2)
Produces the isotropic state of local dimension d
with visibility v
.
Ket.state_psiminus_ket
— Functionstate_psiminus_ket([T=ComplexF64,] d::Integer = 2)
Produces the vector of the maximally entangled state ψ⁻ of local dimension d
.
Ket.state_psiminus
— Functionstate_psiminus([T=ComplexF64,] d::Integer = 2; v::Real = 1)
Produces the maximally entangled state ψ⁻ of local dimension d
with visibility v
.
Ket.state_ghz_ket
— Functionstate_ghz_ket([T=ComplexF64,] d::Integer = 2, N::Integer = 3; coeff = 1/√d)
Produces the vector of the GHZ state local dimension d
.
Ket.state_ghz
— Functionstate_ghz([T=ComplexF64,] d::Integer = 2, N::Integer = 3; v::Real = 1, coeff = 1/√d)
Produces the GHZ state of local dimension d
with visibility v
.
Ket.state_w_ket
— Functionstate_w_ket([T=ComplexF64,] N::Integer = 3; coeff = 1/√d)
Produces the vector of the N
-partite W state.
Ket.state_w
— Functionstate_w([T=ComplexF64,] N::Integer = 3; v::Real = 1, coeff = 1/√d)
Produces the N
-partite W state with visibility v
.
Ket.white_noise
— Functionwhite_noise(rho::AbstractMatrix, v::Real)
Returns v * rho + (1 - v) * id
, where id
is the maximally mixed state.
Ket.white_noise!
— Functionwhite_noise!(rho::AbstractMatrix, v::Real)
Modifies rho
in place to tranform in into v * rho + (1 - v) * id
where id
is the maximally mixed state.
Supermaps
Ket.choi
— Functionchoi(K::Vector{<:AbstractMatrix})
Constructs the Choi-Jamiołkowski representation of the CP map given by the Kraus operators K
. The convention used is that choi(K) = ∑ᵢⱼ |i⟩⟨j|⊗K|i⟩⟨j|K'
Internal functions
Ket._partition
— Functionpartition(n::Integer, k::Integer)
If n ≥ k
partitions the set 1:n
into k
parts as equally sized as possible. Otherwise partitions it into n
parts of size 1.
Ket._fiducial_WH
— Function_fiducial_WH([T=ComplexF64,] d::Integer)
Computes the fiducial Weyl-Heisenberg vector of dimension d
.
Reference: Appleby, Yadsan-Appleby, Zauner, arXiv:1209.1813 http://www.gerhardzauner.at/sicfiducials.html
Ket._idx
— Function_idx(tidx::Vector, dims::Vector)
Converts a tensor index tidx
= [i₁, i₂, ...] with subsystems dimensions dims
to a standard index.
Ket._tidx
— Function_tidx(idx::Integer, dims::Vector)
Converts a standard index idx
to a tensor index [i₁, i₂, ...] with subsystems dimensions dims
.
Ket._idxperm
— Function_idxperm(perm::Vector, dims::Vector)
Computes the index permutation associated with permuting the subsystems of a vector with subsystem dimensions dims
according to perm
.