Author: J. R. Johansson (robert@riken.jp), https://jrjohansson.github.io/
This lecture series was developed by J.R. Johannson. The original lecture notebooks are available here.
This is a slightly modified version of the lectures, to work with the current release of QuTiP. You can find these lectures as a part of the qutip-tutorials repository. This lecture and other tutorial notebooks are indexed at the QuTiP Tutorial webpage.
import matplotlib.pyplot as plt
import numpy as np
from IPython.display import Image
from qutip import (Qobj, about, basis, coherent, coherent_dm, create, destroy,
expect, fock, fock_dm, mesolve, qeye, sigmax, sigmay,
sigmaz, tensor, thermal_dm, anim_matrix_histogram,
anim_fock_distribution)
# set a parameter to see animations in line
from matplotlib import rc
rc('animation', html='jshtml')
%matplotlib inline
QuTiP is a python package for calculations and numerical simulations of quantum systems.
It includes facilities for representing and doing calculations with quantum objects such state vectors (wavefunctions), as bras/kets/density matrices, quantum operators of single and composite systems, and superoperators (useful for defining master equations).
It also includes solvers for a time-evolution of quantum systems, according to: Schrodinger equation, von Neuman equation, master equations, Floquet formalism, Monte-Carlo quantum trajectors, experimental implementations of the stochastic Schrodinger/master equations.
For more information see the project web site at qutip.org, and the QuTiP documentation.
You can install QuTiP directly from pip
by running:
pip install qutip
For further installation details, refer to the GitHub repository.
qobj
¶At the heart of the QuTiP package is the Qobj
class, which is used for representing quantum object such as states and operator.
The Qobj
class contains all the information required to describe a quantum system, such as its matrix representation, composite structure and dimensionality.
Image(filename="images/qobj.png")
We can create a new quantum object using the Qobj
class constructor, like this:
q = Qobj([[1], [0]])
q
Here we passed python list as an argument to the class constructor. The data in this list is used to construct the matrix representation of the quantum objects, and the other properties of the quantum object is by default computed from the same data.
We can inspect the properties of a Qobj
instance using the following class method:
# the dimension, or composite Hilbert state space structure
q.dims
[[2], [1]]
# the shape of the matrix data representation
q.shape
(2, 1)
# the matrix data itself. in sparse matrix format.
q.data
Dense(shape=(2, 1), fortran=True)
# get the dense matrix representation
q.full()
array([[1.+0.j], [0.+0.j]])
# some additional properties
q.isherm, q.type
(False, 'ket')
Qobj
instances for calculations¶With Qobj
instances we can do arithmetic and apply a number of different operations using class methods:
sy = Qobj([[0, -1j], [1j, 0]]) # the sigma-y Pauli operator
sy
sz = Qobj([[1, 0], [0, -1]]) # the sigma-z Pauli operator
sz
# some arithmetic with quantum objects
H = 1.0 * sz + 0.1 * sy
print("Qubit Hamiltonian = \n")
H
Qubit Hamiltonian =
Example of modifying quantum objects using the Qobj
methods:
# The hermitian conjugate
sy.dag()
# The trace
H.tr()
0.0
# Eigen energies
H.eigenenergies()
array([-1.00498756, 1.00498756])
For a complete list of methods and properties of the Qobj
class, see the QuTiP documentation or try help(Qobj)
or dir(Qobj)
.
Normally we do not need to create Qobj
instances from stratch, using its constructor and passing its matrix represantation as argument. Instead we can use functions in QuTiP that generates common states and operators for us. Here are some examples of built-in state functions:
# Fundamental basis states (Fock states of oscillator modes)
N = 2 # number of states in the Hilbert space
n = 1 # the state that will be occupied
basis(N, n) # equivalent to fock(N, n)
fock(4, 2) # another example
# a coherent state
coherent(N=10, alpha=1.0)
# a fock state as density matrix
fock_dm(5, 2) # 5 = hilbert space size, 2 = state that is occupied
# coherent state as density matrix
coherent_dm(N=8, alpha=1.0)
# thermal state
n = 1 # average number of thermal photons
thermal_dm(8, n)
# Pauli sigma x
sigmax()
# Pauli sigma y
sigmay()
# Pauli sigma z
sigmaz()
# annihilation operator
destroy(N=8) # N = number of fock states included in the Hilbert space
# creation operator
create(N=8) # equivalent to destroy(5).dag()
# the position operator is easily constructed from the annihilation operator
a = destroy(8)
x = a + a.dag()
x
Qobj
instances we can check some well known commutation relations:¶def commutator(op1, op2):
return op1 * op2 - op2 * op1
$[a, a^1] = 1$
a = destroy(5)
commutator(a, a.dag())
Ops... The result is not identity! Why? Because we have truncated the Hilbert space. But that's OK as long as the highest Fock state isn't involved in the dynamics in our truncated Hilbert space. If it is, the approximation that the truncation introduces might be a problem.
$[x,p] = i$
x = (a + a.dag()) / np.sqrt(2)
p = -1j * (a - a.dag()) / np.sqrt(2)
commutator(x, p)
Same issue with the truncated Hilbert space, but otherwise OK.
Let's try some Pauli spin inequalities
$[\sigma_x, \sigma_y] = 2i \sigma_z$
commutator(sigmax(), sigmay()) - 2j * sigmaz()
$-i \sigma_x \sigma_y \sigma_z = \mathbf{1}$
-1j * sigmax() * sigmay() * sigmaz()
$\sigma_x^2 = \sigma_y^2 = \sigma_z^2 = \mathbf{1}$
sigmax() ** 2 == sigmay() ** 2 == sigmaz() ** 2 == qeye(2)
True
In most cases we are interested in coupled quantum systems, for example coupled qubits, a qubit coupled to a cavity (oscillator mode), etc.
To define states and operators for such systems in QuTiP, we use the tensor
function to create Qobj
instances for the composite system.
For example, consider a system composed of two qubits. If we want to create a Pauli $\sigma_z$ operator that acts on the first qubit and leaves the second qubit unaffected (i.e., the operator $\sigma_z \otimes \mathbf{1}$), we would do:
sz1 = tensor(sigmaz(), qeye(2))
sz1
We can easily verify that this two-qubit operator does indeed have the desired properties:
psi1 = tensor(basis(N, 1), basis(N, 0)) # excited first qubit
psi2 = tensor(basis(N, 0), basis(N, 1)) # excited second qubit
# this should not be true,
# because sz1 should flip the sign of the excited state of psi1
sz1 * psi1 == psi1
False
# this should be true, because sz1 should leave psi2 unaffected
sz1 * psi2 == psi2
True
Above we used the qeye(N)
function, which generates the identity operator with N
quantum states. If we want to do the same thing for the second qubit we can do:
sz2 = tensor(qeye(2), sigmaz())
sz2
Note the order of the argument to the tensor
function, and the correspondingly different matrix representation of the two operators sz1
and sz2
.
Using the same method we can create coupling terms of the form $\sigma_x \otimes \sigma_x$:
tensor(sigmax(), sigmax())
Now we are ready to create a Qobj
representation of a coupled two-qubit Hamiltonian: $H = \epsilon_1 \sigma_z^{(1)} + \epsilon_2 \sigma_z^{(2)} + g \sigma_x^{(1)}\sigma_x^{(2)}$
epsilon = [1.0, 1.0]
g = 0.1
sz1 = tensor(sigmaz(), qeye(2))
sz2 = tensor(qeye(2), sigmaz())
H = epsilon[0] * sz1 + epsilon[1] * sz2 + g * tensor(sigmax(), sigmax())
H
To create composite systems of different types, all we need to do is to change the operators that we pass to the tensor
function (which can take an arbitrary number of operator for composite systems with many components).
For example, the Jaynes-Cumming Hamiltonian for a qubit-cavity system:
$H = \omega_c a^\dagger a - \frac{1}{2}\omega_a \sigma_z + g (a \sigma_+ + a^\dagger \sigma_-)$
wc = 1.0 # cavity frequency
wa = 1.0 # qubit/atom frenqency
g = 0.1 # coupling strength
# cavity mode operator
a = tensor(destroy(5), qeye(2))
# qubit/atom operators
sz = tensor(qeye(5), sigmaz()) # sigma-z operator
sm = tensor(qeye(5), destroy(2)) # sigma-minus operator
# the Jaynes-Cumming Hamiltonian
H = wc * a.dag() * a - 0.5 * wa * sz + g * (a * sm.dag() + a.dag() * sm)
H
Note that
$a \sigma_+ = (a \otimes \mathbf{1}) (\mathbf{1} \otimes \sigma_+)$
so the following two are identical:
a = tensor(destroy(3), qeye(2))
sp = tensor(qeye(3), create(2))
a * sp
tensor(destroy(3), create(2))
Unitary evolution of a quantum system in QuTiP can be calculated with the mesolve
function.
mesolve
is short for Master-eqaution solve (for dissipative dynamics), but if no collapse operators (which describe the dissipation) are given to the solve it falls back on the unitary evolution of the Schrodinger (for initial states in state vector for) or the von Neuman equation (for initial states in density matrix form).
The evolution solvers in QuTiP returns a class of type Odedata
, which contains the solution to the problem posed to the evolution solver.
For example, considor a qubit with Hamiltonian $H = \sigma_x$ and initial state $\left|1\right>$ (in the sigma-z basis): Its evolution can be calculated as follows:
# Hamiltonian
H = sigmax()
# initial state
psi0 = basis(2, 0)
# list of times for which the solver should store the state vector
tlist = np.linspace(0, 10, 100)
result = mesolve(H, psi0, tlist, [], [])
result
<Result Solver: sesolve Solver stats: method: 'scipy zvode adams' init time: 0.00015783309936523438 preparation time: 0.00011277198791503906 run time: 0.005159139633178711 solver: 'Schrodinger Evolution' Time interval: [0.0, 10.0] (100 steps) Number of e_ops: 0 States saved. >
The result
object contains a list of the wavefunctions at the times requested with the tlist
array.
len(result.states)
100
result.states[-1] # the finial state
You can visualize the time evolution of the state. anim_matrix_histogram
visualizes the elements of it. The animation below shows that the state changes periodically and that the coefficient of $\left|1\right>$ is real and that of $\left|0\right>$ is imaginary.
fig, ani = anim_matrix_histogram(result, limits=[0, 1],
bar_style='abs', color_style='phase')
# close an auto-generated plot and animation
plt.close()
ani
The expectation values of an operator given a state vector or density matrix (or list thereof) can be calculated using the expect
function.
expect(sigmaz(), result.states[-1])
0.40806479408005597
expect(sigmaz(), result.states)
array([ 1. , 0.97966324, 0.9194801 , 0.82189845, 0.69088728, 0.53177528, 0.35103398, 0.15601481, -0.04535002, -0.24487027, -0.43443068, -0.60632131, -0.7535507 , -0.87013055, -0.951319 , -0.99381383, -0.9958867 , -0.9574533 , -0.88007676, -0.76690429, -0.62253905, -0.45285285, -0.26474744, -0.06587387, 0.13567909, 0.3317135 , 0.51425587, 0.67588157, 0.81001677, 0.91120556, 0.97533227, 0.99978876, 0.98358026, 0.92736593, 0.83343229, 0.70559992, 0.54906818, 0.37020385, 0.17628203, -0.0248099 , -0.22489269, -0.41582823, -0.58985057, -0.73988155, -0.85981879, -0.9447841 , -0.9913216 , -0.99753841, -0.96318177, -0.88964906, -0.77993099, -0.6384903 , -0.47107995, -0.284509 , -0.08636599, 0.11528974, 0.31225626, 0.49652225, 0.66059275, 0.79779446, 0.90254706, 0.97058988, 0.99915523, 0.98708135, 0.93485936, 0.84461322, 0.7200135 , 0.56612831, 0.38921658, 0.19647394, -0.00425997, -0.20482058, -0.39705049, -0.57313082, -0.72589982, -0.84914389, -0.93785015, -0.98841059, -0.99876884, -0.9685036 , -0.89884578, -0.7926286 , -0.65417236, -0.48910854, -0.30415083, -0.10682227, 0.09485115, 0.29266668, 0.47857836, 0.64502446, 0.78523515, 0.89350741, 0.96543742, 0.99809965, 0.9901656 , 0.9419579 , 0.85543735, 0.73412306, 0.58294927, 0.40806479])
fig, axes = plt.subplots(1, 1)
axes.plot(tlist, expect(sigmaz(), result.states))
axes.set_xlabel(r"$t$", fontsize=20)
axes.set_ylabel(r"$\left<\sigma_z\right>$", fontsize=20);
If we are only interested in expectation values, we could pass a list of operators to the mesolve
function that we want expectation values for, and have the solver compute then and store the results in the Odedata
class instance that it returns.
For example, to request that the solver calculates the expectation values for the operators $\sigma_x, \sigma_y, \sigma_z$:
result = mesolve(H, psi0, tlist, [], [sigmax(), sigmay(), sigmaz()])
Now the expectation values are available in result.expect[0]
, result.expect[1]
, and result.expect[2]
:
fig, axes = plt.subplots(1, 1)
axes.plot(tlist, result.expect[2], label=r"$\left<\sigma_z\right>$")
axes.plot(tlist, result.expect[1], label=r"$\left<\sigma_y\right>$")
axes.plot(tlist, result.expect[0], label=r"$\left<\sigma_x\right>$")
axes.set_xlabel(r"$t$", fontsize=20)
axes.legend(loc=2);
To add dissipation to a problem, all we need to do is to define a list of collapse operators to the call to the mesolve
solver.
A collapse operator is an operator that describes how the system is interacting with its environment.
For example, consider a quantum harmonic oscillator with Hamiltonian
$H = \hbar\omega a^\dagger a$
and which loses photons to its environment with a relaxation rate $\kappa$. The collapse operator that describes this process is
$\sqrt{\kappa} a$
since $a$ is the photon annihilation operator of the oscillator.
To program this problem in QuTiP:
w = 1.0 # oscillator frequency
kappa = 0.1 # relaxation rate
a = destroy(10) # oscillator annihilation operator
rho0 = fock_dm(10, 5) # initial state, fock state with 5 photons
H = w * a.dag() * a # Hamiltonian
# A list of collapse operators
c_ops = [np.sqrt(kappa) * a]
tlist = np.linspace(0, 50, 100)
# request that the solver return the expectation value
# of the photon number state operator a.dag() * a
result = mesolve(H, rho0, tlist, c_ops, [a.dag() * a],
options={"store_states": True})
fig, axes = plt.subplots(1, 1)
axes.plot(tlist, result.expect[0])
axes.set_xlabel(r"$t$", fontsize=20)
axes.set_ylabel(r"Photon number", fontsize=16);
anim_fock_distribution
enables you to visualize the probability distribution at each time.
fig, ani = anim_fock_distribution(result)
# close an auto-generated plot and animation
plt.close()
ani
about()
QuTiP: Quantum Toolbox in Python ================================ Copyright (c) QuTiP team 2011 and later. Current admin team: Alexander Pitchford, Nathan Shammah, Shahnawaz Ahmed, Neill Lambert, Eric Giguère, Boxi Li, Jake Lishman, Simon Cross, Asier Galicia, Paul Menczel, and Patrick Hopf. Board members: Daniel Burgarth, Robert Johansson, Anton F. Kockum, Franco Nori and Will Zeng. Original developers: R. J. Johansson & P. D. Nation. Previous lead developers: Chris Granade & A. Grimsmo. Currently developed through wide collaboration. See https://github.com/qutip for details. QuTiP Version: 5.1.0.dev0+c874c4a Numpy Version: 1.22.4 Scipy Version: 1.13.1 Cython Version: 3.0.10 Matplotlib Version: 3.5.2 Python Version: 3.10.4 Number of CPUs: 4 BLAS Info: Generic INTEL MKL Ext: False Platform Info: Linux (x86_64) Installation path: /usr/share/miniconda3/envs/test-environment/lib/python3.10/site-packages/qutip ================================================================================ Please cite QuTiP in your publication. ================================================================================ For your convenience a bibtex reference can be easily generated using `qutip.cite()`