forked from qutip/qutip
-
Notifications
You must be signed in to change notification settings - Fork 1
/
bloch_redfield.py
322 lines (248 loc) · 9.86 KB
/
bloch_redfield.py
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
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
# This file is part of QuTiP: Quantum Toolbox in Python.
#
# Copyright (c) 2011 and later, Paul D. Nation and Robert J. Johansson.
# All rights reserved.
#
# Redistribution and use in source and binary forms, with or without
# modification, are permitted provided that the following conditions are
# met:
#
# 1. Redistributions of source code must retain the above copyright notice,
# this list of conditions and the following disclaimer.
#
# 2. Redistributions in binary form must reproduce the above copyright
# notice, this list of conditions and the following disclaimer in the
# documentation and/or other materials provided with the distribution.
#
# 3. Neither the name of the QuTiP: Quantum Toolbox in Python nor the names
# of its contributors may be used to endorse or promote products derived
# from this software without specific prior written permission.
#
# THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
# "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
# LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A
# PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
# HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
# SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
# LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
# DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
# THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
# (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
# OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
###############################################################################
import numpy as np
import scipy.integrate
from qutip.qobj import Qobj
from qutip.superoperator import *
from qutip.expect import expect
from qutip.states import *
from qutip.odeoptions import Odeoptions
from qutip.cy.spmatfuncs import cy_ode_rhs
from qutip.odedata import Odedata
#------------------------------------------------------------------------------
# Solve the Bloch-Redfield master equation
#
#
def brmesolve(H, psi0, tlist, a_ops, e_ops=[], spectra_cb=[],
args={}, options=Odeoptions()):
"""
Solve the dynamics for the system using the Bloch-Redfeild master equation.
.. note::
This solver does not currently support time-dependent Hamiltonian or
collapse operators.
Parameters
----------
H : :class:`qutip.qobj`
System Hamiltonian.
rho0 / psi0: :class:`qutip.qobj`
Initial density matrix or state vector (ket).
tlist : *list* / *array*
List of times for :math:`t`.
a_ops : list of :class:`qutip.qobj`
List of system operators that couple to bath degrees of freedom.
e_ops : list of :class:`qutip.qobj` / callback function
List of operators for which to evaluate expectation values.
args : *dictionary*
Dictionary of parameters for time-dependent Hamiltonians and collapse
operators.
options : :class:`qutip.Qdeoptions`
Options for the ODE solver.
Returns
-------
output: :class:`qutip.odedata`
An instance of the class :class:`qutip.odedata`, which contains either
a list of expectation values, for operators given in e_ops, or a list
of states for the times specified by `tlist`.
"""
if not spectra_cb:
# default to infinite temperature white noise
spectra_cb = [lambda w: 1.0 for _ in a_ops]
R, ekets = bloch_redfield_tensor(H, a_ops, spectra_cb)
output = Odedata()
output.times = tlist
results = bloch_redfield_solve(R, ekets, psi0, tlist, e_ops, options)
if e_ops:
output.expect = results
else:
output.states = results
return output
#------------------------------------------------------------------------------
# Evolution of the Bloch-Redfield master equation given the Bloch-Redfield
# tensor.
#
def bloch_redfield_solve(R, ekets, rho0, tlist, e_ops=[], options=None):
"""
Evolve the ODEs defined by Bloch-Redfield master equation. The
Bloch-Redfield tensor can be calculated by the function
:func:`bloch_redfield_tensor`.
Parameters
----------
R : :class:`qutip.qobj`
Bloch-Redfield tensor.
ekets : array of :class:`qutip.qobj`
Array of kets that make up a basis tranformation for the eigenbasis.
rho0 : :class:`qutip.qobj`
Initial density matrix.
tlist : *list* / *array*
List of times for :math:`t`.
e_ops : list of :class:`qutip.qobj` / callback function
List of operators for which to evaluate expectation values.
options : :class:`qutip.Qdeoptions`
Options for the ODE solver.
Returns
-------
output: :class:`qutip.odedata`
An instance of the class :class:`qutip.odedata`, which contains either
an *array* of expectation values for the times specified by `tlist`.
"""
if options is None:
options = Odeoptions()
if options.tidy:
R.tidyup()
#
# check initial state
#
if isket(rho0):
# Got a wave function as initial state: convert to density matrix.
rho0 = rho0 * rho0.dag()
#
# prepare output array
#
n_e_ops = len(e_ops)
n_tsteps = len(tlist)
dt = tlist[1] - tlist[0]
result_list = []
for op in e_ops:
if op.isherm and rho0.isherm:
result_list.append(np.zeros(n_tsteps))
else:
result_list.append(np.zeros(n_tsteps, dtype=complex))
#
# transform the initial density matrix and the e_ops opterators to the
# eigenbasis
#
rho = rho0.transform(ekets)
e_eb_ops = [e.transform(ekets) for e in e_ops]
#
# setup integrator
#
initial_vector = mat2vec(rho.full())
r = scipy.integrate.ode(cy_ode_rhs)
r.set_f_params(R.data.data, R.data.indices, R.data.indptr)
r.set_integrator('zvode', method=options.method, order=options.order,
atol=options.atol, rtol=options.rtol, nsteps=options.nsteps,
first_step=options.first_step, min_step=options.min_step,
max_step=options.max_step)
r.set_initial_value(initial_vector, tlist[0])
#
# start evolution
#
t_idx = 0
for t in tlist:
if not r.successful():
break
rho.data = vec2mat(r.y)
# calculate all the expectation values, or output rho if no operators
if e_ops:
rho_tmp = Qobj(rho)
for m, e in enumerate(e_eb_ops):
result_list[m][t_idx] = expect(e, rho_tmp)
else:
result_list.append(rho.transform(ekets, True))
r.integrate(r.t + dt)
t_idx += 1
return result_list
# -----------------------------------------------------------------------------
# Functions for calculting the Bloch-Redfield tensor for a time-independent
# system.
#
def bloch_redfield_tensor(H, a_ops, spectra_cb, use_secular=True):
"""
Calculate the Bloch-Redfield tensor for a system given a set of operators
and corresponding spectral functions that describes the system's coupling
to its environment.
Parameters
----------
H : :class:`qutip.qobj`
System Hamiltonian.
a_ops : list of :class:`qutip.qobj`
List of system operators that couple to the environment.
spectra_cb : list of callback functions
List of callback functions that evaluate the noise power spectrum
at a given frequency.
use_secular : bool
Flag (True of False) that indicates if the secular approximation should
be used.
Returns
-------
R, kets: :class:`qutip.qobj`, list of :class:`qutip.qobj`
R is the Bloch-Redfield tensor and kets is a list eigenstates of the
Hamiltonian.
"""
# Sanity checks for input parameters
if not isinstance(H, Qobj):
raise TypeError("H must be an instance of Qobj")
for a in a_ops:
if not isinstance(a, Qobj) or not a.isherm:
raise TypeError("Operators in a_ops must be Hermitian Qobj.")
# default spectrum
if not spectra_cb:
spectra_cb = [lambda w: 1.0 for _ in a_ops]
# use the eigenbasis
evals, ekets = H.eigenstates()
N = len(evals)
K = len(a_ops)
A = np.zeros((K, N, N), dtype=complex) # TODO: use sparse here
W = np.zeros((N, N))
# pre-calculate matrix elements
for n in range(N):
for m in range(N):
W[m, n] = np.real(evals[m] - evals[n])
for k in range(K):
#A[k,n,m] = a_ops[k].matrix_element(ekets[n], ekets[m])
A[k, :, :] = a_ops[k].transform(ekets).full()
dw_min = abs(W[W.nonzero()]).min()
# unitary part
Heb = H.transform(ekets)
R = -1.0j * (spre(Heb) - spost(Heb))
R.data = R.data.tolil()
for I in range(N * N):
a, b = vec2mat_index(N, I)
for J in range(N * N):
c, d = vec2mat_index(N, J)
# unitary part: use spre and spost above, same as this:
#R.data[I,J] = -1j * W[a,b] * (a == c) * (b == d)
if use_secular is False or abs(W[a, b] - W[c, d]) < dw_min / 10.0:
# dissipative part:
for k in range(K):
# for each operator coupling the system to the environment
R.data[I, J] += ((A[k, a, c] * A[k, d, b] / 2) *
(spectra_cb[k](W[c, a]) +
spectra_cb[k](W[d, b])))
s1 = s2 = 0
for n in range(N):
s1 += A[k, a, n] * A[k, n, c] * spectra_cb[k](W[c, n])
s2 += A[k, d, n] * A[k, n, b] * spectra_cb[k](W[d, n])
R.data[I, J] += - (b == d) * s1 / 2 - (a == c) * s2 / 2
R.data = R.data.tocsr()
return R, ekets