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.gellmann
— Functiongellmann([T=ComplexF64,], d::Integer = 3)
Constructs the set G
of generalized Gell-Mann matrices in dimension d
such that G[1] = I
and Tr(G[i]*G[j]) = 2 δ_ij
.
Reference: Generalizations of Pauli matrices
gellmann([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.gellmann!
— Functiongellmann!(res::AbstractMatrix{T}, i::Integer, j::Integer, d::Integer = 3)
In-place version of gellmann
.
Ket.partial_trace
— Functionpartial_trace(X::AbstractMatrix, remove::AbstractVector, dims::AbstractVector = _equal_sizes(X))
Takes the partial trace of matrix X
with subsystem dimensions dims
over the subsystems in remove
. If the argument dims
is omitted two equally-sized subsystems are assumed.
partial_trace(X::AbstractMatrix, remove::Integer, dims::AbstractVector = _equal_sizes(X)))
Takes the partial trace of matrix X
with subsystem dimensions dims
over the subsystem remove
. If the argument dims
is omitted two equally-sized subsystems are assumed.
Ket.partial_transpose
— Functionpartial_transpose(X::AbstractMatrix, transp::AbstractVector, dims::AbstractVector = _equal_sizes(X))
Takes the partial transpose of matrix X
with subsystem dimensions dims
on the subsystems in transp
. If the argument dims
is omitted two equally-sized subsystems are assumed.
partial_transpose(X::AbstractMatrix, transp::Integer, dims::AbstractVector = _equal_sizes(X))
Takes the partial transpose of matrix X
with subsystem dimensions dims
on the subsystem transp
. If the argument dims
is omitted two equally-sized subsystems are assumed.
Ket.permute_systems!
— Functionpermute_systems!(X::AbstractVector, perm::AbstractVector, dims::AbstractVector = _equal_sizes(X))
Permutes the order of the subsystems of vector X
with subsystem dimensions dims
in-place according to the permutation perm
. If the argument dims
is omitted two equally-sized subsystems are assumed.
Ket.permute_systems
— Functionpermute_systems(X::AbstractMatrix, perm::AbstractVector, dims::AbstractVector = _equal_sizes(X))
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
. If the argument dims
is omitted two equally-sized subsystems are assumed.
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,AbstractVector}, perm::AbstractVector)
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
.
Ket.n_body_basis
— Functionn_body_basis(
n::Integer,
n_parties::Integer;
sb::AbstractVector{<:AbstractMatrix} = [pauli(1), pauli(2), pauli(3)],
sparse::Bool = true,
eye::AbstractMatrix = I(size(sb[1], 1))
Return the basis of n
nontrivial operators acting on n_parties
, by default using Pauli matrices.
For example, n_body_basis(2, 3)
generates all products of two Paulis and one identity, so ${X ⊗ X ⊗ 1, X ⊗ 1 ⊗ X, ..., X ⊗ Y ⊗ 1, ..., 1 ⊗ Z ⊗ Z}$.
Instead of Paulis, a basis can be provided by the parameter sb
, and the identity can be changed with eye
. If sparse
is true, the resulting basis will use sparse matrices, otherwise it will agree with sb
.
This function returns a generator, which can then be used e.g. in for loops without fully allocating the entire basis at once. If you need a vector, call collect
on it.
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} = _equal_sizes(ψ))
Produces the Schmidt decomposition of ψ
with subsystem dimensions dims
. If the argument dims
is omitted equally-sized subsystems are assumed. 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} = _equal_sizes(ψ))
Computes the relative entropy of entanglement of a bipartite pure state ψ
with subsystem dimensions dims
. If the argument dims
is omitted equally-sized subsystems are assumed.
entanglement_entropy(ρ::AbstractMatrix, dims::AbstractVector = _equal_sizes(ρ), n::Integer = 1)
Lower bounds the relative entropy of entanglement of a bipartite state ρ
with subsystem dimensions dims
using level n
of the DPS hierarchy. If the argument dims
is omitted equally-sized subsystems are assumed.
Ket.random_robustness
— Functionrandom_robustness(
ρ::AbstractMatrix{T},
dims::AbstractVector{<:Integer} = _equal_sizes(ρ),
n::Integer = 1;
ppt::Bool = true,
verbose::Bool = false,
solver = Hypatia.Optimizer{_solver_type(T)})
Lower bounds the random robustness of state ρ
with subsystem dimensions dims
using level n
of the DPS hierarchy. Argument ppt
indicates whether to include the partial transposition constraints.
Ket.schmidt_number
— Functionschmidt_number(
ρ::AbstractMatrix{T},
s::Integer = 2,
dims::AbstractVector{<:Integer} = _equal_sizes(ρ),
n::Integer = 1;
ppt::Bool = true,
verbose::Bool = false,
solver = Hypatia.Optimizer{_solver_type(T)})
Upper bound on the random robustness of ρ
such that it has a Schmidt number s
.
If a state $ρ$ with local dimensions $d_A$ and $d_B$ has Schmidt number $s$, then there is a PSD matrix $ω$ in the extended space $AA′B′B$, where $A′$ and $B^′$ have dimension $s$, such that $ω / s$ is separable against $AA′|B′B$ and $Π† ω Π = ρ$, where $Π = 1_A ⊗ s ψ^+ ⊗ 1_B$, and $ψ^+$ is a non-normalized maximally entangled state. Separabiity is tested with the DPS hierarchy, with n
controlling the how many copies of the $B′B$ subsystem are used.
References: Hulpke, Bruss, Lewenstein, Sanpera arXiv:quant-ph/0401118Weilenmann, Dive, Trillo, Aguilar, Navascués arXiv:1912.10056
Ket.ppt_mixture
— Functionfunction ppt_mixture(
ρ::AbstractMatrix{T},
dims::AbstractVector{<:Integer};
verbose::Bool = false,
solver = Hypatia.Optimizer{_solver_type(T)})
Lower bound on the white noise such that ρ is still a genuinely multipartite entangled state and a GME witness that detects ρ.
The set of GME states is approximated by the set of PPT mixtures, so the entanglement across the bipartitions is decided with the PPT criterion. If the state is a PPT mixture, returns a 0 matrix instead of a witness.
Reference: Jungnitsch, Moroder, Guehne arXiv:quant-ph/0401118
function ppt_mixture(
ρ::AbstractMatrix{T},
dims::AbstractVector{<:Integer},
obs::AbstractVector{<:AbstractMatrix} = Vector{Matrix}();
verbose::Bool = false,
solver = Hypatia.Optimizer{_solver_type(T)})
Lower bound on the white noise such that ρ is still a genuinely multipartite entangled state that can be detected with a witness using only the operators provided in obs
, and the values of the coefficients defining such a witness.
More precisely, if a list of observables $O_i$ is provided in the parameter obs
, the witness will be of the form $∑_i α_i O_i$ and detects ρ only using these observables. For example, using only two-body operators (and lower order) one can call
julia> two_body_basis = collect(Iterators.flatten(n_body_basis(i, 3) for i ∈ 0:2))
julia> ppt_mixture(state_ghz(2, 3), [2, 2, 2], two_body_basis)
Reference: Jungnitsch, Moroder, Guehne arXiv:quant-ph/0401118
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)
Checks if 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.
Ket.tensor_to_povm
— Functiontensor_to_povm(A::Array{T,4}, o::Vector{Int})
Converts a set of measurements in the common tensor format into a matrix of (hermitian) matrices. By default, the second argument is fixed by the size of A
. It can also contain custom number of outcomes if there are measurements with less outcomes.
Ket.povm_to_tensor
— Functionpovm_to_tensor(Axa::Vector{<:Measurement})
Converts a matrix of (hermitian) matrices into a set of measurements in the common tensor format.
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
— Functiontest_mub(B::Vector{Matrix{<:Number}})
Checks if the input bases are mutually unbiased.
Incompatibility
Ket.incompatibility_robustness
— Functionincompatibility_robustness(A::Vector{Measurement{<:Number}}, measure::String = "g")
Computes the incompatibility robustness of the measurements in the vector A
. Depending on the noise model chosen, the second argument can be "d"
(depolarizing), "r"
(random), "p"
(probabilistic), "jm"
(jointly measurable), or "g"
(generalized), see the corresponding functions below. Returns the parent POVM if return_parent = true
.
Reference: Designolle, Farkas, Kaniewski, arXiv:1906.00448
Ket.incompatibility_robustness_depolarizing
— Functionincompatibility_robustness_depolarizing(A::Vector{Measurement{<:Number}})
Computes the incompatibility depolarizing robustness of the measurements in the vector A
.
Reference: Designolle, Farkas, Kaniewski, arXiv:1906.00448
Ket.incompatibility_robustness_random
— Functionincompatibility_robustness_random(A::Vector{Measurement{<:Number}})
Computes the incompatibility random robustness of the measurements in the vector A
.
Reference: Designolle, Farkas, Kaniewski, arXiv:1906.00448
Ket.incompatibility_robustness_probabilistic
— Functionincompatibility_robustness_probabilistic(A::Vector{Measurement{<:Number}})
Computes the incompatibility probabilistic robustness of the measurements in the vector A
.
Reference: Designolle, Farkas, Kaniewski, arXiv:1906.00448
Ket.incompatibility_robustness_jointly_measurable
— Functionincompatibility_robustness_jointly_measurable(A::Vector{Measurement{<:Number}})
Computes the incompatibility jointly measurable robustness of the measurements in the vector A
.
Reference: Designolle, Farkas, Kaniewski, arXiv:1906.00448
Ket.incompatibility_robustness_generalized
— Functionincompatibility_robustness_generalized(A::Vector{Measurement{<:Number}})
Computes the incompatibility generalized robustness of the measurements in the vector A
.
Reference: Designolle, Farkas, Kaniewski, arXiv:1906.00448
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 = 3)
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.inn22
— Functioninn22([T=Float64,] n::Integer = 3)
inn22 Bell functional in Collins-Gisin notation. Local bound 1.
Reference: Śliwa, arXiv:quant-ph/0305190
Ket.gyni
— Functiongyni([T=Float64,] n::Integer)
Guess your neighbour's input nonlocal game in full probability notation. If T
is an integer type the game is unnormalized.
Reference: Mafalda et al., arXiv:1003.3844.
Ket.local_bound
— Functionlocal_bound(G::Array{T,N}; correlation = N < 4, marg = true)
Computes the local bound of a multipartite Bell functional G
given as an N
-dimensional array. If correlation
is false
, G
is assumed to be written in full probability notation. If correlation
is true
, G
is assumed to be written in correlation notation, with or without marginals depending on marg
.
Reference: Araújo, Hirsch, and Quintino, arXiv:2005.13418.
Ket.tsirelson_bound
— Functiontsirelson_bound(CG::Matrix, scenario::AbstractVecOrTuple, level)
Upper bounds the Tsirelson bound of a multipartite Bell funcional CG
, written in Collins-Gisin notation. scenario
is a tuple detailing the number of inputs and outputs, in the order (oa, ob, ..., ia, ib, ...). level
is an integer or string determining the level of the NPA hierarchy.
tsirelson_bound(FC::Matrix, level::Integer)
Upper bounds the Tsirelson bound of a bipartite Bell funcional FC
, written in full correlation notation. level
is an integer or string determining the level of the NPA hierarchy.
Ket.seesaw
— Functionseesaw(CG::Matrix, scenario::AbstractVecOrTuple, d::Integer)
Maximizes bipartite Bell functional in Collins-Gisin notation CG
using the seesaw heuristic. scenario
is a vector detailing the number of inputs and outputs, in the order [oa, ob, ia, ib]. d
is an integer determining the local dimension of the strategy.
If oa
== ob
== 2 the heuristic reduces to a bunch of eigenvalue problems. Otherwise semidefinite programming is needed and we use the assemblage version of seesaw.
References: Pál and Vértesi, arXiv:1006.3032, section II.B.1 of Tavakoli et al., arXiv:2307.02551
Ket.tensor_probability
— Functiontensor_probability(CG::Array, scenario::AbstractVecOrTuple, behaviour::Bool = false)
Takes a multipartite Bell functional CG
in Collins-Gisin notation and transforms it to full probability notation. scenario
is a tuple detailing the number of inputs and outputs, in the order (oa, ob, ..., ia, ib, ...). If behaviour
is true
do instead the transformation for behaviours. Doesn't assume normalization.
tensor_probability(FC::Matrix, behaviour::Bool = false)
Takes a bipartite Bell functional FC
in full correlator notation and transforms it to full probability notation. If behaviour
is true
do instead the transformation for behaviours. Doesn't assume normalization.
tensor_probability(rho::Hermitian, all_Aax::Vector{Measurement}...)
tensor_probability(rho::Hermitian, Aax::Vector{Measurement}, N::Integer)
Applies N sets of measurements onto a state rho
to form a probability array. If all parties apply the same measurements, use the shorthand notation.
Ket.tensor_collinsgisin
— Functiontensor_collinsgisin(p::Array, behaviour::Bool = false)
Takes a multipartite Bell functional p
in full probability notation and transforms it to Collins-Gisin notation. If behaviour
is true
do instead the transformation for behaviours. Doesn't assume normalization.
Also accepts the arguments of tensor_probability
(state and measurements) for convenience.
Ket.tensor_correlation
— Functiontensor_correlation(p::AbstractArray{T, N2}, behaviour::Bool = false; marg::Bool = true)
Converts a 2x...x2xmx...xm probability array into
- a mx...xm correlation array (no marginals)
- a (m+1)x...x(m+1) correlation array (marginals).
If behaviour
is true
do the transformation for behaviours. Doesn't assume normalization.
Also accepts the arguments of tensor_probability
(state and measurements) for convenience.
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_super_singlet_ket
— Functionstate_super_singlet_ket([T=ComplexF64,] N::Integer = 3; coeff = 1/√N!)
Produces the vector of the N
-partite N
-level singlet state.
Reference: Adán Cabello, arXiv:quant-ph/0203119
Ket.state_super_singlet
— Functionstate_super_singlet([T=ComplexF64,] N::Integer = 3; v::Real = 1, coeff = 1/√d)
Produces the N
-partite N
-level singlet state with visibility v
. This state is invariant under simultaneous rotations on all parties: (U⊗…⊗U)ρ(U⊗…⊗U)'=ρ
.
Reference: Adán Cabello, arXiv:quant-ph/0203119
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 it 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
.