Quantum Systems
using PiccoloQuantumObjects
using SparseArrays # for visualization
⊗ = kron;PiccoloQuantumObjects.QuantumSystems.AbstractQuantumSystem — Type
AbstractQuantumSystemAbstract type for defining systems.
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.QuantumSystem — Type
QuantumSystem <: AbstractQuantumSystemA 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 timeG::Function: The isomorphic generator function: (u, t) -> G(u, t), including the Hamiltonian mapped to superoperator spaceH_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 timedrive_bounds::Vector{Tuple{Float64, Float64}}: Drive amplitude bounds for each control (lower, upper)n_drives::Int: The number of control drives in the systemlevels::Int: The number of levels (dimension) in the systemtime_dependent::Bool: Whether the Hamiltonian has explicit time dependence beyond control modulation
See also OpenQuantumSystem, VariationalQuantumSystem.
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.0imTo extract the drift and drive Hamiltonians from a QuantumSystem, use the get_drift and get_drives functions.
get_drift(system) |> sparse2×2 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 2 stored entries:
1.0+0.0im ⋅
⋅ -1.0+0.0imGet the X drive.
drives = get_drives(system)
drives[1] |> sparse2×2 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 2 stored entries:
⋅ 1.0+0.0im
1.0+0.0im ⋅ And the Y drive.
drives[2] |> sparse2×2 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 2 stored entries:
⋅ 0.0-1.0im
0.0+1.0im ⋅ 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) |> sparse2×2 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 4 stored entries:
1.0+0.0im 1.0+0.0im
1.0+0.0im -1.0+0.0imReachability 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_reachable — Function
is_reachable(gate, hamiltonians; kwargs...)Check if the gate is reachable using the given hamiltonians.
Arguments
gate::AbstractMatrix: target gatehamiltonians::AbstractVector{<:AbstractMatrix}: generators of the Lie algebra
Keyword Arguments
subspace::AbstractVector{<:Int}=1:size(gate, 1): subspace indicescompute_basis::Bool=true: compute the basis or use the Hamiltonians directlyremove_trace::Bool=true: remove trace from generatorsverbose::Bool=true: print information about the operator algebraatol::Float32=eps(Float32): absolute tolerance
See also QuantumSystemUtils.operator_algebra.
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 generatorskwargs...: keyword arguments foris_reachable
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)trueY cannot be reached by X alone.
system = QuantumSystem([PAULIS[:X]], 1.0, [(-1.0, 1.0)])
is_reachable(PAULIS[:Y], system)falseDirect sums
The direct sum of two quantum systems is constructed with the direct_sum function.
PiccoloQuantumObjects.DirectSums.direct_sum — Function
direct_sum(A::AbstractMatrix, B::AbstractMatrix)Returns the direct sum of two matrices.
direct_sum(A::SparseMatrixCSC, B::SparseMatrixCSC)Returns the direct sum of two sparse matrices.
direct_sum(Ã⃗::AbstractVector, B̃⃗::AbstractVector)Returns the direct sum of two iso_vec operators.
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)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) |> sparse4×4 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 4 stored entries:
1.0+0.0im ⋅ ⋅ ⋅
⋅ -1.0+0.0im ⋅ ⋅
⋅ ⋅ 1.0+0.0im ⋅
⋅ ⋅ ⋅ -1.0+0.0imThis page was generated using Literate.jl.