|
| 1 | +--- |
| 2 | +jupyter: |
| 3 | + jupytext: |
| 4 | + text_representation: |
| 5 | + extension: .md |
| 6 | + format_name: markdown |
| 7 | + format_version: '1.3' |
| 8 | + jupytext_version: 1.13.8 |
| 9 | + kernelspec: |
| 10 | + display_name: qutip-tutorials |
| 11 | + language: python |
| 12 | + name: python3 |
| 13 | +--- |
| 14 | + |
| 15 | +<!-- #region --> |
| 16 | +# Excitation-number-restricted states: Jaynes-Cummings Chain |
| 17 | + |
| 18 | +Authors: Robert Johansson ( [email protected]), Neill Lambert ( [email protected]), Maximilian Meyer-Mölleringhof ( [email protected]) |
| 19 | + |
| 20 | +## Introduction |
| 21 | + |
| 22 | +The ENR functions construct a basis set for multipartite systems which contains only states that have an overall number of excitations. |
| 23 | +This is particularly useful for systems where the model conserves excitation number, as in the JC-chain example below. |
| 24 | + |
| 25 | +We can see this by considering a system consisting of 4 modes, each with 5 states. |
| 26 | +The total hilbert space size is $5^4 = 625$. |
| 27 | +If we are only interested in states that contain up to 2 excitations, we only need to include states such as |
| 28 | + |
| 29 | + |
| 30 | + (0, 0, 0, 0) |
| 31 | + (0, 0, 0, 1) |
| 32 | + (0, 0, 0, 2) |
| 33 | + (0, 0, 1, 0) |
| 34 | + (0, 0, 1, 1) |
| 35 | + (0, 0, 2, 0) |
| 36 | + ... |
| 37 | + |
| 38 | +The ENR fucntions create operators and states for the 4 modes that act within this state space. |
| 39 | +For example, |
| 40 | + |
| 41 | +```python |
| 42 | +a1, a2, a3, a4 = enr_destroy([5, 5, 5, 5], excitations=2) |
| 43 | +``` |
| 44 | + |
| 45 | +creates destruction operators for each mode. |
| 46 | +From this point onwards, the annihiltion operators a1, ..., a4 can be used to setup a Hamiltonian, collapse operators and expectation-value operators, etc., following the usual patterne. |
| 47 | + |
| 48 | +In this example we outline the advantage of ENR states by comparing them with the regular qutip implementation. |
| 49 | +For this we calculate the time evolution and the partial trace for each and see consistent results with notable performance improvements. |
| 50 | + |
| 51 | +#### Be aware! |
| 52 | + |
| 53 | +Many default functions in QuTiP will fail on states and operators constructed with this method. |
| 54 | +Additionally, using this formalism, annihilation and creation operators of different sub-systems no longer commute. |
| 55 | +Therefore, when constructing Hamiltonians, annihilation operators must be on the right and creation operators on the left (see the offical publication for QuTiP v5 for more info). |
| 56 | +To find all available functions to work with ENR states see [Energy Restricted Operators in the official documentation](https://qutip.readthedocs.io/en/qutip-5.0.x/apidoc/functions.html#module-qutip.core.energy_restricted). |
| 57 | +<!-- #endregion --> |
| 58 | + |
| 59 | +```python |
| 60 | +import numpy as np |
| 61 | +from qutip import (Options, Qobj, about, basis, destroy, enr_destroy, enr_fock, |
| 62 | + enr_state_dictionaries, identity, liouvillian_ref, mesolve, |
| 63 | + plot_expectation_values, tensor) |
| 64 | + |
| 65 | +%matplotlib inline |
| 66 | +``` |
| 67 | + |
| 68 | +## The Jaynes-Cumming Chain |
| 69 | + |
| 70 | +The general Jaynes-Cumming model describes a single two-level atom interacting with a single electromagnetic cavity mode. |
| 71 | +For this example, we put multiple of these systems in a chain and let them interact with neighbouring systems via their cavities. |
| 72 | +We use $a_i$ ($a^\dag_i$) as annihilation (creation) operators for the cavity $i$ and $s_i$ ($s^\dag_i$) for the atoms. |
| 73 | +We then model the complete Hamiltonian by splitting it into the individual systems: |
| 74 | + |
| 75 | +$H_0 = \sum_{i=0}^{N} a_i^\dag a_i + s_i^\dag s_i$, |
| 76 | + |
| 77 | +the atom-cavity interactions: |
| 78 | + |
| 79 | +$H_{int,AC} = \sum_{i=0}^{N} = \frac{1}{2} (a_i^\dag s_i + s_i^\dag a_i)$, |
| 80 | + |
| 81 | +and the cavity-cavity interactions: |
| 82 | + |
| 83 | +$H_{int,CC} = \sum_{i=0}^{N-1} 0.9 \cdot (a_i^\dag a_{i+1} + a_{i+1}^\dag a_{i})$, |
| 84 | + |
| 85 | +where the interaction strength of $0.9$ was chosen arbitrarily. |
| 86 | + |
| 87 | + |
| 88 | +### Problem paramters |
| 89 | + |
| 90 | +```python |
| 91 | +N = 4 # number of systems |
| 92 | +M = 2 # number of cavity states |
| 93 | +dims = [M, 2] * N # dimensions of JC spin chain |
| 94 | +excite = 1 # total number of excitations |
| 95 | +init_excite = 1 # initial number of excitations |
| 96 | +``` |
| 97 | + |
| 98 | +### Setup to Calculate Time Evolution |
| 99 | + |
| 100 | +```python |
| 101 | +def solve(d, psi0): |
| 102 | + # annihilation operators for cavity modes |
| 103 | + a = d[::2] |
| 104 | + # atomic annihilation operators |
| 105 | + sm = d[1::2] |
| 106 | + |
| 107 | + # notice the ordering of annihilation and creation operators |
| 108 | + H0 = sum([aa.dag() * aa for aa in a]) + sum([s.dag() * s for s in sm]) |
| 109 | + |
| 110 | + # atom-cavity couplings |
| 111 | + Hint_ac = 0 |
| 112 | + for n in range(N): |
| 113 | + Hint_ac += 0.5 * (a[n].dag() * sm[n] + sm[n].dag() * a[n]) |
| 114 | + |
| 115 | + # cavity-cavity couplings |
| 116 | + Hint_cc = 0 |
| 117 | + for n in range(N - 1): |
| 118 | + Hint_cc += 0.9 * (a[n].dag() * a[n + 1] + a[n + 1].dag() * a[n]) |
| 119 | + |
| 120 | + H = H0 + Hint_ac + Hint_cc |
| 121 | + |
| 122 | + e_ops = [x.dag() * x for x in d] |
| 123 | + c_ops = [0.01 * x for x in a] |
| 124 | + |
| 125 | + times = np.linspace(0, 250, 1000) |
| 126 | + L = liouvillian_ref(H, c_ops) |
| 127 | + opt = Options(nsteps=5000, store_states=True) |
| 128 | + result = mesolve(H, psi0, times, c_ops, e_ops, options=opt) |
| 129 | + return result, H, L |
| 130 | +``` |
| 131 | + |
| 132 | +### Regular QuTiP States and Operators |
| 133 | + |
| 134 | +```python |
| 135 | +d = [ |
| 136 | + tensor( |
| 137 | + [ |
| 138 | + destroy(dim1) if idx1 == idx2 else identity(dim1) |
| 139 | + for idx1, dim1 in enumerate(dims) |
| 140 | + ] |
| 141 | + ) |
| 142 | + for idx2, _ in enumerate(dims) |
| 143 | +] |
| 144 | +psi0 = tensor( |
| 145 | + [ |
| 146 | + basis(dim, init_excite) if idx == 1 else basis(dim, 0) |
| 147 | + for idx, dim in enumerate(dims) |
| 148 | + ] |
| 149 | +) |
| 150 | +``` |
| 151 | + |
| 152 | +Regular operators of different systems commute as they belong to different Hilbert spaces. |
| 153 | +Example: |
| 154 | + |
| 155 | +```python |
| 156 | +d[0].dag() * d[1] == d[1] * d[0].dag() |
| 157 | +``` |
| 158 | + |
| 159 | +Solving the time evolution: |
| 160 | + |
| 161 | +```python |
| 162 | +res1, H1, L1 = solve(d, psi0) |
| 163 | +``` |
| 164 | + |
| 165 | +### Using ENR States and Operators |
| 166 | + |
| 167 | +```python |
| 168 | +d_enr = enr_destroy(dims, excite) |
| 169 | +init_enr = [init_excite if n == 1 else 0 for n in range(2 * N)] |
| 170 | +psi0_enr = enr_fock(dims, excite, init_enr) |
| 171 | +``` |
| 172 | + |
| 173 | +Using ENR states forces us to give up on the standard tensor structure of multiple Hilbert spaces. |
| 174 | +Operators for different systems therefore generally no longer commute: |
| 175 | + |
| 176 | +```python |
| 177 | +d_enr[0].dag() * d_enr[1] == d_enr[1] * d_enr[0].dag() |
| 178 | +``` |
| 179 | + |
| 180 | +Solving the time evolution: |
| 181 | + |
| 182 | +```python |
| 183 | +res2, H2, L2 = solve(d_enr, psi0_enr) |
| 184 | +``` |
| 185 | + |
| 186 | +### Comparison of Expectation Values |
| 187 | + |
| 188 | +```python |
| 189 | +fig, axes = plot_expectation_values([res1, res2], figsize=(10, 8)) |
| 190 | +for idx, ax in enumerate(axes[:, 0]): |
| 191 | + if idx % 2: |
| 192 | + ax.set_ylabel(f"Atom {idx//2}") |
| 193 | + else: |
| 194 | + ax.set_ylabel(f"Cavity {idx//2}") |
| 195 | + ax.set_ylim(-0.1, 1.1) |
| 196 | + ax.grid() |
| 197 | +fig.tight_layout() |
| 198 | +``` |
| 199 | + |
| 200 | +### Calculation of Partial Trace |
| 201 | + |
| 202 | +The usage of ENR states makes many standard QuTiP features fail. |
| 203 | +*ptrace* is one of those. |
| 204 | +Below we demonstrate how the partial trace for ENR states can be calculated and show the corresponding result together with the standrad QuTiP approach. |
| 205 | + |
| 206 | +```python |
| 207 | +def ENR_ptrace(rho, sel, excitations): |
| 208 | + if isinstance(sel, int): |
| 209 | + sel = np.array([sel]) |
| 210 | + else: |
| 211 | + sel = np.asarray(sel) |
| 212 | + |
| 213 | + if (sel < 0).any() or (sel >= len(rho.dims[0])).any(): |
| 214 | + raise TypeError("Invalid selection index in ptrace.") |
| 215 | + |
| 216 | + drho = rho.dims[0] |
| 217 | + _, state2idx, idx2state = enr_state_dictionaries(drho, excitations) |
| 218 | + |
| 219 | + dims_short = np.asarray(drho).take(sel).tolist() |
| 220 | + nstates2, state2idx2, _ = enr_state_dictionaries(dims_short, excitations) |
| 221 | + |
| 222 | + # construct new density matrix |
| 223 | + rhout = np.zeros((nstates2, nstates2), dtype=np.complex64) |
| 224 | + # dimensions of traced out system |
| 225 | + rest = np.setdiff1d(np.arange(len(drho)), sel) |
| 226 | + for state in idx2state: |
| 227 | + for state2 in idx2state: |
| 228 | + # add diagonal elements to the new density matrix |
| 229 | + state_red = np.asarray(state).take(rest) |
| 230 | + state2_red = np.asarray(state2).take(rest) |
| 231 | + if np.all(state_red == state2_red): |
| 232 | + rhout[ |
| 233 | + state2idx2[tuple(np.asarray(state).take(sel))], |
| 234 | + state2idx2[tuple(np.asarray(state2).take(sel))], |
| 235 | + ] += rho.data[state2idx[state], state2idx[state2]] |
| 236 | + |
| 237 | + new_dims = np.asarray(drho).take(sel).tolist() |
| 238 | + return Qobj(rhout, dims=[new_dims, new_dims], shape=[nstates2, nstates2]) |
| 239 | +``` |
| 240 | + |
| 241 | +```python |
| 242 | +res1.states[10].ptrace(1) |
| 243 | +``` |
| 244 | + |
| 245 | +```python |
| 246 | +ENR_ptrace(res2.states[10], 1, excite) |
| 247 | +``` |
| 248 | + |
| 249 | +```python |
| 250 | +res1.states[10].ptrace([0, 1, 4]) |
| 251 | +``` |
| 252 | + |
| 253 | +```python |
| 254 | +ENR_ptrace(res2.states[10], [0, 1, 4], excite) |
| 255 | +``` |
| 256 | + |
| 257 | +## About |
| 258 | + |
| 259 | +```python |
| 260 | +about() |
| 261 | +``` |
| 262 | + |
| 263 | +## Testing |
| 264 | + |
| 265 | +```python |
| 266 | +assert np.allclose( |
| 267 | + res1.states[10].ptrace([1]), ENR_ptrace(res2.states[10], [1], excite) |
| 268 | +), "The approaches do not yield the same result." |
| 269 | +``` |
0 commit comments