Skip to content

Commit 6339ca3

Browse files
First version of nm mcsolve notebook
- replicated plots from paper - added introduction and explanations in own words - added general structure for notbeook such as about, etc.
1 parent 522e7b7 commit 6339ca3

File tree

1 file changed

+358
-0
lines changed

1 file changed

+358
-0
lines changed
Lines changed: 358 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,358 @@
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: Python 3
11+
language: python
12+
name: python3
13+
---
14+
15+
# QuTiPv5 Paper Example: Monte Carlo Solver for non-Markovian Baths
16+
17+
Authors: Maximilian Meyer-Mölleringhof ([email protected]), Neill Lambert ([email protected])
18+
19+
## Introduction
20+
21+
When a quantum system experiences non-Markovian effects, it can no longer be described using the standard Lindblad formalism.
22+
To compute the time evolution of the density matrix, we can however use time-convolutionless (TCL) projection operators that lead to a differential equation in time-local form:
23+
24+
$\dot{\rho} (t) = - \dfrac{i}{\hbar} [H(t), \rho(t)] + \sum_n \gamma_n(t) \mathcal{D}_n(t) [\rho(t)]$
25+
26+
with
27+
28+
$\mathcal{D}_n[\rho(t)] = A_n \rho(t) A^\dagger_n - \dfrac{1}{2} [A^\dagger_n A_n \rho(t) + \rho(t) A_n^\dagger A_n]$.
29+
30+
These equations include the system Hamiltonian $H(t)$ and jump operators $A_n$.
31+
Contrary to a Lindblad equation, the coupling rates $\gamma_n(t)$ may be negative here.
32+
33+
In QuTiPv5, the non-Markovian Monte Carlo solver is introduced.
34+
It enables the mapping of the general master equation given above, to a Lindblad equation on the same Hilbert space.
35+
This is achieved by the introduction of so called "influence martingale" which act as trajectory weighting [\[1, 2, 3\]](#References):
36+
37+
$\mu (t) = \exp \left[ \alpha \int_0^t s(\tau) d \tau \right] \Pi_k \dfrac{\gamma_{n_k} (t_k)}{\Gamma_{n_k} (t_k)}$.
38+
39+
Here, the product runs over all jump operators on the trajectory with jump channels $n_k$ and jump times $t_k < t$.
40+
41+
The jump operators are required to fulfill the completeness relation $\sum_n A_n^\dagger A_n = \alpha \mathbb{1}$ for $\alpha > 0$.
42+
This condition is automatically taken care off by QuTiP's `nm_mcsolve()` function.
43+
44+
To finally arrive at the Lindblad form, the shift function
45+
46+
$s(t) = 2 | \min \{ 0, \gamma_1(t), \gamma_2(t), ... \} |$,
47+
48+
such that the shifted rates $\Gamma_n (t) = \gamma_n(t) + s(t)$ non-negative.
49+
Using the regular MCWF method, we obtain the completely positive Lindblad equation
50+
51+
$\dot{\rho}'(t) = - \dfrac{i}{\hbar} [ H(t), \rho'(t) ] + \sum_n \Gamma(t) \mathcal{D}_n[\rho'(t)]$,
52+
53+
so that $\rho'(t) = \mathcal{E} \{\ket{\psi(t)} \bra{\psi(t)}\}$ for the generated trajecotries $\ket{\psi (t)}$.
54+
55+
56+
57+
58+
```python
59+
import matplotlib.pyplot as plt
60+
import numpy as np
61+
import qutip as qt
62+
from qutip import (about, basis, brmesolve, expect, lindblad_dissipator,
63+
liouvillian, mesolve, nm_mcsolve, qeye, sigmam, sigmap,
64+
sigmax, sigmay, sigmaz)
65+
from scipy.interpolate import CubicSpline
66+
from scipy.optimize import root_scalar
67+
```
68+
69+
## Damped Jaynes-Cumming Model
70+
71+
To illustrate the application of `nm_mcsolve()`, we want to look at the damped Jaynes-Cumming model.
72+
It describes a two-level atom coupled to a damped cavity mode.
73+
Such a model can be accurately model as a two-level system coupled to an environment with the power spectrum [\[4\]](#References):
74+
75+
$S(\omega) = \dfrac{\lambda \Gamma^2}{(\omega_0 - \Delta - \omega)^2 + \Gamma^2}$,
76+
77+
where $\lambda$ is the atom-cavity coupling strength, $\omega_0$ is the transition frequency of the atom, $\Delta$ the detuning of the cavity and $\Gamma$ the spectral width.
78+
At zero temperature and after performing the rotating wave approximation, the dynamics of the two-level atom can be described by the master equation
79+
80+
$\dot{\rho}(t) = \dfrac{A(t)}{2i\hbar} [ \sigma_+ \sigma_-, \rho(t) ] + \gamma (t) \mathcal{D}_- [\rho (t)]$,
81+
82+
with the state $\rho(t)$ in the interaction picture, the the ladder operators $\sigma_\pm$, the dissipator $\mathcal{D}_-$ for the Lindblad operator $\sigma_-$, and lastly $\gamma (t)$ and $A(t)$ being the real and imaginary parts of
83+
84+
$\gamma(t) + i A(t) = \dfrac{2 \lambda \Gamma \sinh (\delta t / 2)}{\delta \cosh (\delta t / 2) + (\Gamma - i \Delta) \sinh (\delta t/2)}$,
85+
86+
with $\delta = [(\Gamma - i \Delta)^2 - 2 \lambda \Gamma]^{1/2}$.
87+
88+
```python
89+
H = sigmap() * sigmam() / 2
90+
initial_state = (basis(2, 0) + basis(2, 1)).unit()
91+
tlist = np.linspace(0, 5, 500)
92+
```
93+
94+
```python
95+
# Constants
96+
gamma0 = 1
97+
lamb = 0.3 * gamma0
98+
Delta = 8 * lamb
99+
100+
# Derived Quantities
101+
delta = np.sqrt((lamb - 1j * Delta) ** 2 - 2 * gamma0 * lamb)
102+
deltaR = np.real(delta)
103+
deltaI = np.imag(delta)
104+
deltaSq = deltaR**2 + deltaI**2
105+
```
106+
107+
```python
108+
# calculate gamma and A
109+
def prefac(t):
110+
return (
111+
2
112+
* gamma0
113+
* lamb
114+
/ (
115+
(lamb**2 + Delta**2 - deltaSq) * np.cos(deltaI * t)
116+
- (lamb**2 + Delta**2 + deltaSq) * np.cosh(deltaR * t)
117+
- 2 * (Delta * deltaR + lamb * deltaI) * np.sin(deltaI * t)
118+
+ 2 * (Delta * deltaI - lamb * deltaR) * np.sinh(deltaR * t)
119+
)
120+
)
121+
122+
123+
def cgamma(t):
124+
return prefac(t) * (
125+
lamb * np.cos(deltaI * t)
126+
- lamb * np.cosh(deltaR * t)
127+
- deltaI * np.sin(deltaI * t)
128+
- deltaR * np.sinh(deltaR * t)
129+
)
130+
131+
132+
def cA(t):
133+
return prefac(t) * (
134+
Delta * np.cos(deltaI * t)
135+
- Delta * np.cosh(deltaR * t)
136+
- deltaR * np.sin(deltaI * t)
137+
+ deltaI * np.sinh(deltaR * t)
138+
)
139+
```
140+
141+
```python
142+
_gamma = np.zeros_like(tlist)
143+
_A = np.zeros_like(tlist)
144+
for i in range(len(tlist)):
145+
_gamma[i] = cgamma(tlist[i])
146+
_A[i] = cA(tlist[i])
147+
148+
gamma = CubicSpline(tlist, np.complex128(_gamma))
149+
A = CubicSpline(tlist, np.complex128(_A))
150+
```
151+
152+
```python
153+
unitary_gen = liouvillian(H)
154+
dissipator = lindblad_dissipator(sigmam())
155+
me_sol = mesolve([[unitary_gen, A], [dissipator, gamma]], initial_state, tlist)
156+
```
157+
158+
```python
159+
mc_sol = nm_mcsolve(
160+
[[H, A]],
161+
initial_state,
162+
tlist,
163+
ops_and_rates=[(sigmam(), gamma)],
164+
ntraj=1_000,
165+
options={"map": "parallel"},
166+
seeds=0,
167+
)
168+
```
169+
170+
## Comparison to other methods
171+
172+
We want to compare this computation with the HEOM and the Bloch-Refield solver.
173+
For these methods we directly apply a spin-boson model and the free reservoir auto-correlation function
174+
175+
$C(t) = \dfrac{\lambda \Gamma}{2} e^{-i (\omega - \Delta) t - \lambda |t|}$
176+
177+
as well as the same power specture as we used above.
178+
We use the Hamiltonian $H = \omega_0 \sigma_+ \sigma_-$ in the Schrödinger picture and the coupling operator $Q = \sigma_+ + \sigma_-$.
179+
Here, we chose $\omega_0 \gg \Delta$ to ensure validity of the rotating wave approximation.
180+
181+
```python
182+
omega_c = 100
183+
omega_0 = omega_c + Delta
184+
185+
H = omega_0 * qt.sigmap() * qt.sigmam()
186+
Q = qt.sigmap() + qt.sigmam()
187+
188+
189+
def power_spectrum(w):
190+
return gamma0 * lamb**2 / ((omega_c - w) ** 2 + lamb**2)
191+
```
192+
193+
Firstly, the HEOM solver can directly be applied after separating the auto-correlation function into its real and imaginary parts:
194+
195+
```python
196+
ck_real = [gamma0 * lamb / 4] * 2
197+
vk_real = [lamb - 1j * omega_c, lamb + 1j * omega_c]
198+
ck_imag = np.array([1j, -1j]) * gamma0 * lamb / 4
199+
vk_imag = vk_real
200+
```
201+
202+
```python
203+
heom_bath = qt.heom.BosonicBath(Q, ck_real, vk_real, ck_imag, vk_imag)
204+
heom_sol = qt.heom.heomsolve(H, heom_bath, 10, qt.ket2dm(initial_state), tlist)
205+
```
206+
207+
And secondly, for the Bloch-Redfield solver we can directly use the power spectrum as input:
208+
209+
```python
210+
br_sol = brmesolve(H, initial_state, tlist, a_ops=[(sigmax(), power_spectrum)])
211+
```
212+
213+
Finally, in order to compare these results with the `nm_mesolve()` method, we transform them to the interaction picture:
214+
215+
```python
216+
Us = [(-1j * H * t).expm() for t in tlist]
217+
heom_states = [U * s * U.dag() for (U, s) in zip(Us, heom_sol.states)]
218+
br_states = [U * s * U.dag() for (U, s) in zip(Us, br_sol.states)]
219+
```
220+
221+
### Plotting the Time Evolution
222+
223+
We choose to plot three comparisons to highlight the different solvers' computational results.
224+
In all plots, the gray areas show when $\gamma(t)$ is negative.
225+
226+
```python
227+
root1 = root_scalar(lambda t: cgamma(t), method="bisect", bracket=(1, 2)).root
228+
root2 = root_scalar(lambda t: cgamma(t), method="bisect", bracket=(2, 3)).root
229+
root3 = root_scalar(lambda t: cgamma(t), method="bisect", bracket=(3, 4)).root
230+
root4 = root_scalar(lambda t: cgamma(t), method="bisect", bracket=(4, 5)).root
231+
```
232+
233+
```python
234+
projector = (sigmaz() + qeye(2)) / 2
235+
236+
rho11_me = expect(projector, me_sol.states)
237+
rho11_mc = expect(projector, mc_sol.states[::10])
238+
rho11_br = expect(projector, br_sol.states)
239+
rho11_heom = expect(projector, heom_states)
240+
241+
plt.plot(tlist, rho11_me, "-", color="orange", label="mesolve")
242+
plt.plot(
243+
tlist[::10],
244+
rho11_mc,
245+
"x",
246+
color="blue",
247+
label="nm_mcsolve",
248+
)
249+
plt.plot(tlist, rho11_br, "-.", color="gray", label="brmesolve")
250+
plt.plot(tlist, rho11_heom, "--", color="green", label="heomsolve")
251+
252+
plt.xlabel(r"$t\, /\, \lambda^{-1}$")
253+
plt.xlim((-0.2, 5.2))
254+
plt.xticks([0, 2.5, 5], labels=["0", "2.5", "5"])
255+
plt.title(r"$\rho_{11}$")
256+
plt.ylim((0.4376, 0.5024))
257+
plt.yticks([0.44, 0.46, 0.48, 0.5], labels=["0.44", "0.46", "0.48", "0.50"])
258+
259+
plt.axvspan(root1, root2, color="gray", alpha=0.08, zorder=0)
260+
plt.axvspan(root3, root4, color="gray", alpha=0.08, zorder=0)
261+
262+
plt.legend()
263+
plt.show()
264+
```
265+
266+
```python
267+
me_x = expect(sigmax(), me_sol.states)
268+
mc_x = expect(sigmax(), mc_sol.states[::10])
269+
heom_x = expect(sigmax(), heom_states)
270+
br_x = expect(sigmax(), br_states)
271+
272+
me_y = expect(sigmay(), me_sol.states)
273+
mc_y = expect(sigmay(), mc_sol.states[::10])
274+
heom_y = expect(sigmay(), heom_states)
275+
br_y = expect(sigmay(), br_states)
276+
```
277+
278+
```python
279+
# We smooth the HEOM result because it oscillates quickly and gets hard to see
280+
heom_plot = heom_x * heom_x + heom_y * heom_y
281+
heom_plot = np.convolve(heom_plot, np.array([1 / 11] * 11), mode="valid")
282+
heom_tlist = tlist[5:-5]
283+
```
284+
285+
```python
286+
rho01_me = me_x * me_x + me_y * me_y
287+
rho01_mc = mc_x * mc_x + mc_y * mc_y
288+
rho01_br = br_x * br_x + br_y * br_y
289+
290+
plt.plot(tlist, rho01_me, "-", color="orange", label=r"mesolve")
291+
plt.plot(tlist[::10], rho01_mc, "x", color="blue", label=r"nm_mcsolve")
292+
plt.plot(heom_tlist, heom_plot, "--", color="green", label=r"heomsolve")
293+
plt.plot(tlist, rho01_br, "-.", color="gray", label=r"brmesolve")
294+
295+
plt.xlabel(r"$t\, /\, \lambda^{-1}$")
296+
plt.xlim((-0.2, 5.2))
297+
plt.xticks([0, 2.5, 5], labels=["0", "2.5", "5"])
298+
plt.title(r"$| \rho_{01} |^2$")
299+
plt.ylim((0.8752, 1.0048))
300+
plt.yticks(
301+
[0.88, 0.9, 0.92, 0.94, 0.96, 0.98, 1],
302+
labels=["0.88", "0.90", "0.92", "0.94", "0.96", "0.98", "1"],
303+
)
304+
305+
plt.axvspan(root1, root2, color="gray", alpha=0.08, zorder=0)
306+
plt.axvspan(root3, root4, color="gray", alpha=0.08, zorder=0)
307+
308+
plt.legend()
309+
plt.show()
310+
```
311+
312+
```python
313+
plt.plot(tlist, np.zeros_like(tlist), "-", color="orange", label=r"Zero")
314+
plt.plot(
315+
tlist[::10],
316+
1000 * (mc_sol.trace[::10] - 1),
317+
"x",
318+
color="blue",
319+
label=r"nm_mcsolve",
320+
)
321+
322+
plt.xlabel(r"$t\, /\, \lambda^{-1}$")
323+
plt.xlim((-0.2, 5.2))
324+
plt.xticks([0, 2.5, 5], labels=["0", "2.5", "5"])
325+
plt.title(r"$(\mu - 1)\, /\, 10^{-3}$")
326+
plt.ylim((-5.8, 15.8))
327+
plt.yticks([-5, 0, 5, 10, 15])
328+
329+
plt.axvspan(root1, root2, color="gray", alpha=0.08, zorder=0)
330+
plt.axvspan(root3, root4, color="gray", alpha=0.08, zorder=0)
331+
332+
plt.legend()
333+
plt.show()
334+
```
335+
336+
## References
337+
338+
\[1\] [Donvil and Muratore-Ginanneschi. Nat Commun (2022).](https://www.nature.com/articles/s41467-022-31533-8)
339+
340+
\[2\] [Donvil and Muratore-Ginanneschi. New J. Phys. (2023).](https://dx.doi.org/10.1088/1367-2630/acd4dc)
341+
342+
\[3\] [Donvil and Muratore-Ginanneschi. *Open Systems & Information Dynamics*.](https://www.worldscientific.com/worldscinet/osid)
343+
344+
\[4\] [Breuer and Petruccione *The Theory of Open Quantum Systems*.](https://doi.org/10.1093/acprof:oso/9780199213900.001.0001)
345+
346+
347+
348+
## About
349+
350+
```python
351+
about()
352+
```
353+
354+
## Testing
355+
356+
```python
357+
# TODO
358+
```

0 commit comments

Comments
 (0)