-
Notifications
You must be signed in to change notification settings - Fork 0
/
2023-05-16_map_splitting_R26μ_50pps.jl
205 lines (171 loc) · 5.23 KB
/
2023-05-16_map_splitting_R26μ_50pps.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
# -*- coding: utf-8 -*-
# ---
# jupyter:
# jupytext:
# formats: ipynb,jl:light
# text_representation:
# extension: .jl
# format_name: light
# format_version: '1.5'
# jupytext_version: 1.14.5
# kernelspec:
# display_name: Julia 1.8.5
# language: julia
# name: julia-1.8
# ---
# # Map adiabaticity
# ## Hamiltonian
using QuantumPropagators
using LinearAlgebra
using FFTW
using Serialization
using ProgressMeter
using Revise
using Plots
const 𝕚 = 1im;
const μm = 1;
const μs = 1;
const ns = 1e-3μs;
const cm = 1e4μm;
const met = 1e6μm;
const sec = 1e6μs;
const ms = 1e3μs;
const MHz = 2π;
const Dalton = 1.5746097504353806e+01;
const RUBIDIUM_MASS = 86.91Dalton;
const TAI_RADIUS = 25.46μm
const N_SITES = 8;
const OMEGA_TARGET = 50π / sec;
const EFFECTIVE_MASS = TAI_RADIUS^2 * RUBIDIUM_MASS;
includet("./include/rotating_tai.jl")
includet("./include/split_propagator.jl")
function rotating_tai_hamiltonian(;
tlist,
θ,
ω, # function of time
V₀,
Ω=0.0,
direction=1,
m=N_SITES,
mass=EFFECTIVE_MASS
)
V = Diagonal(V₀ .* cos.(m .* θ))
dθ = θ[2] - θ[1]
nθ = length(θ)
pgrid::Vector{Float64} = 2π * fftfreq(nθ, 1 / dθ)
P::Diagonal{Float64,Vector{Float64}} = Diagonal(pgrid)
K::Diagonal{Float64,Vector{Float64}} = Diagonal(pgrid .^ 2 / (2 * mass))
_Ψ = Array{ComplexF64}(undef, nθ)
fft_op = plan_fft!(_Ψ)
ifft_op = plan_ifft!(_Ψ)
transforms = (Ψ -> fft_op * Ψ, Ψ -> ifft_op * Ψ)
K′::Diagonal{Float64,Vector{Float64}} = K - Ω * P
if ω isa Number
if direction == 1
H = SplitGenerator(K′ + ω * P, V, transforms...)
elseif direction == -1
H = SplitGenerator(K′ - ω * P, V, transforms...)
else
error("direction must be ±1")
end
else
if direction == 1
H = SplitGenerator(hamiltonian(K′, (P, ω)), V, transforms...)
elseif direction == -1
H = SplitGenerator(hamiltonian(K′, (-P, ω)), V, transforms...)
else
error("direction must be ±1")
end
end
end
omega_ramp_up(t; w0=OMEGA_TARGET, t_r) = w0 * sin(π * t / (2t_r))^2;
omega_ramp_down(t; w0=OMEGA_TARGET, t_r) = w0 * cos(π * t / (2t_r))^2;
using QuantumPropagators.Controls: discretize_on_midpoints
function choose_timesteps(separation_time; timesteps_per_microsec=1, minimum_timesteps=1001)
return max(minimum_timesteps, Int(separation_time ÷ μs) * timesteps_per_microsec + 1)
end
choose_timesteps(100ms)
function propagate_splitting(
separation_time,
potential_depth;
ret=:fidelity,
timesteps_per_microsec=1,
minimum_timesteps=1001,
theta_max=0.25π,
theta_steps=1024,
kwargs...
)
nt = choose_timesteps(separation_time; timesteps_per_microsec, minimum_timesteps)
tlist = collect(range(0, separation_time, length=nt))
ω_func(t) = omega_ramp_up(t; w0=OMEGA_TARGET, t_r=separation_time)
θ::Vector{Float64} = collect(range(0, theta_max, length=theta_steps))
Ĥ = rotating_tai_hamiltonian(
tlist=tlist,
V₀=potential_depth,
θ=θ,
ω=discretize_on_midpoints(ω_func, tlist)
)
if ret == :system
return Ĥ, tlist
end
Ĥ₀ = evaluate(Ĥ, tlist, 1)
Ψ₀ = get_ground_state(Ĥ₀, θ, π / 8, d=0.05, steps=10_000)
if ret == :initial_state
return Ψ₀, θ
end
Ĥ_tgt = evaluate(Ĥ, tlist, nt - 1)
Ψ_tgt = get_ground_state(Ĥ_tgt, θ, π / 8, d=0.05, steps=10_000)
if ret == :target
return Ψ_tgt, θ
end
Ψ = propagate(Ψ₀, Ĥ, tlist; method=:splitprop, kwargs...)
if ret == :propagation
return Ψ
end
F = abs2(Ψ ⋅ Ψ_tgt)
if ret == :fidelity
return F
else
error("Invalid ret=$ret")
end
end
# ## Map
# ### V0 = 0.1 - 2.2 MHz; sep time = 10⁻¹ - 10⁵ μs
potential_depth_values = collect(range(0.1MHz, 2.2MHz, length=106))
potential_depth_values ./ MHz
separation_time_orders_of_magnitude = collect(range(-1, 5, length=121))
separation_time_values = [10^x * μs for x in separation_time_orders_of_magnitude]
Threads.nthreads()
using ProgressMeter
function map_fidelity(potential_depth_values, separation_time_values; kwargs...)
N = length(potential_depth_values)
M = length(separation_time_values)
F = zeros(N, M)
progress = Progress(N * M)
Threads.@threads for j = 1:M
for i = 1:N
t_r = separation_time_values[j]
V0 = potential_depth_values[i]
F[i, j] = propagate_splitting(t_r, V0; kwargs...)
next!(progress)
end
end
return F
end
include("./include/workflow.jl")
F = run_or_load("./data/2023-05-16_50pps_map_splitting_fidelity.npz"; force=false) do
map_fidelity(potential_depth_values, separation_time_values)
end
# #### Countour Plot
contourf(
separation_time_values ./ sec,
potential_depth_values ./ MHz,
F,
tick_direction=:out,
xminorticks=9,
yminorticks=2,
xaxis=:log10,
ylabel="V₀ (MHz)",
xlabel="separation time (seconds)",
title=raw"Sep. Fidelity $|⟨Ψ(t_r) | Ψ_{\textrm{tgt}}⟩|^2$ (R=25.46μm, ω=50π/s)",
)