Abstract Quantum Systems

using PiccoloQuantumObjects
using SparseArrays # for visualization
⊗ = kron;

Quantum Systems

The QuantumSystem type is used to represent a quantum system with a drift Hamiltonian and a set of drive Hamiltonians,

\[H = H_{\text{drift}} + \sum_i a_i H_{\text{drives}}^{(i)}\]

PiccoloQuantumObjects.QuantumSystems.QuantumSystemType
QuantumSystem(H_drift::Matrix{<:Number}, H_drives::Vector{Matrix{<:Number}}; kwargs...)
QuantumSystem(H_drift::Matrix{<:Number}; kwargs...)
QuantumSystem(H_drives::Vector{Matrix{<:Number}}; kwargs...)
QuantumSystem(H::Function, n_drives::Int; kwargs...)

Constructs a QuantumSystem object from the drift and drive Hamiltonian terms.

source

QuantumSystem's are containers for quantum dynamics. Internally, they compute the necessary isomorphisms to perform the dynamics in a real vector space.

H_drift = PAULIS[:Z]
H_drives = [PAULIS[:X], PAULIS[:Y]]
system = QuantumSystem(H_drift, H_drives)

a_drives = [1, 0]
system.H(a_drives)
2×2 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 4 stored entries:
 1.0+0.0im   1.0+0.0im
 1.0+0.0im  -1.0+0.0im

To extract the drift and drive Hamiltonians from a QuantumSystem, use the get_drift and get_drives functions.

get_drift(system) |> sparse
2×2 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 2 stored entries:
 1.0+0.0im       ⋅    
     ⋅      -1.0+0.0im

Get the X drive.

drives = get_drives(system)
drives[1] |> sparse
2×2 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 2 stored entries:
     ⋅      1.0+0.0im
 1.0+0.0im      ⋅    

And the Y drive.

drives[2] |> sparse
2×2 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 2 stored entries:
     ⋅      0.0-1.0im
 0.0+1.0im      ⋅    
Note

We can also construct a QuantumSystem directly from a Hamiltonian function. Internally, ForwardDiff.jl is used to compute the drives.

H(a) = PAULIS[:Z] + a[1] * PAULIS[:X] + a[2] * PAULIS[:Y]
system = QuantumSystem(H, 2)
get_drives(system)[1] |> sparse
2×2 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 2 stored entries:
     ⋅      1.0+0.0im
 1.0+0.0im      ⋅    

Create a noise model with a confusion matrix.

function H(a; C::Matrix{Float64}=[1.0 0.0; 0.0 1.0])
    b = C * a
    return b[1] * PAULIS.X + b[2] * PAULIS.Y
end

system = QuantumSystem(a -> H(a, C=[0.99 0.01; -0.01 1.01]), 2)
confused_drives = get_drives(system)
confused_drives[1] |> sparse
2×2 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 2 stored entries:
      ⋅       0.99+0.01im
 0.99-0.01im       ⋅    
confused_drives[2] |> sparse
2×2 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 2 stored entries:
      ⋅       0.01-1.01im
 0.01+1.01im       ⋅    

Open quantum systems

We can also construct an OpenQuantumSystem with Lindblad dynamics, enabling a user to pass a list of dissipation operators.

PiccoloQuantumObjects.QuantumSystems.OpenQuantumSystemType
OpenQuantumSystem(
    H_drift::AbstractMatrix{<:Number},
    H_drives::AbstractVector{<:AbstractMatrix{<:Number}}
    dissipation_operators::AbstractVector{<:AbstractMatrix{<:Number}};
    kwargs...
)
OpenQuantumSystem(
    H_drift::Matrix{<:Number}, H_drives::AbstractVector{Matrix{<:Number}}; 
    dissipation_operators::AbstractVector{<:AbstractMatrix{<:Number}}=Matrix{ComplexF64}[], 
    kwargs...
)
OpenQuantumSystem(H_drift::Matrix{<:Number}; kwargs...)
OpenQuantumSystem(H_drives::Vector{Matrix{<:Number}}; kwargs...)
OpenQuantumSystem(H::Function, n_drives::Int; kwargs...)

Constructs an OpenQuantumSystem object from the drift and drive Hamiltonian terms and dissipation operators.

source

Add a dephasing and annihilation error channel.

H_drives = [PAULIS[:X]]
a = annihilate(2)
dissipation_operators = [a'a, a]
system = OpenQuantumSystem(H_drives, dissipation_operators=dissipation_operators)
system.dissipation_operators[1] |> sparse
2×2 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 1 stored entry:
     ⋅          ⋅    
     ⋅      1.0+0.0im
system.dissipation_operators[2] |> sparse
2×2 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 1 stored entry:
     ⋅      1.0+0.0im
     ⋅          ⋅    
Warning

The Hamiltonian part system.H excludes the Lindblad operators. This is also true for functions that report properties of system.H, such as get_drift, get_drives, and is_reachable.

get_drift(system) |> sparse
2×2 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 0 stored entries:
     ⋅          ⋅    
     ⋅          ⋅    

Composite quantum systems

A CompositeQuantumSystem is constructed from a list of subsystems and their interactions. The interaction, in the form of drift or drive Hamiltonian, acts on the full Hilbert space. The subsystems, with their own drift and drive Hamiltonians, are internally lifted to the full Hilbert space.

system_1 = QuantumSystem([PAULIS[:X]])
system_2 = QuantumSystem([PAULIS[:Y]])
H_drift = PAULIS[:Z] ⊗ PAULIS[:Z]
system = CompositeQuantumSystem(H_drift, [system_1, system_2]);

The drift Hamiltonian is the ZZ coupling.

get_drift(system) |> sparse
4×4 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 4 stored entries:
 1.0+0.0im       ⋅           ⋅          ⋅    
     ⋅      -1.0+0.0im       ⋅          ⋅    
     ⋅           ⋅      -1.0+0.0im      ⋅    
     ⋅           ⋅           ⋅      1.0+0.0im

The drives are the X and Y operators on the first and second subsystems.

drives = get_drives(system)
drives[1] |> sparse
4×4 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 4 stored entries:
     ⋅          ⋅      1.0+0.0im      ⋅    
     ⋅          ⋅          ⋅      1.0+0.0im
 1.0+0.0im      ⋅          ⋅          ⋅    
     ⋅      1.0+0.0im      ⋅          ⋅    
drives[2] |> sparse
4×4 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 4 stored entries:
     ⋅      0.0-1.0im      ⋅          ⋅    
 0.0+1.0im      ⋅          ⋅          ⋅    
     ⋅          ⋅          ⋅      0.0-1.0im
     ⋅          ⋅      0.0+1.0im      ⋅    

The lift function

To lift operators acting on a subsystem into the full Hilbert space, use lift.

PiccoloQuantumObjects.CompositeQuantumSystems.liftFunction
lift(operator::AbstractMatrix{<:Number}, i::Int, subsystem_levels::Vector{Int})
lift(operator::AbstractMatrix{<:Number}, i::Int, n_qubits::Int; kwargs...)
lift(operators::AbstractVector{<:AbstractMatrix{T}}, indices::AbstractVector{Int}, subsystem_levels::Vector{Int})
lift(operators::AbstractVector{<:AbstractMatrix{T}}, indices::AbstractVector{Int}, n_qubits::Int; kwargs...)
lift(operator::AbstractMatrix{T}, indices::AbstractVector{Int}, subsystem_levels::AbstractVector{Int})
lift(operator::AbstractMatrix{T}, indices::AbstractVector{Int}, n_qubits::Int; kwargs...)

Lift an operator acting on the i-th subsystem within subsystem_levels to an operator acting on the entire system spanning subsystem_levels.

source

Create an a + a' operator acting on the 1st subsystem of a qutrit and qubit system.

subspace_levels = [3, 2]
lift(create(3) + annihilate(3), 1, subspace_levels) .|> real |> sparse
6×6 SparseArrays.SparseMatrixCSC{Float64, Int64} with 8 stored entries:
  ⋅    ⋅   1.0       ⋅        ⋅        ⋅ 
  ⋅    ⋅    ⋅       1.0       ⋅        ⋅ 
 1.0   ⋅    ⋅        ⋅       1.41421   ⋅ 
  ⋅   1.0   ⋅        ⋅        ⋅       1.41421
  ⋅    ⋅   1.41421   ⋅        ⋅        ⋅ 
  ⋅    ⋅    ⋅       1.41421   ⋅        ⋅ 

Create IXI operator on the 2nd qubit in a 3-qubit system.

lift(PAULIS[:X], 2, 3) .|> real |> sparse
8×8 SparseArrays.SparseMatrixCSC{Float64, Int64} with 8 stored entries:
  ⋅    ⋅   1.0   ⋅    ⋅    ⋅    ⋅    ⋅ 
  ⋅    ⋅    ⋅   1.0   ⋅    ⋅    ⋅    ⋅ 
 1.0   ⋅    ⋅    ⋅    ⋅    ⋅    ⋅    ⋅ 
  ⋅   1.0   ⋅    ⋅    ⋅    ⋅    ⋅    ⋅ 
  ⋅    ⋅    ⋅    ⋅    ⋅    ⋅   1.0   ⋅ 
  ⋅    ⋅    ⋅    ⋅    ⋅    ⋅    ⋅   1.0
  ⋅    ⋅    ⋅    ⋅   1.0   ⋅    ⋅    ⋅ 
  ⋅    ⋅    ⋅    ⋅    ⋅   1.0   ⋅    ⋅ 

Create an XX operator acting on qubits 3 and 4 in a 4-qubit system.

lift([PAULIS[:X], PAULIS[:X]], [3, 4], 4) .|> real |> sparse
16×16 SparseArrays.SparseMatrixCSC{Float64, Int64} with 16 stored entries:
⎡⡠⠊⠀⠀⠀⠀⠀⠀⎤
⎢⠀⠀⡠⠊⠀⠀⠀⠀⎥
⎢⠀⠀⠀⠀⡠⠊⠀⠀⎥
⎣⠀⠀⠀⠀⠀⠀⡠⠊⎦

We can also lift an operator that entangles different subspaces by passing the indices of the entangled subsystems.

#_Here's another way to create an XX operator acting on qubits 3 and 4 in a 4-qubit system._
lift(kron(PAULIS[:X], PAULIS[:X]), [3, 4], 4) .|> real |> sparse
16×16 SparseArrays.SparseMatrixCSC{Float64, Int64} with 16 stored entries:
⎡⡠⠊⠀⠀⠀⠀⠀⠀⎤
⎢⠀⠀⡠⠊⠀⠀⠀⠀⎥
⎢⠀⠀⠀⠀⡠⠊⠀⠀⎥
⎣⠀⠀⠀⠀⠀⠀⡠⠊⎦

Lift a CX gate acting on the 1st and 3rd qubits in a 3-qubit system.The result is independent of the state of the second qubit.

lift(GATES[:CX], [1, 3], 3) .|> real |> sparse
8×8 SparseArrays.SparseMatrixCSC{Float64, Int64} with 8 stored entries:
 1.0   ⋅    ⋅    ⋅    ⋅    ⋅    ⋅    ⋅ 
  ⋅   1.0   ⋅    ⋅    ⋅    ⋅    ⋅    ⋅ 
  ⋅    ⋅   1.0   ⋅    ⋅    ⋅    ⋅    ⋅ 
  ⋅    ⋅    ⋅   1.0   ⋅    ⋅    ⋅    ⋅ 
  ⋅    ⋅    ⋅    ⋅    ⋅   1.0   ⋅    ⋅ 
  ⋅    ⋅    ⋅    ⋅   1.0   ⋅    ⋅    ⋅ 
  ⋅    ⋅    ⋅    ⋅    ⋅    ⋅    ⋅   1.0
  ⋅    ⋅    ⋅    ⋅    ⋅    ⋅   1.0   ⋅ 

Reachability tests

Whether a quantum system can be used to reach a target state or operator can be tested by computing the dynamical Lie algebra. Access to this calculation is provided by the is_reachable function.

PiccoloQuantumObjects.QuantumSystemUtils.is_reachableFunction
is_reachable(gate, hamiltonians; kwargs...)

Check if the gate is reachable using the given hamiltonians.

Arguments

  • gate::AbstractMatrix: target gate
  • hamiltonians::AbstractVector{<:AbstractMatrix}: generators of the Lie algebra

Keyword Arguments

  • subspace::AbstractVector{<:Int}=1:size(gate, 1): subspace indices
  • compute_basis::Bool=true: compute the basis or use the Hamiltonians directly
  • remove_trace::Bool=true: remove trace from generators
  • verbose::Bool=true: print information about the operator algebra
  • atol::Float32=eps(Float32): absolute tolerance

See also QuantumSystemUtils.operator_algebra.

source
is_reachable(gate::AbstractMatrix{<:Number}, system::AbstractQuantumSystem; kwargs...)

Check if the gate is reachable using the given system.

Keyword Arguments

  • use_drift::Bool=true: include drift Hamiltonian in the generators
  • kwargs...: keyword arguments for is_reachable
source

Y can be reached by commuting Z and X.

system = QuantumSystem(PAULIS[:Z], [PAULIS[:X]])
is_reachable(PAULIS[:Y], system)
true

Y cannot be reached by X alone.

system = QuantumSystem([PAULIS[:X]])
is_reachable(PAULIS[:Y], system)
false

This page was generated using Literate.jl.