Jonathan Zoller (jonathan.zoller@uni-ulm.de)
Example to demonstrate using the control library to determine control pulses using the ctrlpulseoptim.optimize_pulse_unitary function. The CRAB algorithm is used to optimize pulse shapes to minimize the fidelity error, which is equivalent maximising the fidelity to an optimal value of 1.
The system in this example are two qubits, where the interaction can be controlled. The target is to perform a pure state transfer from a down-down state to an up-up state.
The user can experiment with the timeslicing, by means of changing the number of timeslots and/or total time for the evolution. Different initial (starting) pulse types can be tried as well as boundaries on the control and a smooth ramping of the pulse when switching the control on and off (at the beginning and close to the end). The initial and final pulses are displayed in a plot
An in depth discussion of using methods of this type can be found in [1,2]
%matplotlib inline
import numpy as np
import matplotlib.pyplot as plt
import datetime
from qutip import Qobj, identity, sigmax, sigmaz, tensor
import random
import qutip.logging_utils as logging
logger = logging.get_logger()
#Set this to None or logging.WARN for 'quiet' execution
log_level = logging.INFO
#QuTiP control modules
import qutip.control.pulseoptim as cpo
example_name = '2qubitInteract'
The dynamics of the system are governed by the combined Hamiltonian: H(t) = H_d + sum(u1(t)Hc1 + u2(t)Hc2 + ....) That is the time-dependent Hamiltonian has a constant part (called here the drift) and time vary parts, which are the control Hamiltonians scaled by some functions u_j(t) known as control amplitudes In this example we describe an Ising like Hamiltonian, encompassing random coefficients in the drift part and controlling the interaction of the qubits:
$ \hat{H} = \sum_{i=1}^2 \alpha_i \sigma_x^i + \beta_i \sigma_z^i + u(t) \cdot \sigma_z \otimes \sigma_z $
Initial $\newcommand{\ket}[1]{\left|{#1}\right\rangle} \ket{\psi_0} = \text{U_0}$ and target state $\ket{\psi_t} = \text{U_targ}$ are chosen to be:
$ \ket{\psi_0} = \begin{pmatrix} 1 \\ 0 \\ 0 \\ 0 \end{pmatrix}$
$ \ket{\psi_t} = \begin{pmatrix} 0 \\ 0 \\ 0 \\ 1 \end{pmatrix}$
random.seed(20)
alpha = [random.random(),random.random()]
beta = [random.random(),random.random()]
Sx = sigmax()
Sz = sigmaz()
H_d = (alpha[0]*tensor(Sx,identity(2)) +
alpha[1]*tensor(identity(2),Sx) +
beta[0]*tensor(Sz,identity(2)) +
beta[1]*tensor(identity(2),Sz))
H_c = [tensor(Sz,Sz)]
# Number of ctrls
n_ctrls = len(H_c)
q1_0 = q2_0 = Qobj([[1], [0]])
q1_targ = q2_targ = Qobj([[0], [1]])
psi_0 = tensor(q1_0, q2_0)
psi_targ = tensor(q1_targ, q2_targ)
To solve the evolution the control amplitudes are considered constant within piecewise timeslots, hence the evolution during the timeslot can be calculated using U(t_k) = expm(-iH(t_k)dt). Combining these for all the timeslots gives the approximation to the evolution from an initial state $\psi_0$ at t=0 to U(T) at the t=evo_time. The number of timeslots and evo_time have to be chosen such that the timeslot durations (dt) are small compared with the dynamics of the system.
# Number of time slots
n_ts = 100
# Time allowed for the evolution
evo_time = 18
At each iteration the fidelity of the evolution is tested by comparaing the calculated evolution U(T) with the target U_targ. For unitary systems such as this one this is typically: f = normalise(overlap(U(T), U_targ)). The maximum fidelity (for a unitary system) calculated this way would be 1, and hence the error is calculated as fid_err = 1 - fidelity. As such the optimisation is considered completed when the fid_err falls below such a target value.
In some cases the optimisation either gets stuck in some local minima, or the fid_err_targ is just not achievable, therefore some limits are set to the time/effort allowed to find a solution.
The algorithm uses the CRAB algorithm to determine optimized coefficients that lead to a minimal fidelity error. The underlying optimization procedure is set to be the Nelder-Mead downhill simplex. Therefore, when all vertices shrink together, the algorithm will terminate.
# Fidelity error target
fid_err_targ = 1e-3
# Maximum iterations for the optisation algorithm
max_iter = 500
# Maximum (elapsed) time allowed in seconds
max_wall_time = 120
The control amplitudes must be set to some initial values. Typically these are just random values for each control in each timeslot. These do however result in erratic optimised pulses. For this example, a solution will be found for any initial pulse, and so it can be interesting to look at the other initial pulse alternatives.
# pulse type alternatives: RND|ZERO|LIN|SINE|SQUARE|SAW|TRIANGLE|
p_type = 'DEF'
#Set to None to suppress output files
f_ext = "{}_n_ts{}_ptype{}.txt".format(example_name, n_ts, p_type)
In this step, the actual optimization is performed. At each iteration the Nelder-Mead algorithm calculates a new set of coefficients that improves the currently worst set among all set of coefficients. For details see [1,2] and a textbook about static search methods. The algorithm continues until one of the termination conditions defined above has been reached. If undesired results are achieved, rerun the algorithm and/or try to change the number of coefficients to be optimized for, as this is a very crucial parameter.
result = cpo.opt_pulse_crab_unitary(H_d, H_c, psi_0, psi_targ, n_ts, evo_time,
fid_err_targ=fid_err_targ,
max_iter=max_iter, max_wall_time=max_wall_time,
init_coeff_scaling=5.0, num_coeffs=5,
method_params={'xtol':1e-3},
guess_pulse_type=None, guess_pulse_action='modulate',
out_file_ext=f_ext,
log_level=log_level, gen_stats=True)
INFO:qutip.control.dynamics:Setting memory optimisations for level 0 INFO:qutip.control.dynamics:Internal operator data type choosen to be <class 'numpy.ndarray'> INFO:qutip.control.dynamics:phased dynamics generator caching True INFO:qutip.control.dynamics:propagator gradient caching True INFO:qutip.control.dynamics:eigenvector adjoint caching True INFO:qutip.control.dynamics:use sparse eigen decomp False INFO:qutip.control.pulseoptim:System configuration: Drift Hamiltonian: Quantum object: dims = [[2, 2], [2, 2]], shape = (4, 4), type = oper, isherm = True Qobj data = [[ 1.67112549 0.68625416 0.90563968 0. ] [ 0.68625416 -0.13810698 0. 0.90563968] [ 0.90563968 0. 0.13810698 0.68625416] [ 0. 0.90563968 0.68625416 -1.67112549]] Control 1 Hamiltonian: Quantum object: dims = [[2, 2], [2, 2]], shape = (4, 4), type = oper, isherm = True Qobj data = [[ 1. 0. 0. 0.] [ 0. -1. 0. 0.] [ 0. 0. -1. 0.] [ 0. 0. 0. 1.]] Initial state / operator: Quantum object: dims = [[2, 2], [1, 1]], shape = (4, 1), type = ket Qobj data = [[ 1.] [ 0.] [ 0.] [ 0.]] Target state / operator: Quantum object: dims = [[2, 2], [1, 1]], shape = (4, 1), type = ket Qobj data = [[ 0.] [ 0.] [ 0.] [ 1.]] INFO:qutip.control.pulseoptim:Initial amplitudes output to file: ctrl_amps_initial_2qubitInteract_n_ts100_ptypeDEF.txt INFO:qutip.control.optimizer:Optimising pulse(s) using CRAB with 'fmin' (Nelder-Mead) method INFO:qutip.control.pulseoptim:Final amplitudes output to file: ctrl_amps_final_2qubitInteract_n_ts100_ptypeDEF.txt
Firstly the performace statistics are reported, which gives a breakdown of the processing times. In this example it can be seen that the majority of time is spent calculating the propagators, i.e. exponentiating the combined Hamiltonian.
The optimised U(T) is reported as the 'final evolution', which is essentially the string representation of the Qobj that holds the full time evolution at the point when the optimisation is terminated.
The key information is in the summary (given last). Here the final fidelity is reported and the reason for termination of the algorithm.
result.stats.report()
print("Final evolution\n{}\n".format(result.evo_full_final))
print("********* Summary *****************")
print("Final fidelity error {}".format(result.fid_err))
print("Final gradient normal {}".format(result.grad_norm_final))
print("Terminated due to {}".format(result.termination_reason))
print("Number of iterations {}".format(result.num_iter))
print("Completed in {} HH:MM:SS.US".format(
datetime.timedelta(seconds=result.wall_time)))
------------------------------------ ---- Control optimisation stats ---- **** Timings (HH:MM:SS.US) **** Total wall time elapsed during optimisation: 0:00:05.187298 Wall time computing Hamiltonians: 0:00:00.096692 (1.86%) Wall time computing propagators: 0:00:04.988672 (96.17%) Wall time computing forward propagation: 0:00:00.027030 (0.52%) Wall time computing onward propagation: 0:00:00.019223 (0.37%) Wall time computing gradient: 0:00:00 (0.00%) **** Iterations and function calls **** Number of iterations: 125 Number of fidelity function calls: 192 Number of times fidelity is computed: 192 Number of gradient function calls: 0 Number of times gradients are computed: 0 Number of times timeslot evolution is recomputed: 192 **** Control amplitudes **** Number of control amplitude updates: 191 Mean number of updates per iteration: 1.528 Number of timeslot values changed: 19099 Mean number of timeslot changes per update: 99.99476439790575 Number of amplitude values changed: 19099 Mean number of amplitude changes per update: 99.99476439790575 ------------------------------------ Final evolution Quantum object: dims = [[2, 2], [1, 1]], shape = (4, 1), type = ket Qobj data = [[ 0.03046450-0.00243049j] [-0.01434436-0.00097183j] [-0.00583341-0.00069853j] [ 0.97780868-0.20667602j]] ********* Summary ***************** Final fidelity error 0.0005877800057074722 Final gradient normal 0.0 Terminated due to Goal achieved Number of iterations 125 Completed in 0:00:05.187298 HH:MM:SS.US
Here the (random) starting pulse is plotted along with the pulse (control amplitudes) that was found to produce the target gate evolution to within the specified error.
fig1 = plt.figure()
ax1 = fig1.add_subplot(2, 1, 1)
ax1.set_title("Initial Control amps")
ax1.set_ylabel("Control amplitude")
ax1.step(result.time,
np.hstack((result.initial_amps[:, 0], result.initial_amps[-1, 0])),
where='post')
ax2 = fig1.add_subplot(2, 1, 2)
ax2.set_title("Optimised Control Amplitudes")
ax2.set_xlabel("Time")
ax2.set_ylabel("Control amplitude")
ax2.step(result.time,
np.hstack((result.final_amps[:, 0], result.final_amps[-1, 0])),
where='post')
plt.tight_layout()
plt.show()
from qutip.ipynbtools import version_table
version_table()
Software | Version |
---|---|
QuTiP | 4.1.0 |
Numpy | 1.11.3 |
SciPy | 0.18.1 |
matplotlib | 2.0.0 |
Cython | 0.25.2 |
Number of CPUs | 4 |
BLAS Info | INTEL MKL |
IPython | 5.1.0 |
Python | 3.6.0 |Anaconda 4.3.1 (64-bit)| (default, Dec 23 2016, 12:22:00) [GCC 4.4.7 20120313 (Red Hat 4.4.7-1)] |
OS | posix [linux] |
Fri Jul 14 11:59:44 2017 BST |
[1] Doria, P., Calarco, T. & Montangero, S.: Optimal Control Technique for Many-Body Quantum Dynamics. Phys. Rev. Lett. 106, 1–4 (2011).
[2] Caneva, T., Calarco, T. & Montangero, S.: Chopped random-basis quantum optimization. Phys. Rev. A - At. Mol. Opt. Phys. 84, (2011).