Library
Problem Templates
QuantumCollocation.ProblemTemplates.MinimumTimeProblem — Method
MinimumTimeProblem(qcp::QuantumControlProblem; kwargs...)Convert an existing quantum control problem to minimum-time optimization.
IMPORTANT: This function requires an existing QuantumControlProblem (e.g., from SmoothPulseProblem). It cannot be created directly from a quantum trajectory. The workflow is:
- Create base problem with
SmoothPulseProblem(or similar) - Solve base problem to get feasible solution
- Convert to minimum-time with
MinimumTimeProblem
This ensures the problem starts from a good initialization and maintains solution quality through the final fidelity constraint.
Type Dispatch
Automatically handles different quantum trajectory types through the type parameter:
QuantumControlProblem{UnitaryTrajectory}→ UsesFinalUnitaryFidelityConstraintQuantumControlProblem{KetTrajectory}→ UsesFinalKetFidelityConstraintQuantumControlProblem{DensityTrajectory}→ Not yet implemented
The optimization problem is:
\[\begin{aligned} \underset{\vec{\tilde{q}}, u, \Delta t}{\text{minimize}} & \quad J_{\text{original}}(\vec{\tilde{q}}, u) + D \sum_t \Delta t_t \\ \text{ subject to } & \quad \text{original dynamics \& constraints} \\ & F_{\text{final}} \geq F_{\text{threshold}} \\ & \quad \Delta t_{\text{min}} \leq \Delta t_t \leq \Delta t_{\text{max}} \\ \end{aligned}\]
where q represents the quantum state (unitary, ket, or density matrix).
Arguments
qcp::QuantumControlProblem: Existing quantum control problem to convert
Keyword Arguments
final_fidelity::Float64=0.99: Minimum fidelity constraint at final timeD::Float64=100.0: Weight on minimum-time objective ∑Δtpiccolo_options::PiccoloOptions=PiccoloOptions(): Piccolo solver options
Returns
QuantumControlProblem: New problem with minimum-time objective and fidelity constraint
Examples
# Standard workflow
sys = QuantumSystem(H_drift, H_drives, drive_bounds)
pulse = ZeroOrderPulse(0.1 * randn(n_drives, N), collect(range(0.0, T, length=N)))
qtraj = UnitaryTrajectory(sys, pulse, U_goal)
# Step 1: Create and solve base smooth pulse problem (with Δt_bounds for free time)
qcp_smooth = SmoothPulseProblem(qtraj, N; Q=100.0, R=1e-2, Δt_bounds=(0.01, 0.5))
solve!(qcp_smooth; max_iter=100)
# Step 2: Convert to minimum-time
qcp_mintime = MinimumTimeProblem(qcp_smooth; final_fidelity=0.99, D=100.0)
solve!(qcp_mintime; max_iter=100)
# Compare durations
duration_before = sum(get_timesteps(get_trajectory(qcp_smooth)))
duration_after = sum(get_timesteps(get_trajectory(qcp_mintime)))
@assert duration_after <= duration_before
# Nested transformations also work
qcp_final = MinimumTimeProblem(
RobustnessProblem(qcp_smooth); # Future feature
final_fidelity=0.95
)Convenience Constructors
You can also update the goal when creating minimum-time problem:
# Different goal for minimum-time optimization
qcp_mintime = MinimumTimeProblem(qcp_smooth; goal=U_goal_new, final_fidelity=0.98)QuantumCollocation.ProblemTemplates.SamplingProblem — Method
SamplingProblem(qcp::QuantumControlProblem, systems::Vector{<:AbstractQuantumSystem}; kwargs...)Construct a SamplingProblem from an existing QuantumControlProblem and a list of systems.
This creates a robust optimization problem where the controls are shared across all systems, but each system evolves according to its own dynamics. The objective is the weighted sum of fidelity objectives for each system.
Arguments
qcp::QuantumControlProblem: The base problem (defines nominal trajectory, objective, etc.)systems::Vector{<:AbstractQuantumSystem}: List of systems to optimize over
Keyword Arguments
weights::Vector{Float64}=fill(1.0, length(systems)): Weights for each systemQ::Float64=100.0: Weight on infidelity objective (explicit, not extracted from base problem)piccolo_options::PiccoloOptions=PiccoloOptions(): Options for the solver
Returns
QuantumControlProblem{SamplingTrajectory}: A new problem with the sampling trajectory
QuantumCollocation.ProblemTemplates.SmoothPulseProblem — Method
SmoothPulseProblem(qtraj::AbstractQuantumTrajectory{<:ZeroOrderPulse}, N::Int; kwargs...)Construct a QuantumControlProblem for smooth pulse optimization with piecewise constant controls.
Note: This problem template is for ZeroOrderPulse only. For spline-based pulses (LinearSplinePulse, CubicSplinePulse), use SplinePulseProblem instead.
The problem adds discrete derivative variables (du, ddu) that:
- Regularize control changes between timesteps
- Enforce smoothness via
DerivativeIntegratorconstraints
Arguments
qtraj::AbstractQuantumTrajectory{<:ZeroOrderPulse}: Quantum trajectory with piecewise constant pulseN::Int: Number of timesteps for discretization
Keyword Arguments
integrator::Union{Nothing, AbstractIntegrator, Vector{<:AbstractIntegrator}}=nothing: Optional custom integrator(s). If not provided, uses BilinearIntegrator (which does not support global variables). A custom integrator is required whenglobal_namesis specified.global_names::Union{Nothing, Vector{Symbol}}=nothing: Names of global variables to optimize. Requires a custom integrator (e.g., HermitianExponentialIntegrator from Piccolissimo) that supports global variables.global_bounds::Union{Nothing, Dict{Symbol, Union{Float64, Tuple{Float64, Float64}}}}=nothing: Bounds for global variables. Keys are variable names, values are either a scalar (symmetric bounds ±value) or a tuple (lower, upper).du_bound::Float64=Inf: Bound on discrete first derivative (controls jump rate)ddu_bound::Float64=1.0: Bound on discrete second derivative (controls acceleration)Q::Float64=100.0: Weight on infidelity/objectiveR::Float64=1e-2: Weight on regularization terms (u, u̇, ü)R_u::Union{Float64, Vector{Float64}}=R: Weight on control regularizationR_du::Union{Float64, Vector{Float64}}=R: Weight on first derivative regularizationR_ddu::Union{Float64, Vector{Float64}}=R: Weight on second derivative regularizationconstraints::Vector{<:AbstractConstraint}=AbstractConstraint[]: Additional constraintspiccolo_options::PiccoloOptions=PiccoloOptions(): Piccolo solver options
Returns
QuantumControlProblem: Wrapper containing quantum trajectory and optimization problem
Examples
# Unitary gate synthesis with piecewise constant pulse
sys = QuantumSystem(H_drift, H_drives, drive_bounds)
pulse = ZeroOrderPulse(0.1 * randn(n_drives, N), collect(range(0.0, T, length=N)))
qtraj = UnitaryTrajectory(sys, pulse, U_goal)
qcp = SmoothPulseProblem(qtraj, N; Q=100.0, R=1e-2)
solve!(qcp; max_iter=100)
# Quantum state transfer
pulse = ZeroOrderPulse(0.1 * randn(n_drives, N), collect(range(0.0, T, length=N)))
qtraj = KetTrajectory(sys, pulse, ψ_init, ψ_goal)
qcp = SmoothPulseProblem(qtraj, N; Q=50.0, R=1e-3)
solve!(qcp)See also: SplinePulseProblem for spline-based pulses.
QuantumCollocation.ProblemTemplates.SmoothPulseProblem — Method
SmoothPulseProblem(qtraj::MultiKetTrajectory{<:ZeroOrderPulse}, N::Int; kwargs...)Construct a QuantumControlProblem for smooth pulse optimization over an ensemble of ket state transfers with piecewise constant controls.
This handles the case where you want to optimize a single pulse that achieves multiple state transfers simultaneously (e.g., |0⟩→|1⟩ and |1⟩→|0⟩ for an X gate via state transfer).
Note: This problem template is for ZeroOrderPulse only. For spline-based pulses, use SplinePulseProblem instead.
Arguments
qtraj::MultiKetTrajectory{<:ZeroOrderPulse}: Ensemble of ket state transfers with piecewise constant pulseN::Int: Number of timesteps for the discretization
Keyword Arguments
integrator::Union{Nothing, AbstractIntegrator, Vector{<:AbstractIntegrator}}=nothing: Optional custom integrator(s). If not provided, the defaultBilinearIntegratoris used. Whenglobal_namesis specified, you must supply a custom integrator here (i.e., do not rely on the defaultBilinearIntegrator) that supports global variables.global_names::Union{Nothing, Vector{Symbol}}=nothing: Names of global variables to optimize. Requires a custom integrator provided viaintegrator(e.g.,HermitianExponentialIntegratorfrom Piccolissimo) that supports global variables.global_bounds::Union{Nothing, Dict{Symbol, Union{Float64, Tuple{Float64, Float64}}}}=nothing: Bounds for global variables. Keys are variable names, values are either a scalar (symmetric bounds ±value) or a tuple (lower, upper).du_bound::Float64=Inf: Bound on discrete first derivativeddu_bound::Float64=1.0: Bound on discrete second derivativeQ::Float64=100.0: Weight on infidelity/objectiveR::Float64=1e-2: Weight on regularization terms (u, u̇, ü)R_u::Union{Float64, Vector{Float64}}=R: Weight on control regularizationR_du::Union{Float64, Vector{Float64}}=R: Weight on first derivative regularizationR_ddu::Union{Float64, Vector{Float64}}=R: Weight on second derivative regularizationconstraints::Vector{<:AbstractConstraint}=AbstractConstraint[]: Additional constraintspiccolo_options::PiccoloOptions=PiccoloOptions(): Piccolo solver options
Returns
QuantumControlProblem{MultiKetTrajectory}: Wrapper containing ensemble trajectory and optimization problem
Examples
# Create ensemble for X gate via state transfer
sys = QuantumSystem(H_drift, H_drives, drive_bounds)
pulse = ZeroOrderPulse(0.1 * randn(n_drives, N), collect(range(0.0, T, length=N)))
ψ0 = ComplexF64[1.0, 0.0]
ψ1 = ComplexF64[0.0, 1.0]
ensemble_qtraj = MultiKetTrajectory(sys, pulse, [ψ0, ψ1], [ψ1, ψ0])
qcp = SmoothPulseProblem(ensemble_qtraj, N; Q=100.0, R=1e-2)
solve!(qcp; max_iter=100)See also: SplinePulseProblem for spline-based pulses.
QuantumCollocation.ProblemTemplates.SmoothPulseProblem — Method
SmoothPulseProblem(qtraj::AbstractQuantumTrajectory, N::Int; kwargs...)Fallback method that provides helpful error for non-ZeroOrderPulse types.
QuantumCollocation.ProblemTemplates.SplinePulseProblem — Method
SplinePulseProblem(qtraj::AbstractQuantumTrajectory{<:AbstractSplinePulse}, N::Int; kwargs...)Construct a QuantumControlProblem for spline-based pulse optimization.
Unlike SmoothPulseProblem (which uses piecewise constant controls with discrete smoothing variables), this problem template is designed for spline-based pulses where the derivative variables (du) are the actual spline coefficients or slopes.
Pulse Type Semantics
LinearSplinePulse: The du variable represents the slope between knots. A DerivativeIntegrator constraint enforces du[k] = (u[k+1] - u[k]) / Δt, making the slopes consistent with the linear interpolation. This constraint ensures mathematical rigor while allowing slope regularization/bounds.
CubicSplinePulse (Hermite spline): The du variable is the tangent/derivative at each knot point, which is a true independent degree of freedom in Hermite interpolation. No DerivativeIntegrator is added - the optimizer can adjust both :u and :du independently.
Mathematical Notes
- LinearSplinePulse: Always adds
:duandDerivativeIntegratorto enforce slope consistency - CubicSplinePulse:
:duvalues are Hermite tangents (unconstrained, only regularized)
Both pulse types always have :du components in the trajectory, simplifying integrator implementations.
Arguments
qtraj::AbstractQuantumTrajectory{<:AbstractSplinePulse}: Quantum trajectory with spline pulseN::Int: Number of timesteps for the discretization
Keyword Arguments
integrator::Union{Nothing, AbstractIntegrator, Vector{<:AbstractIntegrator}}=nothing: Optional custom integrator(s). If not provided, usesBilinearIntegrator(which does not support global variables). A custom integrator is required whenglobal_namesis specified.global_names::Union{Nothing, Vector{Symbol}}=nothing: Names of global variables to optimize. Requires a custom integrator (e.g.,SplineIntegratorfrom Piccolissimo) that supports global variables.global_bounds::Union{Nothing, Dict{Symbol, Union{Float64, Tuple{Float64, Float64}}}}=nothing: Bounds for global variables. Keys are variable names, values are either a scalar (symmetric bounds ±value) or a tuple (lower, upper).du_bound::Float64=Inf: Bound on derivative (slope) magnitudeQ::Float64=100.0: Weight on infidelity/objectiveR::Float64=1e-2: Weight on regularization termsR_u::Union{Float64, Vector{Float64}}=R: Weight on control regularizationR_du::Union{Float64, Vector{Float64}}=R: Weight on derivative regularizationconstraints::Vector{<:AbstractConstraint}=AbstractConstraint[]: Additional constraintspiccolo_options::PiccoloOptions=PiccoloOptions(): Piccolo solver options
Returns
QuantumControlProblem{<:AbstractQuantumTrajectory}: Wrapper containing trajectory and optimization problem
Examples
# Linear spline pulse
sys = QuantumSystem(H_drift, H_drives, drive_bounds)
pulse = LinearSplinePulse(0.1 * randn(n_drives, N), collect(range(0.0, T, length=N)))
qtraj = UnitaryTrajectory(sys, pulse, U_goal)
qcp = SplinePulseProblem(qtraj, N; Q=100.0, R=1e-2, du_bound=10.0)
solve!(qcp; max_iter=100)See also: SmoothPulseProblem for piecewise constant pulses with discrete smoothing.
QuantumCollocation.ProblemTemplates.SplinePulseProblem — Method
SplinePulseProblem(qtraj::MultiKetTrajectory{<:AbstractSplinePulse}, N; kwargs...)Create a spline-based trajectory optimization problem for ensemble ket state transfers.
Uses coherent fidelity objective (phases must align) for gate implementation.
Arguments
qtraj::MultiKetTrajectory{<:AbstractSplinePulse}: Ensemble trajectory with spline pulseN::Int: Number of timesteps
Keyword Arguments
Same as the base SplinePulseProblem method.
QuantumCollocation.ProblemTemplates.SplinePulseProblem — Method
SplinePulseProblem(qtraj::AbstractQuantumTrajectory, N::Int; kwargs...)Fallback method that provides helpful error for non-spline pulse types.
QuantumCollocation.ProblemTemplates._ensemble_ket_objective — Method
_ensemble_ket_objective(qtraj::MultiKetTrajectory, traj, state_names, weights, goals, Q)Create a coherent fidelity objective for ensemble state transfers.
For ensemble trajectories (implementing a gate via multiple state transfers), we use coherent fidelity: Fcoherent = |1/n ∑ᵢ ⟨ψᵢgoal|ψᵢ⟩|²
This requires all state overlaps to have aligned phases, which is essential for gate implementation (the gate should have a single global phase).
QuantumCollocation.ProblemTemplates._final_fidelity_constraint — Method
_final_fidelity_constraint(qtraj::MultiKetTrajectory, final_fidelity, traj)Create a coherent fidelity constraint for an MultiKetTrajectory.
Uses coherent fidelity: F = |1/n ∑ᵢ ⟨ψᵢ_goal|ψᵢ⟩|²
This enforces that all state transfers have aligned global phases, which is essential when implementing a gate via state transfer (e.g., X gate via |0⟩→|1⟩ and |1⟩→|0⟩).
QuantumCollocation.ProblemTemplates.add_global_bounds_constraints! — Method
add_global_bounds_constraints!(constraints, global_bounds, traj; verbose=false)Add GlobalBoundsConstraint entries for each global variable specified in global_bounds.
Converts bounds from user-friendly formats to the format expected by GlobalBoundsConstraint:
Float64: Symmetric scalar bounds (applied symmetrically to all dimensions)Tuple{Float64, Float64}: Asymmetric scalar bounds (expanded to vectors)VectororTuple{Vector, Vector}: Already in correct format (passed through)
Modifies constraints in place.
QuantumCollocation.ProblemTemplates.extract_regularization — Method
extract_regularization(objective, state_sym::Symbol, new_traj::NamedTrajectory) -> AbstractObjectiveExtract regularization terms (non-state-dependent objectives) from a composite objective, filtering to only include terms for variables that exist in the new trajectory.
Used by SamplingProblem to extract shared regularizers (e.g., control penalty) from the base problem while excluding regularizers for variables that don't exist in the sampling trajectory (e.g., :du, :ddu which are added by SmoothPulseProblem).
QuantumCollocation.ProblemTemplates.sampling_state_objective — Method
sampling_state_objective(qtraj, traj, state_sym, Q)Create the state-dependent objective for a sampling member. Dispatches on quantum trajectory type.
Quantum Objectives
QuantumCollocation.QuantumObjectives.CoherentKetInfidelityObjective — Method
CoherentKetInfidelityObjective(ψ_goals, ψ̃_names, traj; Q=100.0)Create a terminal objective for coherent ket state infidelity across multiple states.
Coherent fidelity is defined as: Fcoherent = |1/n ∑ᵢ ⟨ψᵢgoal|ψᵢ⟩|²
Unlike incoherent fidelity (average of individual |⟨ψᵢ_goal|ψᵢ⟩|²), coherent fidelity requires all state overlaps to have aligned phases. This is essential when implementing a gate via multiple state transfers - the gate should have a single global phase, not independent phases per state.
Arguments
ψ_goals::Vector{<:AbstractVector{<:Complex}}: Target ket statesψ̃_names::Vector{Symbol}: Names of isomorphic state variables in trajectorytraj::NamedTrajectory: The trajectory
Keyword Arguments
Q::Float64=100.0: Weight on the infidelity objective
Example
# For implementing X gate via |0⟩→|1⟩ and |1⟩→|0⟩
goals = [ComplexF64[0, 1], ComplexF64[1, 0]]
names = [:ψ̃1, :ψ̃2]
obj = CoherentKetInfidelityObjective(goals, names, traj; Q=100.0)QuantumCollocation.QuantumObjectives.KetInfidelityObjective — Method
KetInfidelityObjective(ψ_goal, ψ̃_name, traj; Q=100.0)Create a terminal objective for ket state infidelity with an explicit goal state.
This variant is useful for SamplingProblem and EnsembleTrajectory where the goal is shared across multiple state variables that don't have individual goals in traj.goal.
Arguments
ψ_goal::AbstractVector{<:Complex}: The target ket state (complex vector)ψ̃_name::Symbol: Name of the isomorphic state variable in the trajectorytraj::NamedTrajectory: The trajectory
Keyword Arguments
Q::Float64=100.0: Weight on the infidelity objective
QuantumCollocation.QuantumObjectives.KetInfidelityObjective — Method
KetInfidelityObjective(ψ̃_name, traj; Q=100.0)Create a terminal objective for ket state infidelity, using the goal from traj.goal[ψ̃_name].
QuantumCollocation.QuantumObjectives.LeakageObjective — Method
LeakageObjective(indices, name, traj::NamedTrajectory)Construct a KnotPointObjective that penalizes leakage of name at the knot points specified by times at any indices that are outside the computational subspace.
QuantumCollocation.QuantumObjectives.coherent_ket_fidelity — Method
coherent_ket_fidelity(ψ̃s, ψ_goals)Compute coherent fidelity across multiple ket states:
F_coherent = |1/n ∑ᵢ ⟨ψᵢ_goal|ψᵢ⟩|²This requires all overlaps to have consistent phases (global phase alignment), which is necessary for implementing gates via state transfer.
Arguments
ψ̃s::Vector{<:AbstractVector}: List of isomorphic state vectorsψ_goals::Vector{<:AbstractVector{<:Complex}}: List of goal states
Quantum Constraints
QuantumCollocation.QuantumConstraints.FinalCoherentKetFidelityConstraint — Method
FinalCoherentKetFidelityConstraint(ψ_goals, ψ̃_names, final_fidelity, traj)Create a final fidelity constraint using coherent ket fidelity across multiple states.
Coherent fidelity: F = |1/n ∑ᵢ ⟨ψᵢ_goal|ψᵢ⟩|²
This constraint enforces that all state overlaps have aligned phases, which is essential when implementing a gate via multiple state transfers (e.g., MultiKetTrajectory).
Arguments
ψ_goals::Vector{<:AbstractVector{<:Complex}}: Target ket statesψ̃_names::Vector{Symbol}: Names of isomorphic state variables in trajectoryfinal_fidelity::Float64: Minimum fidelity threshold (constraint: F ≥ final_fidelity)traj::NamedTrajectory: The trajectory
Example
# For implementing X gate via |0⟩→|1⟩ and |1⟩→|0⟩
goals = [ComplexF64[0, 1], ComplexF64[1, 0]]
names = [:ψ̃1, :ψ̃2]
constraint = FinalCoherentKetFidelityConstraint(goals, names, 0.99, traj)QuantumCollocation.QuantumConstraints.LeakageConstraint — Method
LeakageConstraint(value, indices, name, traj::NamedTrajectory)Construct a KnotPointConstraint that bounds leakage of name at the knot points specified by times at any indices that are outside the computational subspace.
Quantum Integrators
DirectTrajOpt.Integrators.BilinearIntegrator — Method
BilinearIntegrator(qtraj::DensityTrajectory, N::Int)Create a BilinearIntegrator for density matrix evolution.
DirectTrajOpt.Integrators.BilinearIntegrator — Method
BilinearIntegrator(qtraj::KetTrajectory, N::Int)Create a BilinearIntegrator for ket evolution.
DirectTrajOpt.Integrators.BilinearIntegrator — Method
BilinearIntegrator(qtraj::MultiKetTrajectory, N::Int)Create a vector of BilinearIntegrators for each ket in an MultiKetTrajectory.
DirectTrajOpt.Integrators.BilinearIntegrator — Method
BilinearIntegrator(qtraj::SamplingTrajectory, N::Int)Create a vector of BilinearIntegrators for each system in a SamplingTrajectory.
Each system in the sampling ensemble gets its own dynamics integrator, but they all share the same control variables.
Returns
Vector{BilinearIntegrator}: One integrator per system in the ensemble
DirectTrajOpt.Integrators.BilinearIntegrator — Method
BilinearIntegrator(qtraj::UnitaryTrajectory, N::Int)Create a BilinearIntegrator for unitary evolution.
Options
QuantumCollocation.Options.PiccoloOptions — Type
PiccoloOptionsOptions for the Piccolo quantum optimal control library.
Fields
verbose::Bool = true: Print verbose outputtimesteps_all_equal::Bool = true: Use equal timestepsrollout_integrator::Function = expv: Integrator to use for rolloutgeodesic = true: Use the geodesic to initialize the optimization.zero_initial_and_final_derivative::Bool=false: Zero the initial and final control pulse derivatives.complex_control_norm_constraint_name::Union{Nothing, Symbol} = nothing: Name of the complex control norm constraint.complex_control_norm_constraint_radius::Float64 = 1.0: Radius of the complex control norm constraint.bound_state::Bool = false: Bound the state variables <= 1.0.leakage_constraint::Bool = false: Suppress leakage with constraint and cost.leakage_constraint_value::Float64 = 1e-2: Value for the leakage constraint.leakage_cost::Float64 = 1e-2: Leakage suppression parameter.
Control Problems
QuantumCollocation.QuantumControlProblems.QuantumControlProblem — Type
QuantumControlProblem{QT<:AbstractQuantumTrajectory}Wrapper combining quantum trajectory information with trajectory optimization problem.
This type enables:
- Type-stable dispatch on quantum trajectory type (Unitary, Ket, Density)
- Clean separation of quantum information (system, goal) from optimization details
- Composable problem transformations (e.g., SmoothPulseProblem → MinimumTimeProblem)
Fields
qtraj::QT: Quantum trajectory containing system, goal, and quantum state informationprob::DirectTrajOptProblem: Direct trajectory optimization problem with objective, dynamics, constraints
Construction
Typically created via problem templates:
qtraj = UnitaryTrajectory(sys, U_goal, N)
qcp = SmoothPulseProblem(qtraj; Q=100.0, R=1e-2)Accessors
get_trajectory(qcp): Get the NamedTrajectoryget_system(qcp): Get the QuantumSystemget_goal(qcp): Get the goal state/unitarystate_name(qcp): Get the state variable namedrive_name(qcp): Get the control variable name
Solving
solve!(qcp; max_iter=100, verbose=true)DirectTrajOpt.Solvers.solve! — Method
solve!(qcp::QuantumControlProblem; sync::Bool=true, kwargs...)Solve the quantum control problem by forwarding to the inner DirectTrajOptProblem.
Arguments
sync::Bool=true: If true, callsync_trajectory!after solving to updateqtraj.trajectorywith physical control values. Set to false to skip synchronization (e.g., for debugging).
All other keyword arguments are passed to the DirectTrajOpt solver.
PiccoloQuantumObjects.Pulses.drive_name — Method
drive_name(qcp::QuantumControlProblem)Get the control variable name from the quantum trajectory.
PiccoloQuantumObjects.QuantumTrajectories.get_goal — Method
get_goal(qcp::QuantumControlProblem)Get the goal state/operator from the quantum trajectory.
PiccoloQuantumObjects.QuantumTrajectories.get_system — Method
get_system(qcp::QuantumControlProblem)Get the QuantumSystem from the quantum trajectory.
PiccoloQuantumObjects.QuantumTrajectories.state_name — Method
state_name(qcp::QuantumControlProblem)Get the state variable name from the quantum trajectory.
PiccoloQuantumObjects.Rollouts.fidelity — Method
fidelity(qcp::QuantumControlProblem; kwargs...)Compute the fidelity of the quantum trajectory.
This is a convenience wrapper that forwards to fidelity(qcp.qtraj; kwargs...).
Example
solve!(qcp)
fid = fidelity(qcp) # Equivalent to fidelity(qcp.qtraj)QuantumCollocation.QuantumControlProblems.get_trajectory — Method
get_trajectory(qcp::QuantumControlProblem)Get the NamedTrajectory from the optimization problem.
QuantumCollocation.QuantumControlProblems.sync_trajectory! — Method
sync_trajectory!(qcp::QuantumControlProblem)Update the quantum trajectory in-place from the optimized control values.
After optimization, this function:
- Extracts the optimized controls from
prob.trajectory(unadapting if needed) - Creates a new pulse with those controls via
extract_pulse - Re-solves the ODE to get the updated quantum evolution
- Replaces
qtrajwith the new quantum trajectory
This gives you access to the continuous-time ODE solution with the optimized controls, allowing you to:
- Evaluate the fidelity via
fidelity(qcp.qtraj) - Sample the quantum state at any time via
qcp.qtraj(t) - Get the optimized pulse via
get_pulse(qcp.qtraj)
Example
solve!(qcp; max_iter=100) # Automatically calls sync_trajectory!
fid = fidelity(qcp.qtraj) # Evaluate fidelity with continuous-time solution
pulse = get_pulse(qcp.qtraj) # Get the optimized pulse