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(u, t) = H_{\text{drift}} + \sum_i u_i H_{\text{drives}}^{(i)}\]

where $u$ is the control vector and $t$ is time.

PiccoloQuantumObjects.QuantumSystems.QuantumSystemType
QuantumSystem <: AbstractQuantumSystem

A struct for storing quantum dynamics.

Fields

  • H::Function: The Hamiltonian function: (u, t) -> H(u, t), where u is the control vector and t is time
  • G::Function: The isomorphic generator function: (u, t) -> G(u, t), including the Hamiltonian mapped to superoperator space
  • H_drift::SparseMatrixCSC{ComplexF64, Int}: The drift Hamiltonian (time-independent component)
  • H_drives::Vector{SparseMatrixCSC{ComplexF64, Int}}: The drive Hamiltonians (control-dependent components)
  • T_max::Float64: Maximum evolution time
  • drive_bounds::Vector{Tuple{Float64, Float64}}: Drive amplitude bounds for each control (lower, upper)
  • n_drives::Int: The number of control drives in the system
  • levels::Int: The number of levels (dimension) in the system
  • time_dependent::Bool: Whether the Hamiltonian has explicit time dependence beyond control modulation

See also OpenQuantumSystem, VariationalQuantumSystem.

source

QuantumSystem's are containers for quantum dynamics. Internally, they compute the necessary isomorphisms to perform the dynamics in a real vector space. All systems require explicit specification of T_max (maximum time) and drive_bounds (control bounds).

H_drift = PAULIS[:Z]
H_drives = [PAULIS[:X], PAULIS[:Y]]
T_max = 10.0
drive_bounds = [(-1.0, 1.0), (-1.0, 1.0)]
system = QuantumSystem(H_drift, H_drives, T_max, drive_bounds)

u_controls = [1.0, 0.0]
t = 0.0
system.H(u_controls, t)
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. The function must accept (u, t) arguments where u is the control vector and t is time.

H(u, t) = PAULIS[:Z] + u[1] * PAULIS[:X] + u[2] * PAULIS[:Y]
system = QuantumSystem(H, 10.0, [(-1.0, 1.0), (-1.0, 1.0)])
system.H([1.0, 0.0], 0.0) |> sparse
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

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]], 1.0, [(-1.0, 1.0)])
is_reachable(PAULIS[:Y], system)
true

Y cannot be reached by X alone.

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

Direct sums

The direct sum of two quantum systems is constructed with the direct_sum function.

PiccoloQuantumObjects.DirectSums.direct_sumFunction
direct_sum(A::AbstractMatrix, B::AbstractMatrix)

Returns the direct sum of two matrices.

source
direct_sum(A::SparseMatrixCSC, B::SparseMatrixCSC)

Returns the direct sum of two sparse matrices.

source
direct_sum(Ã⃗::AbstractVector, B̃⃗::AbstractVector)

Returns the direct sum of two iso_vec operators.

source
direct_sum(sys1::QuantumSystem, sys2::QuantumSystem)

Returns the direct sum of two QuantumSystem objects.

Constructs a new system where the Hilbert space is the direct sum of the two input systems: H = H₁ ⊕ H₂ = [H₁ 0 ] [0 H₂]

Both systems must have the same number of drives. The resulting system uses sys1's Tmax and drivebounds.

Example

sys1 = QuantumSystem([PAULIS[:X]], 10.0, [(-1.0, 1.0)])
sys2 = QuantumSystem([PAULIS[:Y]], 10.0, [(-1.0, 1.0)])
sys_combined = direct_sum(sys1, sys2)
source

Create a pair of non-interacting qubits.

system_1 = QuantumSystem(PAULIS[:Z], [PAULIS[:X], PAULIS[:Y]], 1.0, [(-1.0, 1.0), (-1.0, 1.0)])
system_2 = QuantumSystem(PAULIS[:Z], [PAULIS[:X], PAULIS[:Y]], 1.0, [(-1.0, 1.0), (-1.0, 1.0)])
system = direct_sum(system_1, system_2)
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

This page was generated using Literate.jl.