|
| 1 | +--- |
| 2 | +jupyter: |
| 3 | + jupytext: |
| 4 | + text_representation: |
| 5 | + extension: .md |
| 6 | + format_name: markdown |
| 7 | + format_version: '1.3' |
| 8 | + jupytext_version: 1.16.4 |
| 9 | + kernelspec: |
| 10 | + display_name: Python 3 (ipykernel) |
| 11 | + language: python |
| 12 | + name: python3 |
| 13 | +--- |
| 14 | + |
| 15 | +# QuTiPv5 Paper Example: Floquet Speed Test |
| 16 | + |
| 17 | +Authors: Maximilian Meyer-Mölleringhof ( [email protected]), Marc Gali ( [email protected]), Neill Lambert ( [email protected]), Paul Menczel ( [email protected]) |
| 18 | + |
| 19 | +## Introduction |
| 20 | + |
| 21 | +In this example we will discuss Hamiltonians with periodic time-dependence and how we can solve the time evolution of the system in QuTiP. |
| 22 | +A natural approach to this is using the Floquet theorem [\[1, 2\]](#References). |
| 23 | +Similar to the Bloch theorem (for spatial periodicity), this method helps to drastically simplify the problem. |
| 24 | +Let $H$ be a time-dependent and periodic Hamiltonian with period $T$, such that |
| 25 | + |
| 26 | +$H(t) = H(t + T)$ |
| 27 | + |
| 28 | +Following Floquet's theorem, there exist solutions to the Schrödinger equation that take the form |
| 29 | + |
| 30 | +$\ket{\psi_\alpha (t)} = \exp(-i \epsilon_\alpha t / \hbar) \ket{\Psi_\alpha (t)}$, |
| 31 | + |
| 32 | +where $\ket{\Psi_\alpha (t)} = \ket{\Psi_\alpha (t + T)}$ are the Floquet modes and $\epsilon_\alpha$ are the quasi-energies. |
| 33 | +Any solution of the Schrödinger equation can then be written as a linear combination |
| 34 | + |
| 35 | +$\ket{\psi(t)} = \sum_\alpha c_\alpha \ket{\psi_\alpha(t)}$, |
| 36 | + |
| 37 | +where the constants $c_\alpha$ are determined by the initial conditions. |
| 38 | + |
| 39 | +Inserting the Floquet states above into the Schrödinger equation, we can define the Floquet Hamiltonian, |
| 40 | + |
| 41 | +$H_F(t) = H(t) - i \hbar \dfrac{\partial}{\partial t}$, |
| 42 | + |
| 43 | +which renders the problem time-independent, |
| 44 | + |
| 45 | +$H_F(t) \ket{\Psi_\alpha(t)} = \epsilon_\alpha \ket{\Psi_\alpha(t)}$. |
| 46 | + |
| 47 | +Solving this equation then allows us to find $\ket{\psi (t)}$ at any large $t$. |
| 48 | + |
| 49 | +Our goal in this notebook is to show how this procedure is done in QuTiP and how it compares to the standard `sesolve` method - specifically in terms of computation time. |
| 50 | + |
| 51 | +```python |
| 52 | +import time |
| 53 | + |
| 54 | +import matplotlib.pyplot as plt |
| 55 | +import numpy as np |
| 56 | +from qutip import (CoreOptions, FloquetBasis, about, basis, expect, qeye, |
| 57 | + qzero_like, sesolve, sigmax, sigmay, sigmaz, tensor) |
| 58 | + |
| 59 | +%matplotlib inline |
| 60 | +``` |
| 61 | + |
| 62 | +## Two-level system with periodic drive |
| 63 | + |
| 64 | +We consider a two-level system that is exposed to a periodic drive: |
| 65 | + |
| 66 | +$H(t) = - \dfrac{\epsilon}{2} \sigma_z - \dfrac{\Delta}{2} \sigma_x + \dfrac{A}{2} \sin(\omega_d t) \sigma_x$, |
| 67 | + |
| 68 | +where $\epsilon$ is the energy splitting, $\Delta$ the coupling strength and $A$ the drive amplitude. |
| 69 | + |
| 70 | +```python |
| 71 | +delta = 0.2 * 2 * np.pi |
| 72 | +epsilon = 1.0 * 2 * np.pi |
| 73 | +A = 2.5 * 2 * np.pi |
| 74 | +omega = 1.0 * 2 * np.pi |
| 75 | +T = 2 * np.pi / omega |
| 76 | + |
| 77 | +H0 = -1 / 2 * (epsilon * sigmaz() + delta * sigmax()) |
| 78 | +H1 = A / 2 * sigmax() |
| 79 | +args = {"w": omega} |
| 80 | +H = [H0, [H1, lambda t, w: np.sin(w * t)]] |
| 81 | +``` |
| 82 | + |
| 83 | +```python |
| 84 | +psi0 = basis(2, 0) # initial condition |
| 85 | + |
| 86 | +# Numerical parameters |
| 87 | +dt = 0.01 |
| 88 | +N_T = 10 # Number of periods |
| 89 | +tlist = np.arange(0.0, N_T * T, dt) |
| 90 | +``` |
| 91 | + |
| 92 | +```python |
| 93 | +computational_time_floquet = np.zeros_like(tlist) |
| 94 | +computational_time_sesolve = np.zeros_like(tlist) |
| 95 | + |
| 96 | +expect_floquet = np.zeros_like(tlist) |
| 97 | +expect_sesolve = np.zeros_like(tlist) |
| 98 | +``` |
| 99 | + |
| 100 | +```python |
| 101 | +# Running the simulation |
| 102 | +for n, t in enumerate(tlist): |
| 103 | + # Floquet basis |
| 104 | + # -------------------------------- |
| 105 | + tic_f = time.perf_counter() |
| 106 | + floquetbasis = FloquetBasis(H, T, args) |
| 107 | + # Decomposing inital state into Floquet modes |
| 108 | + f_coeff = floquetbasis.to_floquet_basis(psi0) |
| 109 | + # Obtain evolved state in the original basis |
| 110 | + psi_t = floquetbasis.from_floquet_basis(f_coeff, t) |
| 111 | + p_ex = expect(sigmaz(), psi_t) |
| 112 | + toc_f = time.perf_counter() |
| 113 | + |
| 114 | + # Saving data |
| 115 | + computational_time_floquet[n] = toc_f - tic_f |
| 116 | + expect_floquet[n] = p_ex |
| 117 | + |
| 118 | + # sesolve |
| 119 | + # -------------------------------- |
| 120 | + tic_f = time.perf_counter() |
| 121 | + psi_se = sesolve(H, psi0, tlist[: n + 1], e_ops=[sigmaz()], args=args) |
| 122 | + p_ex_ref = psi_se.expect[0] |
| 123 | + toc_f = time.perf_counter() |
| 124 | + |
| 125 | + # Saving data |
| 126 | + computational_time_sesolve[n] = toc_f - tic_f |
| 127 | + expect_sesolve[n] = p_ex_ref[-1] |
| 128 | +``` |
| 129 | + |
| 130 | +```python |
| 131 | +fig, axs = plt.subplots(1, 2, figsize=(13.6, 4.5)) |
| 132 | +axs[0].plot(tlist / T, computational_time_floquet, "-") |
| 133 | +axs[0].plot(tlist / T, computational_time_sesolve, "--") |
| 134 | + |
| 135 | +axs[0].set_yscale("log") |
| 136 | +axs[0].set_yticks(np.logspace(-3, -1, 3)) |
| 137 | +axs[1].plot(tlist / T, np.real(expect_floquet), "-") |
| 138 | +axs[1].plot(tlist / T, np.real(expect_sesolve), "--") |
| 139 | + |
| 140 | +axs[0].set_xlabel(r"$t \, / \, T$") |
| 141 | +axs[1].set_xlabel(r"$t \, / \, T$") |
| 142 | +axs[0].set_ylabel(r"Computational Time [$s$]") |
| 143 | +axs[1].set_ylabel(r"$\langle \sigma_z \rangle$") |
| 144 | +axs[0].legend(("Floquet", "sesolve")) |
| 145 | + |
| 146 | +xticks = np.rint(np.linspace(0, N_T, N_T + 1, endpoint=True)) |
| 147 | +axs[0].set_xticks(xticks) |
| 148 | +axs[0].set_xlim([0, N_T]) |
| 149 | +axs[1].set_xticks(xticks) |
| 150 | +axs[1].set_xlim([0, N_T]) |
| 151 | + |
| 152 | +plt.show() |
| 153 | +``` |
| 154 | + |
| 155 | +## One Dimensional Ising Chain |
| 156 | + |
| 157 | +As a second use case for the Floquet method, we look at the one-dimensional Ising chain with a periodic drive. |
| 158 | +We describe such a system using the Hamiltonian |
| 159 | + |
| 160 | +$H(t) = g_0 \sum_{n=1}^N \sigma_z^{(n)} - J_0 \sum_{n=1}^{N - 1} \sigma_x^{(n)} \sigma_x^{(n+1)} + A \sin (\omega_d t) \sum_{n=1}^N \sigma_x^{(n)}$, |
| 161 | + |
| 162 | +where $g_0$ is the level splitting, $J_0$ is the nearest-neighbour coupling constant and $A$ is the drive amplitude. |
| 163 | + |
| 164 | +As outlined in the QuTiPv5 paper, it's educational to study the dimensional scaling of the Floquet method compared to the standard `sesolve`. |
| 165 | +Specifically how the crossing time (when the Floquet method becomes faster) depends on the dimensions $N$. |
| 166 | +Although we will not reproduce the whole plot for computation time reasons, we implement the solution for one specific choice of $N$. |
| 167 | +Feel free to play around with the parameters! |
| 168 | + |
| 169 | +```python |
| 170 | +N = 4 # number of spins |
| 171 | +g0 = 1 # energy-splitting |
| 172 | +J0 = 1.4 # coupling strength |
| 173 | +A = 1.0 # drive strength |
| 174 | +omega = 1.0 * 2 * np.pi # drive frequency |
| 175 | +T = 2 * np.pi / omega # drive period |
| 176 | +``` |
| 177 | + |
| 178 | +```python |
| 179 | +# For Hamiltonian Setup |
| 180 | +def setup_Ising_drive(N, g0, J0, A, omega, data_type="CSR"): |
| 181 | + """ |
| 182 | + # N : number of spins |
| 183 | + # g0 : splitting, |
| 184 | + # J0 : couplings |
| 185 | + # A : drive amplitude |
| 186 | + # omega: drive frequency |
| 187 | + """ |
| 188 | + with CoreOptions(default_dtype=data_type): |
| 189 | + |
| 190 | + sx_list, sy_list, sz_list = [], [], [] |
| 191 | + for i in range(N): |
| 192 | + op_list = [qeye(2)] * N |
| 193 | + op_list[i] = sigmax().to(data_type) |
| 194 | + sx_list.append(tensor(op_list)) |
| 195 | + op_list[i] = sigmay().to(data_type) |
| 196 | + sy_list.append(tensor(op_list)) |
| 197 | + op_list[i] = sigmaz().to(data_type) |
| 198 | + sz_list.append(tensor(op_list)) |
| 199 | + |
| 200 | + # Hamiltonian - Energy splitting terms |
| 201 | + H_0 = 0.0 |
| 202 | + for i in range(N): |
| 203 | + H_0 += g0 * sz_list[i] |
| 204 | + |
| 205 | + # Interaction terms |
| 206 | + H_1 = qzero_like(H_0) |
| 207 | + for n in range(N - 1): |
| 208 | + H_1 += -J0 * sx_list[n] * sx_list[n + 1] |
| 209 | + |
| 210 | + # Driving terms |
| 211 | + if A > 0: |
| 212 | + H_d = 0.0 |
| 213 | + for i in range(N): |
| 214 | + H_d += A * sx_list[i] |
| 215 | + args = {"w": omega} |
| 216 | + H = [H_0, H_1, [H_d, lambda t, w: np.sin(w * t)]] |
| 217 | + else: |
| 218 | + args = {} |
| 219 | + H = [H_0, H_1] |
| 220 | + |
| 221 | + # Defining initial conditions |
| 222 | + state_list = [basis(2, 1)] * (N - 1) |
| 223 | + state_list.append(basis(2, 0)) |
| 224 | + psi0 = tensor(state_list) |
| 225 | + |
| 226 | + # Defining expectation operator |
| 227 | + e_ops = sz_list |
| 228 | + return H, psi0, e_ops, args |
| 229 | +``` |
| 230 | + |
| 231 | +```python |
| 232 | +H, psi0, e_ops, args = setup_Ising_drive(N, g0, J0, A, omega) |
| 233 | +``` |
| 234 | + |
| 235 | +```python |
| 236 | +# Simulation parameters |
| 237 | +N_T = 10 |
| 238 | +dt = 0.01 |
| 239 | +tlist = np.arange(0, N_T * T, dt) |
| 240 | +tlist_0 = np.arange(0, T, dt) # One period tlist |
| 241 | + |
| 242 | +options = {"progress_bar": False, "store_floquet_states": True} |
| 243 | +``` |
| 244 | + |
| 245 | +```python |
| 246 | +computational_time_floquet = np.ones(tlist.shape) * np.nan |
| 247 | +computational_time_sesolve = np.ones(tlist.shape) * np.nan |
| 248 | + |
| 249 | +expect_floquet = np.zeros((N, len(tlist))) |
| 250 | +expect_sesolve = np.zeros((N, len(tlist))) |
| 251 | +``` |
| 252 | + |
| 253 | +```python |
| 254 | +# Running the simulation |
| 255 | +for n, t in enumerate(tlist): |
| 256 | + # Floquet basis |
| 257 | + # -------------------------------- |
| 258 | + tic_f = time.perf_counter() |
| 259 | + if t < T: |
| 260 | + # find the floquet modes for the time-dependent hamiltonian |
| 261 | + floquetbasis = FloquetBasis(H, T, args) |
| 262 | + else: |
| 263 | + floquetbasis = FloquetBasis(H, T, args, precompute=tlist_0) |
| 264 | + |
| 265 | + # Decomposing inital state into Floquet modes |
| 266 | + f_coeff = floquetbasis.to_floquet_basis(psi0) |
| 267 | + # Obtain evolved state in the original basis |
| 268 | + psi_t = floquetbasis.from_floquet_basis(f_coeff, t) |
| 269 | + p_ex = expect(e_ops, psi_t) |
| 270 | + toc_f = time.perf_counter() |
| 271 | + |
| 272 | + # Saving data |
| 273 | + computational_time_floquet[n] = toc_f - tic_f |
| 274 | + |
| 275 | + # sesolve |
| 276 | + # -------------------------------- |
| 277 | + tic_f = time.perf_counter() |
| 278 | + output = sesolve(H, psi0, tlist[: n + 1], e_ops=e_ops, args=args) |
| 279 | + p_ex_r = output.expect |
| 280 | + toc_f = time.perf_counter() |
| 281 | + |
| 282 | + # Saving data |
| 283 | + computational_time_sesolve[n] = toc_f - tic_f |
| 284 | + |
| 285 | + for i in range(N): |
| 286 | + expect_floquet[i, n] = p_ex[i] |
| 287 | + expect_sesolve[i, n] = p_ex_r[i][-1] |
| 288 | +``` |
| 289 | + |
| 290 | +```python |
| 291 | +fig, axs = plt.subplots(1, 2, figsize=(12, 5)) |
| 292 | +ax = axs[0] |
| 293 | +axs[0].plot(tlist / T, computational_time_floquet, "-") |
| 294 | +axs[0].plot(tlist / T, computational_time_sesolve, "--") |
| 295 | + |
| 296 | +axs[0].set_yscale("log") |
| 297 | +axs[0].set_yticks(np.logspace(-3, -1, 3)) |
| 298 | +lg = [] |
| 299 | +for i in range(N): |
| 300 | + axs[1].plot(tlist / T, np.real(expect_floquet[i, :]), "-") |
| 301 | + lg.append(f"n={i+1}") |
| 302 | +for i in range(N): |
| 303 | + axs[1].plot(tlist / T, np.real(expect_sesolve[i, :]), "--", color="black") |
| 304 | +lg.append("sesolve") |
| 305 | + |
| 306 | +axs[0].set_xlabel(r"$t \, / \, T$") |
| 307 | +axs[1].set_xlabel(r"$t \, / \, T$") |
| 308 | +axs[0].set_ylabel(r"$Computational \; Time,\; [s]$") |
| 309 | +axs[1].set_ylabel(r"$\langle \sigma_z^{{(n)}} \rangle$") |
| 310 | +axs[0].legend(("Floquet", "sesolve")) |
| 311 | +axs[1].legend(lg) |
| 312 | +xticks = np.rint(np.linspace(0, N_T, N_T + 1, endpoint=True)) |
| 313 | +axs[0].set_xticks(xticks) |
| 314 | +axs[0].set_xlim([0, N_T]) |
| 315 | +axs[1].set_xticks(xticks) |
| 316 | +axs[1].set_xlim([0, N_T]) |
| 317 | +txt = f"$N={N}$, $g_0={g0}$, $J_0={J0}$, $A={A}$" |
| 318 | +axs[0].text(0.0, 1.05, txt, transform=axs[0].transAxes) |
| 319 | + |
| 320 | +plt.show() |
| 321 | +``` |
| 322 | + |
| 323 | +## References |
| 324 | + |
| 325 | +[1] [Floquet, Annales scientifiques de l’École Normale Supérieure (1883)](http://www.numdam.org/articles/10.24033/asens.220/) |
| 326 | + |
| 327 | +[2] [Shirley, Phys.Rev. (1965)](https://link.aps.org/doi/10.1103/PhysRev.138.B979) |
| 328 | + |
| 329 | +[3] [QuTiP 5: The Quantum Toolbox in Python](https://arxiv.org/abs/2412.04705) |
| 330 | + |
| 331 | + |
| 332 | +## About |
| 333 | + |
| 334 | +```python |
| 335 | +about() |
| 336 | +``` |
| 337 | + |
| 338 | +## Testing |
| 339 | + |
| 340 | +```python |
| 341 | +for i in range(N): |
| 342 | + assert np.allclose( |
| 343 | + np.real(expect_floquet[i, :]), np.real(expect_sesolve[i, :]), atol=1e-5 |
| 344 | + ), f"floquet and sesolve solutions for Ising chain element {i} deviate." |
| 345 | +``` |
0 commit comments