Quantum Systems
PiccoloQuantumObjects.QuantumSystems.AbstractQuantumSystem — TypeAbstractQuantumSystemAbstract type for defining systems.
PiccoloQuantumObjects.QuantumSystems.CompositeQuantumSystem — TypeCompositeQuantumSystem <: AbstractQuantumSystemA 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 couplingsH_drives::Vector{SparseMatrixCSC{ComplexF64, Int}}: All drive Hamiltonians (coupling drives + subsystem drives)T_max::Float64: Maximum evolution timedrive_bounds::Vector{Tuple{Float64, Float64}}: Drive amplitude bounds for each controln_drives::Int: Total number of control driveslevels::Int: Total dimension of the composite system (product of subsystem dimensions)subsystem_levels::Vector{Int}: Dimensions of each subsystemsubsystems::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[])PiccoloQuantumObjects.QuantumSystems.CompositeQuantumSystem — MethodCompositeQuantumSystem(
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 Hamiltonianssubsystems::AbstractVector{<:QuantumSystem}: Vector of subsystems to composeT_max::Float64: Maximum evolution timedrive_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)
- Tuples
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[])PiccoloQuantumObjects.QuantumSystems.CompositeQuantumSystem — MethodCompositeQuantumSystem(
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 composeT_max::Float64: Maximum evolution timedrive_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)
- Tuples
Example
sys1 = QuantumSystem([PAULIS[:X]], 10.0, [1.0])
sys2 = QuantumSystem([PAULIS[:Y]], 10.0, [1.0])
csys = CompositeQuantumSystem([sys1, sys2], 10.0, Float64[])PiccoloQuantumObjects.QuantumSystems.CompositeQuantumSystem — MethodCompositeQuantumSystem(
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 Hamiltoniansubsystems::AbstractVector{<:QuantumSystem}: Vector of subsystems to composeT_max::Float64: Maximum evolution timedrive_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)
- Tuples
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[])PiccoloQuantumObjects.QuantumSystems.CompositeQuantumSystem — MethodCompositeQuantumSystem(
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 Hamiltonianssubsystems::AbstractVector{<:QuantumSystem}: Vector of subsystems to composeT_max::Float64: Maximum evolution timedrive_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)
- Tuples
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 boundPiccoloQuantumObjects.QuantumSystems.DriveBounds — TypeDriveBoundsType 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)]PiccoloQuantumObjects.QuantumSystems.OpenQuantumSystem — TypeOpenQuantumSystem <: AbstractQuantumSystemA 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 HamiltonianH_drives::Vector{SparseMatrixCSC{ComplexF64, Int}}: The drive HamiltoniansT_max::Float64: Maximum evolution timedrive_bounds::Vector{Tuple{Float64, Float64}}: Drive amplitude boundsn_drives::Int: The number of control driveslevels::Int: The number of levels in the systemdissipation_operators::Vector{SparseMatrixCSC{ComplexF64, Int}}: The dissipation operators
See also QuantumSystem.
PiccoloQuantumObjects.QuantumSystems.OpenQuantumSystem — MethodOpenQuantumSystem(
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)
PiccoloQuantumObjects.QuantumSystems.QuantumSystem — TypeQuantumSystem <: 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 system
See also OpenQuantumSystem, VariationalQuantumSystem.
PiccoloQuantumObjects.QuantumSystems.QuantumSystem — MethodQuantumSystem(
H_drift::AbstractMatrix{<:Number},
H_drives::Vector{<:AbstractMatrix{<:Number}},
T_max::Float64,
drive_bounds::DriveBounds
)Construct a QuantumSystem from drift and drive Hamiltonian terms.
Arguments
H_drift::AbstractMatrix: The drift (time-independent) HamiltonianH_drives::Vector{<:AbstractMatrix}: Vector of drive Hamiltonians, one for each controlT_max::Float64: Maximum evolution timedrive_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)
- Tuples
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)]
)PiccoloQuantumObjects.QuantumSystems.QuantumSystem — MethodQuantumSystem(H::Function, T_max::Float64, drive_bounds::DriveBounds)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 timeT_max::Float64: Maximum evolution timedrive_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)
- Tuples
Example
# Define a time-dependent Hamiltonian
H = (u, t) -> PAULIS[:Z] + u[1] * PAULIS[:X] + u[2] * PAULIS[:Y]
# Using symmetric bounds (scalars)
sys = QuantumSystem(H, 10.0, [1.0, 1.0])
# Equivalent to: [(-1.0, 1.0), (-1.0, 1.0)]PiccoloQuantumObjects.QuantumSystems.QuantumSystem — MethodQuantumSystem(H_drift::AbstractMatrix, T_max::Float64)Convenience constructor for a system with only a drift Hamiltonian (no drives).
Example
sys = QuantumSystem(PAULIS[:Z], 10.0)PiccoloQuantumObjects.QuantumSystems.QuantumSystem — MethodQuantumSystem(H_drives::Vector{<:AbstractMatrix}, T_max::Float64, drive_bounds::DriveBounds)Convenience constructor for a system with no drift Hamiltonian (H_drift = 0).
Arguments
H_drives::Vector{<:AbstractMatrix}: Vector of drive HamiltoniansT_max::Float64: Maximum evolution timedrive_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)
- Tuples
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)]PiccoloQuantumObjects.QuantumSystems.VariationalQuantumSystem — TypeVariationalQuantumSystem <: AbstractQuantumSystemA 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 directionT_max::Float64: Maximum evolution timedrive_bounds::Vector{Tuple{Float64, Float64}}: Drive amplitude boundsn_drives::Int: The number of control drives in the systemlevels::Int: The number of levels (dimension) in the system
See also QuantumSystem, OpenQuantumSystem.
PiccoloQuantumObjects.QuantumSystems.VariationalQuantumSystem — MethodVariationalQuantumSystem(
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) HamiltonianH_drives::AbstractVector{<:AbstractMatrix}: Vector of drive Hamiltonians for controlH_vars::AbstractVector{<:AbstractMatrix}: Vector of variational Hamiltonians representing perturbation directionsT_max::Float64: Maximum evolution timedrive_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)
- Tuples
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)]
)PiccoloQuantumObjects.QuantumSystems.get_drift — Methodget_drift(sys::AbstractQuantumSystem)Returns the drift Hamiltonian of the system.
PiccoloQuantumObjects.QuantumSystems.get_drives — Methodget_drives(sys::AbstractQuantumSystem)Returns the drive Hamiltonians of the system.
PiccoloQuantumObjects.QuantumSystems.lift_operator — Functionlift_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.
PiccoloQuantumObjects.QuantumSystems.normalize_drive_bounds — Methodnormalize_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)]Gates
PiccoloQuantumObjects.Gates.GATES — ConstantA 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.
PiccoloQuantumObjects.Gates.PAULIS — ConstantThe 2×2 Pauli matrics and identity.
Embedded Operators
PiccoloQuantumObjects.EmbeddedOperators.AbstractPiccoloOperator — TypeAbstractPiccoloOperatorUnion type for operators.
PiccoloQuantumObjects.EmbeddedOperators.EmbeddedOperator — TypeEmbeddedOperatorEmbedded operator type to represent an operator embedded in a subspace of a larger quantum system.
Fields
operator::Matrix{<:Number}: Embedded operator of sizeprod(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.
PiccoloQuantumObjects.EmbeddedOperators.EmbeddedOperator — MethodEmbeddedOperator(
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.
PiccoloQuantumObjects.EmbeddedOperators.EmbeddedOperator — MethodEmbeddedOperator(
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.
PiccoloQuantumObjects.EmbeddedOperators.EmbeddedOperator — MethodEmbeddedOperator(subspace_operator::AbstractMatrix{<:Number}, system::QuantumSystem; kwargs...)Embed the subspace_operator into a quantum system.
PiccoloQuantumObjects.EmbeddedOperators.EmbeddedOperator — MethodEmbeddedOperator(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.
PiccoloQuantumObjects.EmbeddedOperators.embed — Methodembed(subspace_operator::AbstractMatrix{<:Number}, embedded_operator::EmbeddedOperator)Embed the subspace_operator in the subspace of a larger embedded_operator.
PiccoloQuantumObjects.EmbeddedOperators.embed — Methodembed(operator::AbstractMatrix{<:Number}, subspace::AbstractVector{Int}, levels::Int)Embed an operator in the subspace of a larger matrix of size levels x levels.
PiccoloQuantumObjects.EmbeddedOperators.get_enr_subspace_indices — Methodget_enr_subspace_indices(excitation_restriction::Int, subsystem_levels::AbstractVector{Int})Get the indices for the subspace of the quantum system with an excitation restriction.
PiccoloQuantumObjects.EmbeddedOperators.get_iso_vec_leakage_indices — Functionget_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.
PiccoloQuantumObjects.EmbeddedOperators.get_iso_vec_subspace_indices — Functionget_iso_vec_subspace_indices(subspace::AbstractVector{Int}, subsystem_levels::AbstractVector{Int})
get_iso_vec_subspace_indices(op::EmbeddedOperator)Get the indices for the subspace in the isomorphic vector space for operators.
PiccoloQuantumObjects.EmbeddedOperators.get_leakage_indices — Functionget_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.
PiccoloQuantumObjects.EmbeddedOperators.get_subspace_indices — Functionget_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.
PiccoloQuantumObjects.EmbeddedOperators.unembed — Methodunembed(matrix::AbstractMatrix{<:Number}, subspace::AbstractVector{Int})Unembed a subspace operator from the matrix. This is equivalent to calling matrix[subspace, subspace].
PiccoloQuantumObjects.EmbeddedOperators.unembed — Methodunembed(op::AbstractMatrix, embedded_op::EmbeddedOperator)Unembed a sub-matrix from the op at the subspace defined by embedded_op.
PiccoloQuantumObjects.EmbeddedOperators.unembed — Methodunembed(embedded_op::EmbeddedOperator)Unembed an embedded operator, returning the original operator.
Isomorphisims
PiccoloQuantumObjects.Isomorphisms.G — MethodG(H::AbstractMatrix)::Matrix{Float64}Returns the isomorphism of $-iH$, i.e. $G(H) = \text{iso}(-iH)$.
See also Isomorphisms.iso, Isomorphisms.H.
PiccoloQuantumObjects.Isomorphisms.H — MethodH(G::AbstractMatrix{<:Real})Returns the inverse of $G(H) = iso(-iH)$, i.e. returns H.
See also Isomorphisms.iso, Isomorphisms.G.
PiccoloQuantumObjects.Isomorphisms.ad_vec — Methodad_vec(H::AbstractMatrix{ℂ}; anti::Bool=false) where ℂ <: NumberReturns 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^*\]
PiccoloQuantumObjects.Isomorphisms.density_to_iso_vec — Methoddensity_to_iso_vec(ρ::AbstractMatrix{<:Number})Returns the isomorphism ρ⃗̃ = ket_to_iso(vec(ρ)) of a density matrix ρ
PiccoloQuantumObjects.Isomorphisms.iso — Methodiso(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.
PiccoloQuantumObjects.Isomorphisms.iso_D — Methodiso_D(L::AbstractMatrix{ℂ}) where ℂ <: NumberReturns the isomorphic representation of the Lindblad dissipator L.
PiccoloQuantumObjects.Isomorphisms.iso_operator_to_iso_vec — Methodiso_operator_to_iso_vec(Ũ::AbstractMatrix{ℝ}) where ℝ <: RealConvert a real matrix Ũ representing an isomorphism operator into a real vector.
PiccoloQuantumObjects.Isomorphisms.iso_operator_to_operator — Methodiso_operator_to_operator(Ũ)PiccoloQuantumObjects.Isomorphisms.iso_to_ket — Methodiso_to_ket(ψ̃::AbstractVector{<:Real})Convert a real isomorphism vector ψ̃ into a ket vector.
PiccoloQuantumObjects.Isomorphisms.iso_vec_to_density — Methodiso_vec_to_density(ρ⃗̃::AbstractVector{<:Real})Returns the density matrix ρ from its isomorphism ρ⃗̃
PiccoloQuantumObjects.Isomorphisms.iso_vec_to_iso_operator — Methodiso_vec_to_iso_operator(Ũ⃗::AbstractVector{ℝ}) where ℝ <: RealConvert a real vector Ũ⃗ into a real matrix representing an isomorphism operator.
PiccoloQuantumObjects.Isomorphisms.iso_vec_to_operator — Methodiso_vec_to_operator(Ũ⃗::AbstractVector{ℝ}) where ℝ <: RealConvert a real vector Ũ⃗ into a complex matrix representing an operator.
PiccoloQuantumObjects.Isomorphisms.ket_to_iso — Methodket_to_iso(ψ::AbstractVector{<:Number})Convert a ket vector ψ into a complex vector with real and imaginary parts.
PiccoloQuantumObjects.Isomorphisms.mat — Methodmat(x::AbstractVector)Convert a vector x into a square matrix. The length of x must be a perfect square.
PiccoloQuantumObjects.Isomorphisms.operator_to_iso_operator — Methodoperator_to_iso_operator(U)PiccoloQuantumObjects.Isomorphisms.operator_to_iso_vec — Methodoperator_to_iso_vec(U::AbstractMatrix{ℂ}) where ℂ <: NumberConvert a complex matrix U representing an operator into a real vector.
PiccoloQuantumObjects.Isomorphisms.var_G — Methodvar_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.
Quantum Object Utilities
PiccoloQuantumObjects.QuantumObjectUtils.annihilate — Methodannihilate(levels::Int)Get the annihilation operator for a system with levels.
PiccoloQuantumObjects.QuantumObjectUtils.create — Methodcreate(levels::Int)Get the creation operator for a system with levels.
PiccoloQuantumObjects.QuantumObjectUtils.haar_identity — Methodhaar_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.
PiccoloQuantumObjects.QuantumObjectUtils.haar_random — Methodhaar_random(n::Int)Generate a random unitary matrix using the Haar measure for an n-dimensional system.
PiccoloQuantumObjects.QuantumObjectUtils.ket_from_bitstring — Methodket_from_bitstring(ket::String)Get the state vector for a qubit system given a ket string ket of 0s and 1s.
PiccoloQuantumObjects.QuantumObjectUtils.ket_from_string — Methodket_from_string(
ket::String,
levels::Vector{Int};
level_dict=Dict(:g => 0, :e => 1, :f => 2, :h => 3, :i => 4, :j => 5, :k => 6, :l => 7),
return_states=false
)Construct a quantum state from a string ket representation.
PiccoloQuantumObjects.QuantumObjectUtils.operator_from_string — Methodoperator_from_string(operator::String; lookup=PAULIS)Reduce the string (each character is one key) via operators from a dictionary.
Quantum System Utilities
PiccoloQuantumObjects.QuantumSystemUtils.is_linearly_dependent — Methodis_linearly_dependent(M::AbstractMatrix; eps=eps(Float32), verbose=true)Check if the columns of the matrix M are linearly dependent.
PiccoloQuantumObjects.QuantumSystemUtils.is_reachable — Methodis_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
PiccoloQuantumObjects.QuantumSystemUtils.is_reachable — Methodis_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.
PiccoloQuantumObjects.QuantumSystemUtils.operator_algebra — Methodoperator_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 layersnormalize::Bool=false: normalize the basisverbose::Bool=false: print informationremove_trace::Bool=true: remove trace from generators
PiccoloQuantumObjects.DirectSums.direct_sum — Methoddirect_sum(A::AbstractMatrix, B::AbstractMatrix)Returns the direct sum of two matrices.
PiccoloQuantumObjects.DirectSums.direct_sum — Methoddirect_sum(Ã⃗::AbstractVector, B̃⃗::AbstractVector)Returns the direct sum of two iso_vec operators.
PiccoloQuantumObjects.DirectSums.direct_sum — Methoddirect_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)PiccoloQuantumObjects.DirectSums.direct_sum — Methoddirect_sum(A::SparseMatrixCSC, B::SparseMatrixCSC)Returns the direct sum of two sparse matrices.
Rollouts
PiccoloQuantumObjects.Rollouts.fidelity — Methodfidelity(ρ::AbstractMatrix{<:Number}, ρ_goal::AbstractMatrix{<:Number})Calculate the fidelity between two density matrices ρ and ρ_goal.
PiccoloQuantumObjects.Rollouts.fidelity — Methodfidelity(ψ::AbstractVector{<:Number}, ψ_goal::AbstractVector{<:Number})Calculate the fidelity between two quantum states ψ and ψ_goal.
PiccoloQuantumObjects.Rollouts.free_phase — Methodfree_phase(phases::AbstractVector{<:Real}, phase_operators::AbstractVector{<:AbstractMatrix{<:ℂ}})Rotate the phase_operators by the phases and return the Kronecker product.
PiccoloQuantumObjects.Rollouts.infer_is_evp — Methodinfer_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.
PiccoloQuantumObjects.Rollouts.open_rollout — Functionopen_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 vectorcontrols::AbstractMatrix{<:Real}: Control matrixΔt::AbstractVector: Time stepssystem::OpenQuantumSystem: Quantum system
Keyword Arguments
show_progress::Bool=false: Show progress barintegrator::Function=expv: Integrator functionexp_vector_product::Bool: Infer whether the integrator is an exponential-vector product
PiccoloQuantumObjects.Rollouts.open_rollout — Methodopen_rollout(
ρ₁::AbstractMatrix{<:Complex},
controls::AbstractMatrix{<:Real},
Δt::AbstractVector,
system::OpenQuantumSystem;
kwargs...
)Rollout a density matrix ρ₁ under the control controls and timesteps Δt
PiccoloQuantumObjects.Rollouts.open_rollout_fidelity — Functionopen_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.
PiccoloQuantumObjects.Rollouts.rollout — Functionrollout(
ψ̃_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.
PiccoloQuantumObjects.Rollouts.rollout_fidelity — Functionrollout_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.
PiccoloQuantumObjects.Rollouts.unitary_fidelity — Methodunitary_fidelity(U::AbstractMatrix{<:Number}, U_goal::AbstractMatrix{<:Number})Calculate the fidelity between unitary operators U and U_goal in the subspace.
PiccoloQuantumObjects.Rollouts.unitary_free_phase_fidelity — Methodunitary_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.
PiccoloQuantumObjects.Rollouts.unitary_rollout — Functionunitary_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 vectorcontrols::AbstractMatrix{<:Real}: Control matrixΔt::AbstractVector: Time stepssystem::AbstractQuantumSystem: Quantum system
Keyword Arguments
show_progress::Bool=false: Show progress barintegrator::Function=expv: Integrator functionexp_vector_product::Bool: Infer whether the integrator is an exponential-vector product
PiccoloQuantumObjects.Rollouts.unitary_rollout_fidelity — Functionunitary_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.
PiccoloQuantumObjects.Rollouts.variational_rollout — Functionvariational_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.
PiccoloQuantumObjects.Rollouts.variational_unitary_rollout — Functionvariational_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.