Measurements
Distance measures
PastaQ.fidelity
— Functionfidelity(ψ::MPS, ϕ::MPS; kwargs...)
Quantum state fidelity between two wavefunctions:
fidelity(ψ::MPS, ρ::MPO; kwargs...)
fidelity(ρ::MPO, ψ::MPS; kwargs...)
Quantum state fidelity between an MPS wavefunction and an MPO density operator:
fidelity(Ψ::MPS, ϱ::LPDO{MPO}; cutoff::Float64 = 1e-15)
fidelity(ϱ::LPDO{MPO}, ψ::MPS; kwargs...)
Quantum state fidelity between an MPS wavefunction and a LPDO density operator $\varrho=XX^\dagger$
fidelity(A::MPO, B::MPO; process::Bool = false, cutoff::Float64 = 1e-15)
Fidelity $F$ between two MPOs $A$ and $B$. Implements the following:
- If $A$ and $B$ are density operators, $F$ is the full quantum state fidelity
Note: this scales exponentially with the number $n$ of qubits, as it involves a full diagonalization.
- If $A$ and $B$ are unitary operators (i.e. rank-1 channels) and
process = true
,
$F$ is the process fidelity
where $|\Phi_j\rangle = |j\rangle\rangle$ is the MPS corresponding to the vectorization of the unitary operator.
- If $A$ is a Choi matrix and
B
is a unitary operator (or viceversa), return the process fidelity
- If $A$ and $B$ are both Choi matrices, return the full process fidelity
which, as above, scale exponentially with $n$.
fidelity(A::ITensor, B::ITensor)
Compute the quantum fidelity between two ITensors. Wrap each one in the ITensorState and the ITensorOperator according to the index structure to allow dispatch.
fidelity(A::ITensorState, B::ITensorState; kwargs...)
State fidelity between two wavefunctions
fidelity(A::ITensorState, B::ITensorOperator; kwargs...)
fidelity(A::ITensorOperator, B::ITensorState; kwargs...)
State fidelity between a wavefunction and a density operator
fidelity(A::ITensorOperator, B::ITensorOperator; kwargs...)
Fidelity between two operators. Wrap the Choi type as for the TN (above)
PastaQ.frobenius_distance
— Functionfrobenius_distance(ρ::Union{MPO, LPDO}, σ::Union{MPO, LPDO})
Compute the trace norm of the difference between two LPDOs or MPOs:
Projective measurements
PastaQ.fullpreparations
— Functionfullpreparations(n::Int; local_input_states = "Pauli")
Generate the full set of $n$-qubit input states built out of a collection of $D=M^n$ states out of $M$ single-qubit states. Predefined options:
local_input_states = "Pauli"
: $D=6^n$ Pauli eigenstateslocal_input_states = "Tetra"
: $D=4^n$ 1-qubit SIC-POVM
PastaQ.fullbases
— Functionfullbases(n::Int; local_basis = "Pauli")
Generate the full set of measurement bases for a choice of local single-qubit basis set. Predefined option:
local_basis = "Pauli"
: set of $3^n$ Pauli measurement bases
PastaQ.randompreparations
— Functionrandompreparations(n::Int, npreps::Int;
local_input_state = "Pauli")
Generate npreps
random input states to a quantum circuit. Each $n$-qubit input state is selected according to the following options:
local_input_states = "Pauli"
: $D=6^n$ Pauli eigenstateslocal_input_states = "Tetra"
: $D=4^n$ 1-qubit SIC-POVM
The input states can also be set to a user-defined set, local_input_states = ["A","B","C",...]
, assuming single-qubit states $|A\rangle$, $|B\rangle$ have been properly defined.
PastaQ.randombases
— Functionrandombases(n::Int, nbases::Int; local_basis = "Pauli")
Generate nbases
measurement bases composed by $n$ single-qubit bases. By default, each local basis is randomly selected between Pauli bases ["X","Y","Z"]
, with "Z"
being the default basis where the quantum state is written.
The measurement bases can also be defined by the user, local_basis = [O1, O2,...]
, as long as the single-qubit Hermitian operators $O_j$ are defined.
PastaQ.getsamples
— Functiongetsamples(M::Union{MPS,MPO}, nshots::Int; kwargs...)
getsamples(T::ITensor, nshots::Int)
Perform nshots
projective measurements of a wavefunction $|\psi\rangle$ or density operator $\rho$ in the MPS/MPO reference basis. Each measurement consists of a binary vector $\sigma = (\sigma_1,\sigma_2,\dots)$, drawn from the probabilty distributions:
- $P(\sigma) = |\langle\sigma|\psi\rangle|^2$, if $M = |\psi\rangle$ is an
MPS
. - $P(\sigma) = \langle\sigma|\rho|\sigma\rangle$ : if $M = \rho$ is an
MPO
.
getsamples(M::Union{MPS,MPO,ITensor}, bases::Matrix::Matrix{<:String}, nshots::int; kwargs...)
getsamples(M::Union{MPS,MPO,ITensor}, bases::Vector{<:Vector}, nshots::Int; kwargs...) =
Generate a set of measurements acccording to a set of input bases
, performing nshots
measurements per basis. For a single measurement, a depth-1 unitary $U$ is applied to the input state $M$ according to the basis
. The probability of recording outcome $\sigma = (\sigma_1,\sigma_2,\dots)$ in the basis defined by $U$ is
- $P(\sigma) = |\langle\sigma|U\psi\rangle|^2$, if $M = |\psi\rangle$ is an
MPS
. - $P(\sigma) = \langle\sigma|U\rho U^\dagger|\sigma\rangle$, if $M = \rho$ is an
MPO
.
getsamples(
M::Union{LPDO,MPO,ITensor},
preps::Matrix,
bases::Matrix,
nshots::Int;
kwargs...
)
Generate a set of process measurement data acccording to a set of input states preps
and measurement bases
, performing nshots
measurements per configuration. For a single measurement, the input state $|\xi\rangle=\otimes_j|\xi_j\rangle$ is evolved according to the channel $M$, and then measured in a given measurement basis. The probability that the final state returns outcome $\sigma = (\sigma_1,\sigma_2,\dots)$ in the basis defined by $U$ is given by
where $\rho_\xi = |\xi\rangle\langle\xi|$ and $\Lambda_M$ is the Choi matrix.