Solving Problems with Time-dependent Hamiltonians¶
Methods for Writing Time-Dependent Operators¶
In the previous examples of quantum evolution,
we assumed that the systems under consideration were described by time-independent Hamiltonians.
However, many systems have explicit time dependence in either the Hamiltonian,
or the collapse operators describing coupling to the environment, and sometimes both components might depend on time.
The time-evolutions solvers
qutip.mesolve
, qutip.mcsolve
, qutip.sesolve
, qutip.brmesolve
qutip.ssesolve
, qutip.photocurrent_sesolve
, qutip.smesolve
, and qutip.photocurrent_mesolve
are all capable of handling time-dependent Hamiltonians and collapse terms.
There are, in general, three different ways to implement time-dependent problems in QuTiP:
Function based: Hamiltonian / collapse operators expressed using [qobj, func] pairs, where the time-dependent coefficients of the Hamiltonian (or collapse operators) are expressed using Python functions.
String (Cython) based: The Hamiltonian and/or collapse operators are expressed as a list of [qobj, string] pairs, where the time-dependent coefficients are represented as strings. The resulting Hamiltonian is then compiled into C code using Cython and executed.
Array Based: The Hamiltonian and/or collapse operators are expressed as a list of [qobj, np.array] pairs. The arrays are 1 dimensional and dtype are complex or float. They must contain one value for each time in the tlist given to the solver. Cubic spline interpolation will be used between the given times.
Hamiltonian function (outdated): The Hamiltonian is itself a Python function with time-dependence. Collapse operators must be time independent using this input format.
Give the multiple choices of input style, the first question that arrises is which option to choose?
In short, the function based method (option #1) is the most general,
allowing for essentially arbitrary coefficients expressed via user defined functions.
However, by automatically compiling your system into C++ code,
the second option (string based) tends to be more efficient and will run faster
[This is also the only format that is supported in the qutip.brmesolve
solver].
Of course, for small system sizes and evolution times, the difference will be minor.
Although this method does not support all time-dependent coefficients that one can think of,
it does support essentially all problems that one would typically encounter.
Time-dependent coefficients using any of the following functions,
or combinations thereof (including constants) can be compiled directly into C++-code:
'abs', 'acos', 'acosh', 'arg', 'asin', 'asinh', 'atan', 'atanh', 'conj',
'cos', 'cosh','exp', 'erf', 'zerf', 'imag', 'log', 'log10', 'norm', 'pi',
'proj', 'real', 'sin', 'sinh', 'sqrt', 'tan', 'tanh'
In addition, QuTiP supports cubic spline based interpolation functions [Modeling Non-Analytic and/or Experimental Time-Dependent Parameters using Interpolating Functions].
If you require mathematical functions other than those listed above,
it is possible to call any of the functions in the NumPy library using the prefix np.
before the function name in the string, i.e 'np.sin(t)'
and scipy.special
imported as spe
.
This includes a wide range of functionality, but comes with a small overhead created by going from C++->Python->C++.
Finally option #4, expressing the Hamiltonian as a Python function, is the original method for time dependence in QuTiP 1.x. However, this method is somewhat less efficient then the previously mentioned methods. However, in contrast to the other options this method can be used in implementing time-dependent Hamiltonians that cannot be expressed as a function of constant operators with time-dependent coefficients.
A collection of examples demonstrating the simulation of time-dependent problems can be found on the tutorials web page.
Function Based Time Dependence¶
A very general way to write a time-dependent Hamiltonian or collapse operator is by using Python functions as the time-dependent coefficients. To accomplish this, we need to write a Python function that returns the time-dependent coefficient. Additionally, we need to tell QuTiP that a given Hamiltonian or collapse operator should be associated with a given Python function. To do this, one needs to specify operator-function pairs in list format: [Op, py_coeff]
, where Op
is a given Hamiltonian or collapse operator and py_coeff
is the name of the Python function representing the coefficient. With this format, the form of the Hamiltonian for both mesolve
and mcsolve
is:
>>> H = [H0, [H1, py_coeff1], [H2, py_coeff2], ...]
where H0
is a time-independent Hamiltonian, while H1
,``H2``, are time dependent. The same format can be used for collapse operators:
>>> c_ops = [[C0, py_coeff0], C1, [C2, py_coeff2], ...]
Here we have demonstrated that the ordering of time-dependent and time-independent terms does not matter. In addition, any or all of the collapse operators may be time dependent.
Note
While, in general, you can arrange time-dependent and time-independent terms in any order you like, it is best to place all time-independent terms first.
As an example, we will look at an example that has a time-dependent Hamiltonian of the form \(H=H_{0}-f(t)H_{1}\) where \(f(t)\) is the time-dependent driving strength given as \(f(t)=A\exp\left[-\left( t/\sigma \right)^{2}\right]\). The follow code sets up the problem
In [1]: ustate = basis(3, 0)
In [2]: excited = basis(3, 1)
In [3]: ground = basis(3, 2)
In [4]: N = 2 # Set where to truncate Fock state for cavity
In [5]: sigma_ge = tensor(qeye(N), ground * excited.dag()) # |g><e|
In [6]: sigma_ue = tensor(qeye(N), ustate * excited.dag()) # |u><e|
In [7]: a = tensor(destroy(N), qeye(3))
In [8]: ada = tensor(num(N), qeye(3))
In [9]: c_ops = [] # Build collapse operators
In [10]: kappa = 1.5 # Cavity decay rate
In [11]: c_ops.append(np.sqrt(kappa) * a)
In [12]: gamma = 6 # Atomic decay rate
In [13]: c_ops.append(np.sqrt(5*gamma/9) * sigma_ue) # Use Rb branching ratio of 5/9 e->u
In [14]: c_ops.append(np.sqrt(4*gamma/9) * sigma_ge) # 4/9 e->g
In [15]: t = np.linspace(-15, 15, 100) # Define time vector
In [16]: psi0 = tensor(basis(N, 0), ustate) # Define initial state
In [17]: state_GG = tensor(basis(N, 1), ground) # Define states onto which to project
In [18]: sigma_GG = state_GG * state_GG.dag()
In [19]: state_UU = tensor(basis(N, 0), ustate)
In [20]: sigma_UU = state_UU * state_UU.dag()
In [21]: g = 5 # coupling strength
In [22]: H0 = -g * (sigma_ge.dag() * a + a.dag() * sigma_ge) # time-independent term
In [23]: H1 = (sigma_ue.dag() + sigma_ue) # time-dependent term
Given that we have a single time-dependent Hamiltonian term, and constant collapse terms, we need to specify a single Python function for the coefficient \(f(t)\). In this case, one can simply do
In [24]: def H1_coeff(t, args):
....: return 9 * np.exp(-(t / 5.) ** 2)
....:
In this case, the return value dependents only on time. However, when specifying Python functions for coefficients, the function must have (t,args) as the input variables, in that order. Having specified our coefficient function, we can now specify the Hamiltonian in list format and call the solver (in this case qutip.mesolve
)
In [25]: H = [H0,[H1,H1_coeff]]
In [26]: output = mesolve(H, psi0, t, c_ops, [ada, sigma_UU, sigma_GG])
We can call the Monte Carlo solver in the exact same way (if using the default ntraj=500
):
In [27]: output = mcsolve(H, psi0, t, c_ops, [ada, sigma_UU, sigma_GG])
10.0%. Run time: 0.25s. Est. time left: 00:00:00:02
20.0%. Run time: 0.48s. Est. time left: 00:00:00:01
30.0%. Run time: 0.75s. Est. time left: 00:00:00:01
40.0%. Run time: 0.97s. Est. time left: 00:00:00:01
50.0%. Run time: 1.21s. Est. time left: 00:00:00:01
60.0%. Run time: 1.46s. Est. time left: 00:00:00:00
70.0%. Run time: 1.70s. Est. time left: 00:00:00:00
80.0%. Run time: 1.95s. Est. time left: 00:00:00:00
90.0%. Run time: 2.21s. Est. time left: 00:00:00:00
100.0%. Run time: 2.46s. Est. time left: 00:00:00:00
Total run time: 2.47s
The output from the master equation solver is identical to that shown in the examples, the Monte Carlo however will be noticeably off, suggesting we should increase the number of trajectories for this example. In addition, we can also consider the decay of a simple Harmonic oscillator with time-varying decay rate
In [28]: kappa = 0.5
In [29]: def col_coeff(t, args): # coefficient function
....: return np.sqrt(kappa * np.exp(-t))
....:
In [30]: N = 10 # number of basis states
In [31]: a = destroy(N)
In [32]: H = a.dag() * a # simple HO
In [33]: psi0 = basis(N, 9) # initial state
In [34]: c_ops = [[a, col_coeff]] # time-dependent collapse term
In [35]: times = np.linspace(0, 10, 100)
In [36]: output = mesolve(H, psi0, times, c_ops, [a.dag() * a])
Using the args variable¶
In the previous example we hardcoded all of the variables, driving amplitude \(A\) and width \(\sigma\), with their numerical values. This is fine for problems that are specialized, or that we only want to run once. However, in many cases, we would like to change the parameters of the problem in only one location (usually at the top of the script), and not have to worry about manually changing the values on each run. QuTiP allows you to accomplish this using the keyword args
as an input to the solvers. For instance, instead of explicitly writing 9 for the amplitude and 5 for the width of the gaussian driving term, we can make us of the args variable
In [37]: def H1_coeff(t, args):
....: return args['A'] * np.exp(-(t/args['sigma'])**2)
....:
or equivalently,
In [38]: def H1_coeff(t, args):
....: A = args['A']
....: sig = args['sigma']
....: return A * np.exp(-(t / sig) ** 2)
....:
where args is a Python dictionary of key: value
pairs args = {'A': a, 'sigma': b}
where a
and b
are the two parameters for the amplitude and width, respectively. Of course, we can always hardcode the values in the dictionary as well args = {'A': 9, 'sigma': 5}
, but there is much more flexibility by using variables in args
. To let the solvers know that we have a set of args to pass we append the args
to the end of the solver input:
In [39]: output = mesolve(H, psi0, times, c_ops, [a.dag() * a], args={'A': 9, 'sigma': 5})
or to keep things looking pretty
In [40]: args = {'A': 9, 'sigma': 5}
In [41]: output = mesolve(H, psi0, times, c_ops, [a.dag() * a], args=args)
Once again, the Monte Carlo solver qutip.mcsolve
works in an identical manner.
String Format Method¶
Note
You must have Cython installed on your computer to use this format. See Installation for instructions on installing Cython.
The string-based time-dependent format works in a similar manner as the previously discussed Python function method. That being said, the underlying code does something completely different. When using this format, the strings used to represent the time-dependent coefficients, as well as Hamiltonian and collapse operators, are rewritten as Cython code using a code generator class and then compiled into C code. The details of this meta-programming will be published in due course. however, in short, this can lead to a substantial reduction in time for complex time-dependent problems, or when simulating over long intervals.
Like the previous method, the string-based format uses a list pair format [Op, str]
where str
is now a string representing the time-dependent coefficient. For our first example, this string would be '9 * exp(-(t / 5.) ** 2)'
. The Hamiltonian in this format would take the form:
In [42]: H = [H0, [H1, '9 * exp(-(t / 5) ** 2)']]
Notice that this is a valid Hamiltonian for the string-based format as exp
is included in the above list of suitable functions. Calling the solvers is the same as before:
In [43]: output = mesolve(H, psi0, t, c_ops, [a.dag() * a])
We can also use the args
variable in the same manner as before, however we must rewrite our string term to read: 'A * exp(-(t / sig) ** 2)'
In [44]: H = [H0, [H1, 'A * exp(-(t / sig) ** 2)']]
In [45]: args = {'A': 9, 'sig': 5}
In [46]: output = mesolve(H, psi0, times, c_ops, [a.dag()*a], args=args)
Important
Naming your args
variables exp
, sin
, pi
etc. will cause errors when using the string-based format.
Collapse operators are handled in the exact same way.
Modeling Non-Analytic and/or Experimental Time-Dependent Parameters using Interpolating Functions¶
Note
New in QuTiP 4.1
Sometimes it is necessary to model a system where the time-dependent parameters are non-analytic functions, or are derived from experimental data (i.e. a collection of data points). In these situations, one can use interpolating functions as an approximate functional form for input into a time-dependent solver. QuTiP includes it own custom cubic spline interpolation class qutip.interpolate.Cubic_Spline
to provide this functionality. To see how this works, lets first generate some noisy data:
In [47]: t = np.linspace(-15, 15, 100)
In [48]: func = lambda t: 9*np.exp(-(t / 5)** 2)
In [49]: noisy_func = lambda t: func(t)+(0.05*func(t))*np.random.randn(t.shape[0])
In [50]: noisy_data = noisy_func(t)
In [51]: plt.figure()
Out[51]: <Figure size 640x480 with 0 Axes>
In [52]: plt.plot(t, func(t))