API reference
EntanglementDetection.AlternatingSeparableLMO
— TypeAlternatingSeparableLMO{T, N, MB <: AbstractMatrix{Complex{T}}} <: SeparableLMO{T, N}
AlternatingSeparableLMO
implements compute_extreme_point(lmo, direction)
which returns a pure product state used in Frank-Wolfe algorithms. The method used is an alternating algorithm starting from random pure states on each party and alternatively optimizing each reduced state via an eigendecomposition.
Type parameters:
T
: element type of the correlation tensorN
: number of partiesMB
: type of the matrix basis
Fields:
dims
: dimensions of the reduced state on each partymatrix_basis
: matrix basis of the correlation tensormax_iter
: maximum number of alternation stepsthreshold
: threshold to stop the alternationnb
: number of random rounds to find the possible global optimalworkspace
: contains fields pre-allocated performance-critical functionstmp
: temporary vector for fast scalar products of bipartite tensors
EntanglementDetection.CrossPolytopeSubdivision
— TypeCrossPolytopeSubdivision{T} <: ApproximationQuantumStates{T}
Approximates the complex sphere in dimension d
via an edgewise subdivision of the cross-polytope on the real sphere in dimension 2d-1
.
EntanglementDetection.EnumeratingSeparableLMO
— TypeEnumeratingSeparableLMO{T, N, MB <: AbstractMatrix{Complex{T}}} <: SeparableLMO{T, N}
EnumeratingSeparableLMO
implements compute_extreme_point(lmo, direction)
which returns a pure product state used in Frank-Wolfe algorithms. The method used is the enumeration of a fine approximation of the sets of states on all parties except the first one, for which an eigendecomposition suffices.
Type parameters:
T
: element type of the correlation tensorN
: number of partiesMB
: type of the matrix basis
Fields:
dims
: dimensions of the reduced state on each partymatrix_basis
: matrix basis of the correlation tensorworkspace
: contains fields pre-allocated performance-critical functionstmp
: temporary vector for fast scalar products of bipartite tensors
EntanglementDetection.KSeparableLMO
— TypeKSeparableLMO{T, N, LMO} <: SeparableLMO{T, N}
KSeparableLMO implements compute_extreme_point(lmo::LMO, direction)
which returns a pure product state by considering all K-partite partitions by iterating over all LMO{T, K}.
Type parameters:
T
: element type of the correlation tensorN
: number of partiesMB
: type of the matrix basisLMO
: type of the LMO for k-partite separability
Fields:
dims
: dimensions of the reduced state on each partymatrix_basis
: matrix basis of the correlation tensorlmos
: list of LMOs for all K-partitionspartitions
: list of K-partite partitions
EntanglementDetection.PolytopePhasedComplexSphere
— TypePolytopePhasedComplexSphere{T} <: ApproximationQuantumStates{T}
Approximates the complex sphere in dimension d
via a polytope on the real sphere in dimension d
and d-1
phases.
EntanglementDetection.PolytopeReducedComplexSphere
— TypePolytopeReducedComplexSphere{T} <: ApproximationQuantumStates{T}
Approximates the complex sphere in dimension d
via a polytope on the real sphere in dimension 2d-1
.
EntanglementDetection.PureState
— TypePureState{T, N} <: AbstractArray{T, N}
Represents a pure product state. Each subsystem is a pure state stored as a tensor PureState.tensors[n].
EntanglementDetection.Workspace
— TypeWorkspace{T, N}
Structure for initial pre-allocation of performance-critical functions.
EntanglementDetection._correlation_tensor_ket!
— Method_correlation_tensor_ket!(tensor::Vector{T}, φ::Vector{Complex{T}}, matrix_basis)
Convert a ket φ
to a correlation tensor Vector{T}, with same subspace dimension d
.
EntanglementDetection._eigmin!
— Method_eigmin!(ket::Vector, matrix::Matrix)
Computes the minimal real eigenvalue and updates ket
in place The variable matrix
of size d × d also gets overwritten. For BLAS-compatible types, uses LAPACK.syev!
for d ≤ 5 and LAPACK.syevr!
otherwise. For other types, falls back to eigen!
.
EntanglementDetection._reduced_tensor!
— Method_reduced_tensor!(tensor::Vector{T}, pure_tensors::NTuple{N, Vector{T}}, dir::Array{T, N}, j::Int, s::Vector{Int} = setdiff(1:N, j)) where {T <: Real, N}
Computes the correlation tensor of the j
-th subsystem (the tensor-version of the partial trace). When N=2, j=1, computes ⟨ϕ2|dir|ϕ2⟩ for dir ∈ H₁ ⊗ H₂
EntanglementDetection.correlation_tensor
— Methodcorrelation_tensor(ρ::Matrix{T}, dims::NTuple{N, Int})
Convert a dimension-(a)symmetry density matrix ρ
to a correlation tensor Array{T, N}, with subspace dimensions dims
.
EntanglementDetection.density_matrix
— Methodfunction density_matrix(tensor::Vector{T}, dims::NTuple{N,Int}, matrix_basis = _gellmann(T, dims)) where {T <: Real} where {N}
Convert tensors of pure states to density matrix (for eigendecomposition).
EntanglementDetection.entangled_bound_with_white_noise
— Methodentangled_bound_with_white_noise(ρ::AbstractMatrix{CT}, dims::NTuple{N, Int}; atol = 1e-4) where {CT <: Number, N}
Given a quantum state ρ
, using PPT criterion calculate the upper bound for the amount of white noise it can tolerate such that ρ = (1-p) * ρ + p * I/d
is entangled. Note this is corresponding to the fully separability. Returns the upper bound p
.
EntanglementDetection.entanglement_detection
— Methodentanglement_detection(ρ::AbstractMatrix{T}, dims::NTuple{N, Int}; measure, algorithm, kwargs...) where {T <: Number, N}
Given a quantum state ρ
with subsystem dimensions as dims
, make a judgment about whether it is entangled or separable. If the argument dims
is omitted equally-sized subsystems are assumed, which is solving on the symmetry bipartite separable space.
Returns a named tuple (ent, atoms, witness)
with:
ent
true/false/nothing, the judgment as entangled/separable/can't tellwitness
the witness for the entangledρ
,nothing
for separableρ
decompose
a tuple with weights and pure separable states, the best separable decomposition for separableρ
or for the closest separable state of entangledρ
EntanglementDetection.entanglement_robustness
— Methodentanglement_robustness(ρ_p::Function, p_list::Array{Vector{T}}, dims::NTuple{N, Int}; monotone, robust_monitor, kwargs...)
Decide the entanglement and separable regions for a family of quantum state ρ_p
.
Inputs:
ρ_p
the function for a family of quantum statesp_list
the parameter region forρ_p(p)
monotone
the monotony of the function if with only one parameter (default: true)robust_monitor
the monitor for the calculation of the robustness problem
Returns a named tuple (ent_range, nan_range, sep_range)
with:
ent_range
the parameter range forp
such thatρ_p(p)
is entangled.nan_range
the parameter range forp
such thatρ_p(p)
can not be decided.sep_range
the parameter range forp
such thatρ_p(p)
is separable.
EntanglementDetection.entanglement_witness
— Methodentanglement_witness(ρ::AbstractMatrix, σ::AbstractMatrix; dims, kwargs...)
Given a direction from a state σ
in the separable space to a entangled state ρ
, find the best witness that detect the entanglement of the state along this direction:
W = (σ-ρ + Tr[σ(ρ-σ)]*I)/||ρ-σ||
which satisfies ∀σ ∈ SEP, Tr(Wσ) ≥ 0, and ∃ρ, Tr(Wρ) < 0.
Returns a named tuple (W, α, ε, σ, ϕ)
with:
W
the witness operatorα
= Tr[ϕ(ρ-σ)]/||ρ-σ||, the overlapε
the ε-net, normalized by ||ρ-σ||β
= Tr[σ(ρ-σ)]/||ρ-σ||, the overlapσ
decide the direction fromσ
toρ
ϕ
the closest pure separable state
EntanglementDetection.entanglement_witness
— Methodentanglement_witness(ρ::AbstractMatrix{T}; dims, kwargs...)
Given a quantum state ρ
, find the best witness that detect the entanglement of the state, defined by
W = (σ-ρ + Tr[σ(ρ-σ)]*I)/||ρ-σ||
which satisfies ∀σ ∈ SEP, Tr(Wσ) ≥ 0, and ∃ρ, Tr(Wρ) < 0.
Returns a tuple (W, σ, α)
with:
W
the witness operatorσ
the closest pure separable stateα
= Tr[σ(ρ-σ)], the overlap
EntanglementDetection.separable_ball_criterion
— Methodseparable_ball_criterion(ρ, σ, dims; kwargs...)
Given a quantum state ρ
and a separable state σ
, define a quantum state
ρ' = (1+t)/t * ρ - 1/t * σ, s.t. ρ = t/(1-t) * ρ' + 1/(1-t) * σ.
By confirming ρ' inside a separable ball to confirm the separability of ρ
.
Returns a named tuple (sep, atoms, witness)
with:
sep
true/false, the judgment of separabilityp
= t/(1-t), the weight forρ'
EntanglementDetection.separable_ball_criterion
— Methodseparable_ball_criterion(ρ, dims; kwargs...)
Given a quantum state ρ
with subsystem dimensions as dims
, define a quantum state
ρ' = (1+t)/t * ρ - 1/t * σ, s.t. ρ = t/(1-t) * ρ' + 1/(1-t) * σ.
where σ
is the closest separable pure state. By confirming ρ'
inside a separable ball to confirm the separability of ρ
.
Returns a named tuple (sep, atoms, witness)
with:
sep
true/false, the judgment of separabilityp
= t/(1-t), the weight forρ'
EntanglementDetection.separable_ball_membership
— Methodseparable_ball_membership(ρ::AbstractMatrix{T}, dims::NTuple{N, Int}; kwargs...)
Judges whether a given quantum state ρ
in a separable ball with subsystem dimensions as dims
.
EntanglementDetection.separable_ball_radius
— Methodseparable_ball_radius(::Type{T}, dims::NTuple{N, Int})
Lower bound of the radius of a separable ball, centered around the maximally mixed state.
Reference:
EntanglementDetection.separable_bound_with_white_noise
— Methodseparable_bound_with_white_noise(::Type{T}, dims::NTuple{N, Int}, ent_bound::T, distance::T) where {T <: Real, N}
Consider a quantum state mixed with white noise as ρ(p) = (1-p) * ρ + p * I/d
for the separability problem.
Given
dims
for the structure of the separable spaceent_bound
the noise level such that ρ(ent_bound
) could be entangled (not necessary)distance
= min{σ ∈ SEP}||ρ(`entbound`) - σ||, the distance between the state and the separable space
Return the noise level such that ρ(sep_bound
) must be separable.
EntanglementDetection.separable_distance
— Methodseparable_distance(ρ::AbstractMatrix{CT}, dims::NTuple{N, Int}, lmo::LMO=AlternatingSeparableLMO(float(real(CT)), dims); measure, fw_algorithm, kwargs...)
separable_distance(C::Array{T, N}; matrix_basis, measure, fw_algorithm, kwargs...)
Computes the distance between the quantum density matrix ρ
and the separable space under a specific measure
via a specific fw_algorithm
:
f(ρ) = min_{σ ∈ SEP} g(ρ,σ)
For the density matrix ρ
, if the argument dims
is omitted equally-sized subsystems are assumed, which is solving on the symmetry bipartite separable space. For the correlation tensor C
, if the argument matrix_basis
is omitted, Gell-Mann matrix is assumed, which is Pauli basis for qubit systems.
The quantum state can also be given by a correlation tensor C
corresponding to the experimental data from a set of (over-)completed matrix_basis
. If the argument matrix_basis
is omitted the generalized Gell-Mann basis are assumed.
The measure
g(ρ,σ) can be set as
"2-norm"
,"relative-entropy"
,- squared
"Bures metric"
.
The fw_algorithm
can be used as
FrankWolfe.frank_wolfe
FrankWolfe.lazified_conditional_gradient
FrankWolfe.away_frank_wolfe
FrankWolfe.blended_pairwise_conditional_gradient
Returns a named tuple (σ, v, primal)
with:
σ
the closest density matrix in the separable spacev
the closest pure separable state ket on the boundary of the separable spaceprimal
primal value f(x), the distance to the separable spaceactive_set
all the pure separable states, which combined to the closest separable stateσ
lmo
the structure for related computation
EntanglementDetection.white_noise_robustness
— Methodfunction whitenoiserobustness(ρ::AbstractMatrix{CT}, dims::NTuple{N, Int}, lmo::LMO; entproof, sepsearch, noise_atol, kwargs...) where {CT <: Number, N, LMO <: SeparableLMO}
Given a quantum state ρ
, calculate the upper and lower bounds for the amount of white noise it can tolerate to be separable.
ent_proof
:true
/false
use the enumerate method to guarantee the entanglementsep_search
:true
/false
use the linear search method to find a possible better separable boundnoise_atol
: the numerical accuracy for the noise level
Returns a tuple (ent_bound, sep_bound)
.