Skip to content
Merged
Show file tree
Hide file tree
Changes from 1 commit
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
Prev Previous commit
Next Next commit
Tests for free phases
  • Loading branch information
andgoldschmidt committed Nov 9, 2024
commit f7605bfa491bb6186dc9666914bd6dd6ea7bebbe
30 changes: 29 additions & 1 deletion src/isomorphisms.jl
Original file line number Diff line number Diff line change
Expand Up @@ -76,7 +76,7 @@ function iso_vec_to_operator(Ũ⃗::AbstractVector{R}) where R
N = Int(sqrt(Ũ⃗_dim))
U = Matrix{complex(eltype(Ũ⃗))}(undef, N, N)
for i=0:N-1
U[:, i+1] .= @view(Ũ⃗[i * 2N .+ (1:N)]) + one(R) * im * @view(Ũ⃗[i * 2N .+ (N+1:2N)])
U[:, i+1] .= @view(Ũ⃗[i * 2N .+ (1:N)]) + one(eltype(Ũ⃗)) * im * @view(Ũ⃗[i * 2N .+ (N+1:2N)])
end
return U
end
Expand Down Expand Up @@ -228,4 +228,32 @@ iso_dm(ρ::AbstractMatrix) = ket_to_iso(vec(ρ))
@test iso_operator_to_iso_vec(iso_vec_to_iso_operator(iso_vec)) == iso_vec
end

@testitem "Test isomorphism type promotion" begin
# Check that generic types are converted to Complex
X = Any[0 1; 1 0]
@test operator_to_iso_vec(X) == [0; 1; 0; 0; 1; 0; 0; 0]

X = Number[0 1; 1 0]
@test operator_to_iso_vec(X) == [0; 1; 0; 0; 1; 0; 0; 0]


# Check that Any types are converted to Real
x = Any[0; 1; 0; 0; 1; 0; 0; 0]
@test iso_vec_to_operator(x) == [0 1; 1 0]

# Check that valid Complex can be converted to Real
x = Complex[0; 1; 0; 0; 1; 0; 0; 0]
@test iso_vec_to_operator(x) == [0 1; 1 0]

# Check that actual Complex throws an error
x = Complex[0; im; 0; 0; 1; 0; 0; 0]
@test_throws TypeError iso_vec_to_operator(x)

# Check the non-converting Isomorphisms
x = Any[0; 1; 0; 0; 1; 0; 0; 0]
X = Any[0 1 0 0; 1 0 0 0; 0 0 0 1; 0 0 1 0]
@test iso_vec_to_iso_operator(x) == X
@test iso_operator_to_iso_vec(X) == x
end

end
36 changes: 36 additions & 0 deletions src/losses/unitary_infidelity_loss.jl
Original file line number Diff line number Diff line change
Expand Up @@ -419,6 +419,11 @@ end
Y = [0 -im; im 0]
@test unitary_fidelity(X, X) ≈ 1
@test unitary_fidelity(X, Y) ≈ 0

# Tr undefined on Type{Any}
X = Any[0 1; 1 0]
Y = Complex[0 -im; im 0]
@test unitary_fidelity(X, Y) ≈ 0
end

@testitem "Isovec Unitary Fidelity" begin
Expand Down Expand Up @@ -513,5 +518,36 @@ end


@testitem "Isovec Unitary Fidelity Gradient" begin
@test_skip true
end

@testitem "Free phase fidelity" begin
n_levels = 3
phase_data = [1.9, 2.7]
phase_operators = [PAULIS[:Z], PAULIS[:Z]]
subspace = get_subspace_indices([1:2, 1:2], [n_levels, n_levels])

R = Losses.free_phase(phase_data, phase_operators)
@test R'R ≈ [1 0 0 0; 0 1 0 0; 0 0 1 0; 0 0 0 1]
@test size(R) == (2^2, 2^2)

U_goal = EmbeddedOperator(GATES[:CZ], subspace, [n_levels, n_levels]).operator
U_final = EmbeddedOperator(R'GATES[:CZ], subspace, [n_levels, n_levels]).operator
# Value is ~0.3 without phases
@test unitary_fidelity(U_final, U_goal, subspace=subspace) < 0.5
@test unitary_free_phase_fidelity(
U_final,
U_goal,
phase_data,
phase_operators,
subspace=subspace
) ≈ 1

# Forgot subspace
@test_throws DimensionMismatch iso_vec_unitary_free_phase_fidelity(
operator_to_iso_vec(U_final),
operator_to_iso_vec(U_goal),
phase_data,
phase_operators,
)
end
30 changes: 28 additions & 2 deletions src/problem_templates/unitary_minimum_time_problem.jl
Original file line number Diff line number Diff line change
Expand Up @@ -133,8 +133,8 @@ end
@testitem "Minimum time Hadamard gate" begin
using NamedTrajectories

H_drift = GATES[:Z]
H_drives = [GATES[:X], GATES[:Y]]
H_drift = PAULIS[:Z]
H_drives = [PAULIS[:X], PAULIS[:Y]]
U_goal = GATES[:H]
T = 51
Δt = 0.2
Expand Down Expand Up @@ -165,3 +165,29 @@ end
UnitaryMinimumTimeProblem(prob)

end

@testitem "Minimum time free phase" begin
using NamedTrajectories

phase_operators = [PAULIS[:Z]]

prob = UnitarySmoothPulseProblem(
QuantumSystem([PAULIS[:X]]), GATES[:H], 51, 0.2,
ipopt_options=IpoptOptions(print_level=1),
piccolo_options=PiccoloOptions(verbose=false),
phase_operators=phase_operators
)

# Soft fidelity constraint
final_fidelity = minimum([0.99, unitary_fidelity(prob)])
mintime_prob = UnitaryMinimumTimeProblem(
prob,
final_fidelity=final_fidelity,
phase_operators=phase_operators
)
solve!(mintime_prob; max_iter=100)

duration_after = sum(get_timesteps(mintime_prob.trajectory))
duration_before = sum(get_timesteps(prob.trajectory))
@test duration_after < duration_before
end
29 changes: 29 additions & 0 deletions src/problem_templates/unitary_smooth_pulse_problem.jl
Original file line number Diff line number Diff line change
Expand Up @@ -308,3 +308,32 @@ end
final = unitary_fidelity(prob, subspace=U_goal.subspace_indices)
@test final > initial
end

@testitem "Free phase Y gate using X" begin
using Random
Random.seed!(1234)
phase_name = :ϕ
phase_operators = [PAULIS[:Z]]

prob = UnitarySmoothPulseProblem(
QuantumSystem([PAULIS[:X]]), GATES[:Y], 51, 0.2;
phase_operators=phase_operators,
phase_name=phase_name,
ipopt_options=IpoptOptions(print_level=1),
piccolo_options=PiccoloOptions(verbose=false, free_time=false)
)

before = prob.trajectory.global_data[phase_name]
solve!(prob, max_iter=20)
after = prob.trajectory.global_data[phase_name]

@test before ≠ after

@test unitary_fidelity(
prob,
phases=prob.trajectory.global_data[phase_name],
phase_operators=phase_operators
) > 0.9

@test unitary_fidelity(prob) < 0.9
end
3 changes: 3 additions & 0 deletions src/rollouts.jl
Original file line number Diff line number Diff line change
Expand Up @@ -474,6 +474,9 @@ end
@test unitary_fidelity(prob) ≈ 1
@test unitary_fidelity(embedded_U_goal, as, ts, sys) ≈ 1

# Free phase unitary
@test unitary_fidelity(prob, phases=[0.0], phase_operators=[PAULIS[:Z]]) ≈ 1

# Expv explicit
# State fidelity
@test fidelity(ψ, ψ_goal, as, ts, sys, integrator=expv) ≈ 1
Expand Down