Quantum Systems

PiccoloQuantumObjects.QuantumSystems.CompositeQuantumSystemType
CompositeQuantumSystem <: AbstractQuantumSystem

A composite quantum system consisting of multiple subsystems with optional coupling terms.

Composite systems represent multiple quantum subsystems (e.g., multiple qubits or oscillators) that may be coupled together. Each subsystem's Hamiltonians are automatically lifted to the full tensor product space, and subsystem drives are appended to any coupling drives.

Fields

  • H::Function: The total Hamiltonian function: (u, t) -> H(u, t)
  • G::Function: The isomorphic generator function: (u, t) -> G(u, t)
  • H_drift::SparseMatrixCSC{ComplexF64, Int}: The total drift Hamiltonian including subsystem drifts and couplings
  • H_drives::Vector{SparseMatrixCSC{ComplexF64, Int}}: All drive Hamiltonians (coupling drives + subsystem drives)
  • T_max::Float64: Maximum evolution time
  • drive_bounds::Vector{Tuple{Float64, Float64}}: Drive amplitude bounds for each control
  • n_drives::Int: Total number of control drives
  • levels::Int: Total dimension of the composite system (product of subsystem dimensions)
  • subsystem_levels::Vector{Int}: Dimensions of each subsystem
  • subsystems::Vector{QuantumSystem}: The individual quantum subsystems

See also QuantumSystem, lift_operator.

Example

# Two qubits with ZZ coupling
sys1 = QuantumSystem([PAULIS[:X]], 10.0, [(-1.0, 1.0)])
sys2 = QuantumSystem([PAULIS[:Y]], 10.0, [(-1.0, 1.0)])
H_coupling = 0.1 * kron(PAULIS[:Z], PAULIS[:Z])
csys = CompositeQuantumSystem(H_coupling, [sys1, sys2], 10.0, Float64[])
source
PiccoloQuantumObjects.QuantumSystems.CompositeQuantumSystemMethod
CompositeQuantumSystem(
    H_drift::AbstractMatrix,
    H_drives::AbstractVector{<:AbstractMatrix},
    subsystems::AbstractVector{<:QuantumSystem},
    T_max::Float64,
    drive_bounds::DriveBounds
)

Construct a CompositeQuantumSystem with coupling drift and drive terms.

Arguments

  • H_drift::AbstractMatrix: Coupling drift Hamiltonian (in full tensor product space)
  • H_drives::AbstractVector{<:AbstractMatrix}: Coupling drive Hamiltonians
  • subsystems::AbstractVector{<:QuantumSystem}: Vector of subsystems to compose
  • T_max::Float64: Maximum evolution time
  • drive_bounds::DriveBounds: Drive bounds for the coupling drives (subsystem bounds are inherited). Can be:
    • Tuples (lower, upper) for asymmetric bounds
    • Scalars which are interpreted as symmetric bounds (-value, value)

The total drift includes both the coupling drift and all subsystem drifts (automatically lifted). The total drives include coupling drives followed by all subsystem drives (automatically lifted).

Example

sys1 = QuantumSystem(PAULIS[:Z], [PAULIS[:X]], 10.0, [1.0])
sys2 = QuantumSystem([PAULIS[:Y]], 10.0, [1.0])
g12 = 0.1 * kron(PAULIS[:X], PAULIS[:X])  # coupling drift
csys = CompositeQuantumSystem(g12, Matrix{ComplexF64}[], [sys1, sys2], 10.0, Float64[])
source
PiccoloQuantumObjects.QuantumSystems.CompositeQuantumSystemMethod
CompositeQuantumSystem(
    subsystems::AbstractVector{<:QuantumSystem},
    T_max::Float64,
    drive_bounds::DriveBounds
)

Convenience constructor for a composite system with no coupling terms (neither drift nor drives).

Use this when you have independent subsystems that you want to represent in a single composite space, but without any direct coupling between them.

Arguments

  • subsystems::AbstractVector{<:QuantumSystem}: Vector of subsystems to compose
  • T_max::Float64: Maximum evolution time
  • drive_bounds::DriveBounds: Drive bounds for the coupling drives (typically empty). Can be:
    • Tuples (lower, upper) for asymmetric bounds
    • Scalars which are interpreted as symmetric bounds (-value, value)

Example

sys1 = QuantumSystem([PAULIS[:X]], 10.0, [1.0])
sys2 = QuantumSystem([PAULIS[:Y]], 10.0, [1.0])
csys = CompositeQuantumSystem([sys1, sys2], 10.0, Float64[])
source
PiccoloQuantumObjects.QuantumSystems.CompositeQuantumSystemMethod
CompositeQuantumSystem(
    H_drift::AbstractMatrix,
    subsystems::AbstractVector{<:QuantumSystem},
    T_max::Float64,
    drive_bounds::DriveBounds
)

Convenience constructor for a composite system with coupling drift but no coupling drives.

Arguments

  • H_drift::AbstractMatrix: Coupling drift Hamiltonian
  • subsystems::AbstractVector{<:QuantumSystem}: Vector of subsystems to compose
  • T_max::Float64: Maximum evolution time
  • drive_bounds::DriveBounds: Drive bounds for the coupling drives (typically empty). Can be:
    • Tuples (lower, upper) for asymmetric bounds
    • Scalars which are interpreted as symmetric bounds (-value, value)

Example

sys1 = QuantumSystem([PAULIS[:X]], 10.0, [1.0])
sys2 = QuantumSystem([PAULIS[:Y]], 10.0, [1.0])
H_coupling = 0.1 * kron(PAULIS[:Z], PAULIS[:Z])  # coupling drift
csys = CompositeQuantumSystem(H_coupling, [sys1, sys2], 10.0, Float64[])
source
PiccoloQuantumObjects.QuantumSystems.CompositeQuantumSystemMethod
CompositeQuantumSystem(
    H_drives::AbstractVector{<:AbstractMatrix},
    subsystems::AbstractVector{<:QuantumSystem},
    T_max::Float64,
    drive_bounds::DriveBounds
)

Convenience constructor for a composite system with coupling drives but no coupling drift.

Arguments

  • H_drives::AbstractVector{<:AbstractMatrix}: Coupling drive Hamiltonians
  • subsystems::AbstractVector{<:QuantumSystem}: Vector of subsystems to compose
  • T_max::Float64: Maximum evolution time
  • drive_bounds::DriveBounds: Drive bounds for the coupling drives. Can be:
    • Tuples (lower, upper) for asymmetric bounds
    • Scalars which are interpreted as symmetric bounds (-value, value)

Example

sys1 = QuantumSystem([PAULIS[:X]], 10.0, [1.0])
sys2 = QuantumSystem([PAULIS[:Y]], 10.0, [1.0])
g12 = 0.1 * kron(PAULIS[:X], PAULIS[:X])  # coupling drive
csys = CompositeQuantumSystem([g12], [sys1, sys2], 10.0, [1.0])  # symmetric bound
source
PiccoloQuantumObjects.QuantumSystems.DriveBoundsType
DriveBounds

Type alias for drive amplitude bounds input. Bounds can be specified as:

  • A tuple (lower, upper) for asymmetric bounds
  • A scalar value which is interpreted as symmetric bounds (-value, value)

Examples

drive_bounds = [(-1.0, 1.0), 0.5, (-0.3, 0.7)]
# Interpreted as: [(-1.0, 1.0), (-0.5, 0.5), (-0.3, 0.7)]
source
PiccoloQuantumObjects.QuantumSystems.OpenQuantumSystemType
OpenQuantumSystem <: AbstractQuantumSystem

A struct for storing open quantum dynamics.

Fields

  • H::Function: The Hamiltonian function: (u, t) -> H(u, t)
  • 𝒢::Function: The Lindbladian generator function: u -> 𝒢(u)
  • H_drift::SparseMatrixCSC{ComplexF64, Int}: The drift Hamiltonian
  • H_drives::Vector{SparseMatrixCSC{ComplexF64, Int}}: The drive Hamiltonians
  • T_max::Float64: Maximum evolution time
  • drive_bounds::Vector{Tuple{Float64, Float64}}: Drive amplitude bounds
  • n_drives::Int: The number of control drives
  • levels::Int: The number of levels in the system
  • dissipation_operators::Vector{SparseMatrixCSC{ComplexF64, Int}}: The dissipation operators
  • time_dependent::Bool: Whether the Hamiltonian has explicit time dependence

See also QuantumSystem.

source
PiccoloQuantumObjects.QuantumSystems.OpenQuantumSystemMethod
OpenQuantumSystem(
    H_drift::AbstractMatrix{<:Number},
    H_drives::AbstractVector{<:AbstractMatrix{<:Number}},
    T_max::Float64,
    drive_bounds::DriveBounds;
    dissipation_operators::AbstractVector{<:AbstractMatrix{<:Number}}=Matrix{ComplexF64}[]
)
OpenQuantumSystem(
    H_drift::AbstractMatrix{<:Number}, 
    T_max::Float64;
    dissipation_operators::AbstractVector{<:AbstractMatrix{<:Number}}=Matrix{ComplexF64}[]
)
OpenQuantumSystem(
    H_drives::Vector{<:AbstractMatrix{<:Number}},
    T_max::Float64, 
    drive_bounds::DriveBounds;
    dissipation_operators::AbstractVector{<:AbstractMatrix{<:Number}}=Matrix{ComplexF64}[]
)
OpenQuantumSystem(
    H::Function, 
    T_max::Float64,
    drive_bounds::DriveBounds;
    dissipation_operators::Vector{<:AbstractMatrix{<:Number}}=Matrix{ComplexF64}[]
)
OpenQuantumSystem(
    system::QuantumSystem; 
    dissipation_operators::Vector{<:AbstractMatrix{<:Number}}=Matrix{ComplexF64}[]
)

Constructs an OpenQuantumSystem object from the drift and drive Hamiltonian terms and dissipation operators. All constructors require Tmax (maximum time) and drivebounds (control bounds for each drive) to be explicitly specified.

Drive Bounds

The drive_bounds parameter can be:

  • Tuples (lower, upper) for asymmetric bounds
  • Scalars which are interpreted as symmetric bounds (-value, value)
source
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
PiccoloQuantumObjects.QuantumSystems.QuantumSystemMethod
QuantumSystem(
    H_drift::AbstractMatrix{<:Number},
    H_drives::Vector{<:AbstractMatrix{<:Number}},
    T_max::Float64,
    drive_bounds::Vector{<:Union{Tuple{Float64, Float64}, Float64}};
    time_dependent::Bool=false
)

Construct a QuantumSystem from drift and drive Hamiltonian terms.

Arguments

  • H_drift::AbstractMatrix: The drift (time-independent) Hamiltonian
  • H_drives::Vector{<:AbstractMatrix}: Vector of drive Hamiltonians, one for each control
  • T_max::Float64: Maximum evolution time
  • drive_bounds::DriveBounds: Drive amplitude bounds for each control. Can be:
    • Tuples (lower, upper) for asymmetric bounds
    • Scalars which are interpreted as symmetric bounds (-value, value)

Keyword Arguments

  • time_dependent::Bool=false: Set to true if using time-dependent modulation (typically handled at a higher level)

The resulting Hamiltonian is: H(u, t) = Hdrift + Σᵢ uᵢ * Hdrives[i]

Example

sys = QuantumSystem(
    PAULIS[:Z],                    # drift
    [PAULIS[:X], PAULIS[:Y]],      # drives
    10.0,                          # T_max
    [1.0, 1.0]                     # symmetric bounds: [(-1.0, 1.0), (-1.0, 1.0)]
)
source
PiccoloQuantumObjects.QuantumSystems.QuantumSystemMethod
QuantumSystem(H::Function, T_max::Float64, drive_bounds::Vector{<:Union{Tuple{Float64, Float64}, Float64}}; time_dependent::Bool=false)

Construct a QuantumSystem from a Hamiltonian function.

Arguments

  • H::Function: Hamiltonian function with signature (u, t) -> H(u, t) where u is the control vector and t is time
  • T_max::Float64: Maximum evolution time
  • drive_bounds::DriveBounds: Drive amplitude bounds for each control. Can be:
    • Tuples (lower, upper) for asymmetric bounds
    • Scalars which are interpreted as symmetric bounds (-value, value)

Keyword Arguments

  • time_dependent::Bool=false: Set to true if the Hamiltonian has explicit time dependence (e.g., cos(ωt) modulation)

Example

# Define a time-dependent Hamiltonian
H = (u, t) -> PAULIS[:Z] + u[1] * cos(ω * t) * PAULIS[:X]
sys = QuantumSystem(H, 10.0, [(-1.0, 1.0)]; time_dependent=true)
source
PiccoloQuantumObjects.QuantumSystems.QuantumSystemMethod
QuantumSystem(H_drives::Vector{<:AbstractMatrix}, T_max::Float64, drive_bounds::Vector; time_dependent::Bool=false)

Convenience constructor for a system with no drift Hamiltonian (H_drift = 0).

Arguments

  • H_drives::Vector{<:AbstractMatrix}: Vector of drive Hamiltonians
  • T_max::Float64: Maximum evolution time
  • drive_bounds::DriveBounds: Drive amplitude bounds for each control. Can be:
    • Tuples (lower, upper) for asymmetric bounds
    • Scalars which are interpreted as symmetric bounds (-value, value)

Example

# Using scalars for symmetric bounds
sys = QuantumSystem([PAULIS[:X], PAULIS[:Y]], 10.0, [1.0, 1.0])
# Equivalent to: drive_bounds = [(-1.0, 1.0), (-1.0, 1.0)]
source
PiccoloQuantumObjects.QuantumSystems.VariationalQuantumSystemType
VariationalQuantumSystem <: AbstractQuantumSystem

A struct for storing variational quantum dynamics, used for sensitivity and robustness analysis.

Variational systems allow exploring how the dynamics change under perturbations to the Hamiltonian. The variational operators represent directions of uncertainty or perturbation in the system.

Fields

  • H::Function: The Hamiltonian function: (u, t) -> H(u, t)
  • G::Function: The isomorphic generator function: (u, t) -> G(u, t)
  • G_vars::AbstractVector{<:Function}: Variational generator functions, one for each perturbation direction
  • T_max::Float64: Maximum evolution time
  • drive_bounds::Vector{Tuple{Float64, Float64}}: Drive amplitude bounds
  • n_drives::Int: The number of control drives in the system
  • levels::Int: The number of levels (dimension) in the system

See also QuantumSystem, OpenQuantumSystem.

source
PiccoloQuantumObjects.QuantumSystems.VariationalQuantumSystemMethod
VariationalQuantumSystem(
    H_drift::AbstractMatrix,
    H_drives::AbstractVector{<:AbstractMatrix},
    H_vars::AbstractVector{<:AbstractMatrix},
    T_max::Float64,
    drive_bounds::DriveBounds
)

Construct a VariationalQuantumSystem from drift, drive, and variational Hamiltonian terms.

Arguments

  • H_drift::AbstractMatrix: The drift (time-independent) Hamiltonian
  • H_drives::AbstractVector{<:AbstractMatrix}: Vector of drive Hamiltonians for control
  • H_vars::AbstractVector{<:AbstractMatrix}: Vector of variational Hamiltonians representing perturbation directions
  • T_max::Float64: Maximum evolution time
  • drive_bounds::DriveBounds: Drive amplitude bounds for each control. Can be:
    • Tuples (lower, upper) for asymmetric bounds
    • Scalars which are interpreted as symmetric bounds (-value, value)

The variational operators allow sensitivity analysis by exploring how dynamics change under perturbations: Hperturbed = H + Σᵢ εᵢ * Hvars[i]

Example

varsys = VariationalQuantumSystem(
    PAULIS[:Z],                    # drift
    [PAULIS[:X], PAULIS[:Y]],      # drives
    [PAULIS[:X]],                  # variational perturbations
    10.0,                          # T_max
    [1.0, 1.0]                     # symmetric bounds: [(-1.0, 1.0), (-1.0, 1.0)]
)
source
PiccoloQuantumObjects.QuantumSystems.lift_operatorFunction
lift_operator(operator::AbstractMatrix{<:Number}, i::Int, subsystem_levels::Vector{Int})
lift_operator(operator::AbstractMatrix{<:Number}, i::Int, n_qubits::Int; kwargs...)
lift_operator(operators::AbstractVector{<:AbstractMatrix{T}}, indices::AbstractVector{Int}, subsystem_levels::Vector{Int})
lift_operator(operators::AbstractVector{<:AbstractMatrix{T}}, indices::AbstractVector{Int}, n_qubits::Int; kwargs...)
lift_operator(operator::AbstractMatrix{T}, indices::AbstractVector{Int}, subsystem_levels::AbstractVector{Int})
lift_operator(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
PiccoloQuantumObjects.QuantumSystems.normalize_drive_boundsMethod
normalize_drive_bounds(bounds::DriveBounds)

Convert drive bounds to a consistent tuple format. Scalar values are converted to symmetric bounds around zero: b becomes (-b, b).

Arguments

  • bounds::DriveBounds: Input bounds, can be tuples or scalars

Returns

  • Vector{Tuple{Float64, Float64}}: Normalized bounds as tuples

Examples

# All scalars (symmetric bounds)
normalize_drive_bounds([1.0, 1.5, 0.5])
# Returns: [(-1.0, 1.0), (-1.5, 1.5), (-0.5, 0.5)]

# All tuples (asymmetric bounds)
normalize_drive_bounds([(-2.0, 3.0), (-1.0, 1.0)])
# Returns: [(-2.0, 3.0), (-1.0, 1.0)]

# Mixed types (requires explicit type annotation)
normalize_drive_bounds(Union{Float64, Tuple{Float64,Float64}}[1.0, (-2.0, 3.0)])
# Returns: [(-1.0, 1.0), (-2.0, 3.0)]
source

Gates

PiccoloQuantumObjects.Gates.GATESConstant

A constant dictionary GATES containing common quantum gate matrices as complex-valued matrices. Each gate is represented by its unitary matrix.

  • GATES[:I] - Identity: Leaves the state unchanged.
  • GATES[:X] - Pauli-X (NOT): Flips the qubit state.
  • GATES[:Y] - Pauli-Y: Rotates the qubit state around the Y-axis of the Bloch sphere.
  • GATES[:Z] - Pauli-Z: Flips the phase of the qubit state.
  • GATES[:H] - Hadamard: Creates superposition by transforming basis states.
  • GATES[:CX] - Controlled-X (CNOT): Flips the 2nd qubit (target) if the first qubit (control) is |1⟩.
  • GATES[:CZ] - Controlled-Z (CZ): Flips the phase of the 2nd qubit (target) if the 1st qubit (control) is |1⟩.
  • GATES[:XI] - Complex: A gate for complex operations.
  • GATES[:sqrtiSWAP] - Square root of iSWAP: Partially swaps two qubits with a phase.
source

Embedded Operators

PiccoloQuantumObjects.EmbeddedOperators.EmbeddedOperatorType
EmbeddedOperator

Embedded operator type to represent an operator embedded in a subspace of a larger quantum system.

Fields

  • operator::Matrix{<:Number}: Embedded operator of size prod(subsystem_levels) x prod(subsystem_levels).
  • subspace::Vector{Int}: Indices of the subspace the operator is embedded in.
  • subsystem_levels::Vector{Int}: Levels of the subsystems in the composite system.
source
PiccoloQuantumObjects.EmbeddedOperators.EmbeddedOperatorMethod
EmbeddedOperator(
    subspace_operator::AbstractMatrix{<:Number},
    subsystem_indices::AbstractVector{Int},
    subspaces::AbstractVector{<:AbstractVector{Int}},
    subsystem_levels::AbstractVector{Int}
)

Embed the subspace_operator into the provided subspaces of a composite system, where the subsystem_indices list the subspaces at which the operator is defined, and the subsystem_levels list the levels of the subsystems in which the operator is embedded.

source
PiccoloQuantumObjects.EmbeddedOperators.EmbeddedOperatorMethod
EmbeddedOperator(
    subspace_operator::AbstractMatrix{<:Number},
    subsystem_indices::AbstractVector{Int},
    subspaces::AbstractVector{<:AbstractVector{Int}},
    composite_system::CompositeQuantumSystem
)

Embed the subspace_operator into the provided subspaces of a composite system.

source
PiccoloQuantumObjects.EmbeddedOperators.EmbeddedOperatorMethod
EmbeddedOperator(subspace_operator::Matrix{<:Number}, subspace::AbstractVector{Int}, subsystem_levels::AbstractVector{Int})

Create an embedded operator. The operator is embedded at the subspace of the system spanned by the subsystem_levels.

source
PiccoloQuantumObjects.EmbeddedOperators.get_iso_vec_leakage_indicesFunction
get_iso_vec_leakage_indices(subspace::AbstractVector{Int}, levels::Int)
get_iso_vec_leakage_indices(subspaces::AbstractVector{<:AbstractVector{Int}}, subsystem_levels::AbstractVector{Int})
get_iso_vec_leakage_indices(subsystem_levels::AbstractVector{Int}; subspace=1:2)
get_iso_vec_leakage_indices(op::EmbeddedOperator)

Get the indices for the leakage in the isomorphic vector space for operators.

source
PiccoloQuantumObjects.EmbeddedOperators.get_leakage_indicesFunction
get_leakage_indices(subspace::AbstractVector{Int}, levels::Int)
get_leakage_indices(subspaces::AbstractVector{<:AbstractVector{Int}}, subsystem_levels::AbstractVector{Int})
get_leakage_indices(subsystem_levels::AbstractVector{Int}; subspace=1:2)
get_leakage_indices(op::EmbeddedOperator)

Get the indices for the states that are outside of the provided subspace of the quantum system.

source
PiccoloQuantumObjects.EmbeddedOperators.get_subspace_indicesFunction
get_subspace_indices(subspace::AbstractVector{Int}, levels::Int)
get_subspace_indices(subspaces::Vector{<:AbstractVector{Int}}, subsystem_levels::AbstractVector{Int})
get_subspace_indices(subsystem_levels::AbstractVector{Int}; subspace=1:2)
get_subspace_indices(op::EmbeddedOperator)

Get the indices for the provided subspace of the quantum system.

source

Isomorphisims

PiccoloQuantumObjects.Isomorphisms.ad_vecMethod
ad_vec(H::AbstractMatrix{ℂ}; anti::Bool=false) where ℂ <: Number

Returns the vectorized adjoint action of a matrix H:

\[\text{ad_vec}(H) = \mqty(1 & 0 \\ 0 & 1) \otimes H - (-1)^{\text{anti}} \mqty(0 & 1 \\ 1 & 0) \otimes H^*\]

source
PiccoloQuantumObjects.Isomorphisms.isoMethod
iso(H::AbstractMatrix{<:Number})

Returns the isomorphism of $H$:

\[iso(H) = \widetilde{H} = \mqty(1 & 0 \\ 0 & 1) \otimes \Re(H) + \mqty(0 & -1 \\ 1 & 0) \otimes \Im(H)\]

where $\Im(H)$ and $\Re(H)$ are the imaginary and real parts of $H$ and the tilde indicates the standard isomorphism of a complex valued matrix:

\[\widetilde{H} = \mqty(1 & 0 \\ 0 & 1) \otimes \Re(H) + \mqty(0 & -1 \\ 1 & 0) \otimes \Im(H)\]

See also Isomorphisms.G, Isomorphisms.H.

source
PiccoloQuantumObjects.Isomorphisms.var_GMethod
var_G(G::AbstractMatrix{<:Real}, G_vars::AbstractVector{<:AbstractMatrix{<:Real}})

Returns the variational generator of G with variational derivatives, G_vars.

The variational generator is

\[\text{var}_G(G, [G_a, G_b]) = \mqty( G & 0 & 0 \\ G_a & G & 0 \\ G_b & 0 & G )\]

where G is the isomorphism of a Hamiltonian and G_a and G_b are the variational derivatives of G for parameters a and b, respectively.

source

Quantum Object Utilities

PiccoloQuantumObjects.QuantumObjectUtils.haar_identityMethod
haar_identity(n::Int, radius::Number)

Generate a random unitary matrix close to the identity matrix using the Haar measure for an n-dimensional system with a given radius. The smaller the radius, the closer the matrix will be to the identity.

source

Quantum System Utilities

PiccoloQuantumObjects.QuantumSystemUtils.is_reachableMethod
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
PiccoloQuantumObjects.QuantumSystemUtils.is_reachableMethod
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
PiccoloQuantumObjects.QuantumSystemUtils.operator_algebraMethod
operator_algebra(generators; kwargs...)

Compute the Lie algebra basis for the given generators.

Arguments

  • generators::Vector{<:AbstractMatrix}: generators of the Lie algebra

Keyword Arguments

  • return_layers::Bool=false: return the Lie tree layers
  • normalize::Bool=false: normalize the basis
  • verbose::Bool=false: print information
  • remove_trace::Bool=true: remove trace from generators
source
PiccoloQuantumObjects.DirectSums.direct_sumMethod
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

Rollouts

PiccoloQuantumObjects.Rollouts.infer_is_evpMethod
infer_is_evp(integrator::Function)

Infer whether the integrator is a exponential-vector product (EVP) function.

If true, the integrator is expected to have a signature like the exponential action, expv. Otherwise, it is expected to have a signature like exp.

source
PiccoloQuantumObjects.Rollouts.open_rolloutFunction
open_rollout(
    ρ⃗₁::AbstractVector{<:Complex},
    controls::AbstractMatrix{<:Real},
    Δt::AbstractVector,
    system::OpenQuantumSystem;
    kwargs...
)

Rollout a quantum state ρ⃗₁ under the control controls for a time Δt

Arguments

  • ρ⃗₁::AbstractVector{<:Complex}: Initial state vector
  • controls::AbstractMatrix{<:Real}: Control matrix
  • Δt::AbstractVector: Time steps
  • system::OpenQuantumSystem: Quantum system

Keyword Arguments

  • show_progress::Bool=false: Show progress bar
  • integrator::Function=expv: Integrator function
  • exp_vector_product::Bool: Infer whether the integrator is an exponential-vector product
source
PiccoloQuantumObjects.Rollouts.open_rolloutMethod
open_rollout(
    ρ₁::AbstractMatrix{<:Complex},
    controls::AbstractMatrix{<:Real},
    Δt::AbstractVector,
    system::OpenQuantumSystem;
    kwargs...
)

Rollout a density matrix ρ₁ under the control controls and timesteps Δt

source
PiccoloQuantumObjects.Rollouts.open_rollout_fidelityFunction
open_rollout_fidelity(
    ρ⃗₁::AbstractVector{<:Complex},
    ρ⃗₂::AbstractVector{<:Complex},
    controls::AbstractMatrix{<:Real},
    Δt::AbstractVector,
    system::OpenQuantumSystem
)
open_rollout_fidelity(
    ρ₁::AbstractMatrix{<:Complex},
    ρ₂::AbstractMatrix{<:Complex},
    controls::AbstractMatrix{<:Real},
    Δt::AbstractVector,
    system::OpenQuantumSystem
)
open_rollout_fidelity(
    traj::NamedTrajectory,
    system::OpenQuantumSystem;
    state_name::Symbol=:ρ⃗̃,
    control_name::Symbol=:a,
    kwargs...
)

Calculate the fidelity between the final state of an open quantum system rollout and a goal state.

source
PiccoloQuantumObjects.Rollouts.rolloutFunction
rollout(
    ψ̃_init::AbstractVector{<:Real},
    controls::AbstractMatrix{<:Real},
    Δt::AbstractVector,
    system::AbstractQuantumSystem
)
rollout(
    ψ_init::AbstractVector{<:Complex},
    controls::AbstractMatrix{<:Real},
    Δt::AbstractVector,
    system::AbstractQuantumSystem
)
rollout(
    inits::AbstractVector{<:AbstractVector},
    controls::AbstractMatrix{<:Real},
    Δt::AbstractVector,
    system::AbstractQuantumSystem
)

Rollout a quantum state ψ̃_init under the control controls for a time Δt using the system system.

If exp_vector_product is true, the integrator is expected to have a signature like the exponential action, expv. Otherwise, it is expected to have a signature like exp.

Types should allow for autodifferentiable controls and times.

source
PiccoloQuantumObjects.Rollouts.rollout_fidelityFunction
rollout_fidelity(
    ψ̃_init::AbstractVector{<:Real},
    ψ̃_goal::AbstractVector{<:Real},
    controls::AbstractMatrix{<:Real},
    Δt::AbstractVector,
    system::AbstractQuantumSystem
)
rollout_fidelity(
    ψ_init::AbstractVector{<:Complex},
    ψ_goal::AbstractVector{<:Complex},
    controls::AbstractMatrix{<:Real},
    Δt::AbstractVector,
    system::AbstractQuantumSystem
)
rollout_fidelity(
    trajectory::NamedTrajectory,
    system::AbstractQuantumSystem
)

Calculate the fidelity between the final state of a rollout and a goal state.

source
PiccoloQuantumObjects.Rollouts.unitary_free_phase_fidelityMethod
unitary_free_phase_fidelity(
    U::AbstractMatrix,
    U_goal::AbstractMatrix,
    phases::AbstractVector{<:Real},
    phase_operators::AbstractVector{<:AbstractMatrix};
    subspace::AbstractVector{Int}=axes(U, 1)
)

Calculate the fidelity between unitary operators U and U_goal in the subspace, including the phase rotations about the phase_operators.

source
PiccoloQuantumObjects.Rollouts.unitary_rolloutFunction
unitary_rollout(
    Ũ⃗_init::AbstractVector{<:Real},
    controls::AbstractMatrix{<:Real},
    Δt::AbstractVector,
    system::AbstractQuantumSystem;
    kwargs...
)

Rollout a isomorphic unitary operator Ũ⃗_init under the control controls for a time Δt using the system system.

Arguments

  • Ũ⃗_init::AbstractVector{<:Real}: Initial unitary vector
  • controls::AbstractMatrix{<:Real}: Control matrix
  • Δt::AbstractVector: Time steps
  • system::AbstractQuantumSystem: Quantum system

Keyword Arguments

  • show_progress::Bool=false: Show progress bar
  • integrator::Function=expv: Integrator function
  • exp_vector_product::Bool: Infer whether the integrator is an exponential-vector product
source
PiccoloQuantumObjects.Rollouts.unitary_rollout_fidelityFunction
unitary_rollout_fidelity(
    Ũ⃗_init::AbstractVector{<:Real},
    Ũ⃗_goal::AbstractVector{<:Real},
    controls::AbstractMatrix{<:Real},
    Δt::AbstractVector,
    system::AbstractQuantumSystem;
    kwargs...
)
unitary_rollout_fidelity(
    Ũ⃗_goal::AbstractVector{<:Real},
    controls::AbstractMatrix{<:Real},
    Δt::AbstractVector,
    system::AbstractQuantumSystem;
    kwargs...
)
unitary_rollout_fidelity(
    U_init::AbstractMatrix{<:Complex},
    U_goal::AbstractMatrix{<:Complex},
    controls::AbstractMatrix{<:Real},
    Δt::AbstractVector,
    system::AbstractQuantumSystem;
    kwargs...
)
unitary_rollout_fidelity(
    U_goal::AbstractMatrix{<:Complex},
    controls::AbstractMatrix{<:Real},
    Δt::AbstractVector,
    system::AbstractQuantumSystem;
    kwargs...
)
unitary_rollout_fidelity(
    U_goal::EmbeddedOperator,
    controls::AbstractMatrix{<:Real},
    Δt::AbstractVector,
    system::AbstractQuantumSystem;
    subspace::AbstractVector{Int}=U_goal.subspace,
    kwargs...
)
unitary_rollout_fidelity(
    traj::NamedTrajectory,
    sys::AbstractQuantumSystem;
    kwargs...
)

Calculate the fidelity between the final state of a unitary rollout and a goal state. If the initial unitary is not provided, the identity operator is assumed. If phases and phase_operators are provided, the free phase unitary fidelity is calculated.

source
PiccoloQuantumObjects.Rollouts.variational_rolloutFunction
variational_rollout(
    ψ̃_init::AbstractVector{<:Real},
    controls::AbstractMatrix{<:Real},
    Δt::AbstractVector{<:Real},
    system::VariationalQuantumSystem;
    show_progress::Bool=false,
    integrator::Function=expv,
    exp_vector_product::Bool=infer_is_evp(integrator)
)
variational_rollout(ψ::Vector{<:Complex}, args...; kwargs...)
variational_rollout(inits::AbstractVector{<:AbstractVector}, args...; kwargs...)
variational_rollout(
    traj::NamedTrajectory, 
    system::AbstractQuantumSystem; 
    state_name::Symbol=:ψ̃,
    drive_name::Symbol=:a,
    kwargs...
)

Simulates the variational evolution of a quantum state under a given control trajectory.

Returns

  • Ψ̃::Matrix{<:Real}: The evolved quantum state at each timestep.
  • Ψ̃_vars::Vector{<:Matrix{<:Real}}: The variational derivatives of the quantum state with respect to the variational parameters.

Notes

This function computes the variational evolution of a quantum state using the variational generators of the system. It supports autodifferentiable controls and timesteps, making it suitable for optimization tasks. The variational derivatives are computed alongside the state evolution, enabling sensitivity analysis and gradient-based optimization.

source
PiccoloQuantumObjects.Rollouts.variational_unitary_rolloutFunction
variational_unitary_rollout(
    Ũ⃗_init::AbstractVector{<:Real},
    controls::AbstractMatrix{<:Real},
    Δt::AbstractVector{<:Real},
    system::VariationalQuantumSystem;
    show_progress::Bool=false,
    integrator::Function=expv,
    exp_vector_product::Bool=infer_is_evp(integrator)
)
variational_unitary_rollout(
    controls::AbstractMatrix{<:Real},
    Δt::AbstractVector,
    system::VariationalQuantumSystem;
    kwargs...
)
variational_unitary_rollout(
    traj::NamedTrajectory,
    system::VariationalQuantumSystem;
    unitary_name::Symbol=:Ũ⃗,
    drive_name::Symbol=:a,
    kwargs...
)

Simulates the variational evolution of a quantum state under a given control trajectory.

Returns

  • Ũ⃗::Matrix{<:Real}: The evolved unitary at each timestep.
  • Ũ⃗_vars::Vector{<:Matrix{<:Real}}: The variational derivatives of the unitary with respect to the variational parameters.

Notes

This function computes the variational evolution of a unitary using the variational generators of the system. It supports autodifferentiable controls and timesteps, making it suitable for optimization tasks. The variational derivatives are computed alongside the state evolution, enabling sensitivity analysis and gradient-based optimization.

source