Calculation of control fields for QFT gate on two qubits using L-BFGS-B algorithm

Alexander Pitchford ([email protected])

Example to demonstrate using the control library to determine control pulses using the ctrlpulseoptim.create_pulse_optimizer function to generate an Optimizer object, through which the configuration can be manipulated before running the optmisation algorithm. In this case it is demonstrated by modifying the initial ctrl pulses. Also re-uses objects in repeated runs with different total evolution times.

The (default) L-BFGS-B algorithm is used to optimise the pulse to minimise the fidelity error, which is equivalent maximising the fidelity to optimal value of 1.

The system in this example is two qubits in constant fields in x, y and z with variable independant controls fields in x and y acting on each qubit The target evolution is the QFT gate. The user can experiment with the different:

  • evolution times - evo_times list values, try anything
  • phase options - phase_option = SU or PSU
  • propagtor computer type prop_type = DIAG or FRECHET
  • fidelity measures - fid_type = UNIT or TRACEDIFF

The user can experiment with the timeslicing, by means of changing the timeslots durations. Different initial (starting) pulse types can be tried. The initial and final pulses are displayed in a plot

This example assumes that the example-control-pulseoptim-Hadamard has already been tried, and hence explanations in that notebook are not repeated here.

In [1]:
%matplotlib inline
import numpy as np
import matplotlib.pyplot as plt
import datetime
In [2]:
from qutip import Qobj, identity, sigmax, sigmay, sigmaz, tensor
from qutip.qip.algorithms import qft
import qutip.logging 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
import qutip.control.pulsegen as pulsegen

example_name = 'QFT'

Defining the physics

Note here that there are two controls acting on each qubit.

In [3]:
Sx = sigmax()
Sy = sigmay()
Sz = sigmaz()
Si = 0.5*identity(2)

# Drift Hamiltonian
H_d = 0.5*(tensor(Sx, Sx) + tensor(Sy, Sy) + tensor(Sz, Sz))
# The (four) control Hamiltonians
H_c = [tensor(Sx, Si), tensor(Sy, Si), tensor(Si, Sx), tensor(Si, Sy)]
n_ctrls = len(H_c)
# start point for the gate evolution
U_0 = identity(4)
# Target for the gate evolution - Quantum Fourier Transform gate
U_targ = qft.qft(2)

Defining the time evolution parameters

Multiple total evolution times will be tried. Using this approach, the minimum evolution time required to achieve the target fidelity could be determined (iteratively).

Note that the timeslot duration dt is fixed, and so the number of timeslots depends on the evo_time

In [4]:
# Duration of each timeslot
dt = 0.02
# List of evolution times to try
evo_times = [1, 3, 6]
n_evo_times = len(evo_times)
evo_time = evo_times[0]
n_ts = int(float(evo_time) / dt)
#Empty list that will hold the results for each evolution time
results = list()

Set the conditions which will cause the pulse optimisation to terminate

In [5]:
# Fidelity error target
fid_err_targ = 1e-10
# Maximum iterations for the optisation algorithm
max_iter = 200
# Maximum (elapsed) time allowed in seconds
max_wall_time = 120
# Minimum gradient (sum of gradients squared)
# as this tends to 0 -> local minima has been found
min_grad = 1e-20

Set the initial pulse type

Here the linear initial pulse type is used, simply because it results in smooth final pulses

In [6]:
# pulse type alternatives: RND|ZERO|LIN|SINE|SQUARE|SAW|TRIANGLE|
p_type = 'LIN'

Give an extension for output files

In [7]:
#Set to None to suppress output files
f_ext = "{}_n_ts{}_ptype{}.txt".format(example_name, n_ts, p_type)

Create the optimiser objects

Here is the main difference between this and the Hadamard example. In this case we use a different pulseoptim function that just creates the objects that can be used to set the physics and configure the optimisation algorithm. This gives greater flexibility (shown here by seting different initial pulse parameters for each control) and is also more efficient when running multiple optimisations on the same system.

In [8]:
optim = cpo.create_pulse_optimizer(H_d, H_c, U_0, U_targ, n_ts, evo_time, 
                amp_lbound=-5.0, amp_ubound=5.0, 
                fid_err_targ=fid_err_targ, min_grad=min_grad, 
                max_iter=max_iter, max_wall_time=max_wall_time, 
                optim_alg='LBFGSB', max_metric_corr=20, accuracy_factor=1e8,
                dyn_type='UNIT', 
                prop_type='DIAG', fid_type='UNIT', phase_option='PSU', 
                init_pulse_type=p_type, 
                log_level=log_level, gen_stats=True)

# **** get handles to the other objects ****
optim.test_out_files = 0
dyn = optim.dynamics
dyn.test_out_files = 0
p_gen = optim.pulse_generator
                

Optimise the pulse for each of the different evolution times

Here a loop is used to perform the optimisation for each of the evo_times given in the list above. The first optimisation is completed using the timeslot parameters passed when the optimisation objects are created. For the subsequent runs, the Dynamics object 'dyn' is used to set the timeslot parameters before the initial pulses are generated and optimisation is completed. Note that using this method, the dyn.initialize_controls method must be called with an array of the initial amplitudes before the optim.run_optimization method is called.

In [9]:
for i in range(n_evo_times):
    # Generate the tau (duration) and time (cumulative) arrays
    # so that it can be used to create the pulse generator
    # with matching timeslots
    dyn.init_time_slots()
    if i > 0:
        # Create a new pulse generator for the new dynamics
        p_gen = pulsegen.create_pulse_gen(p_type, dyn)
        
    #Generate different initial pulses for each of the controls
    init_amps = np.zeros([n_ts, n_ctrls])
    if (p_gen.periodic):
        phase_diff = np.pi / n_ctrls
        for j in range(n_ctrls):
            init_amps[:, j] = p_gen.gen_pulse(start_phase=phase_diff*j)
    elif (isinstance(p_gen, pulsegen.PulseGenLinear)):
        for j in range(n_ctrls):
            p_gen.scaling = float(j) - float(n_ctrls - 1)/2
            init_amps[:, j] = p_gen.gen_pulse()
    elif (isinstance(p_gen, pulsegen.PulseGenZero)):
        for j in range(n_ctrls):
            p_gen.offset = sf = float(j) - float(n_ctrls - 1)/2
            init_amps[:, j] = p_gen.gen_pulse()
    else:
        # Should be random pulse
        for j in range(n_ctrls):
            init_amps[:, j] = p_gen.gen_pulse()
    
    dyn.initialize_controls(init_amps)
    
    # Save initial amplitudes to a text file
    if f_ext is not None:
        pulsefile = "ctrl_amps_initial_" + f_ext
        dyn.save_amps(pulsefile)
        print("Initial amplitudes output to file: " + pulsefile)

    print("***********************************")
    print("Starting pulse optimisation for T={}".format(evo_time))
    result = optim.run_optimization()
    results.append(result)

    # Save final amplitudes to a text file
    if f_ext is not None:
        pulsefile = "ctrl_amps_final_" + f_ext
        dyn.save_amps(pulsefile)
        print("Final amplitudes output to file: " + pulsefile)
    
    # Report the results
    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)))
    
    if i+1 < len(evo_times):
        # reconfigure the dynamics for the next evo time
        evo_time = evo_times[i+1]
        n_ts = int(float(evo_time) / dt)
        dyn.tau = None
        dyn.evo_time = evo_time
        dyn.num_tslots = n_ts
INFO:qutip.control.optimizer:Optimising pulse using L-BFGS-B
Initial amplitudes output to file: ctrl_amps_initial_QFT_n_ts50_ptypeLIN.txt
***********************************
Starting pulse optimisation for T=1
Final amplitudes output to file: ctrl_amps_final_QFT_n_ts50_ptypeLIN.txt
INFO:qutip.control.optimizer:Optimising pulse using L-BFGS-B

------------------------------------
---- Control optimisation stats ----
**** Timings (HH:MM:SS.US) ****
Total wall time elapsed during optimisation: 0:00:01.842925
Wall time computing Hamiltonians: 0:00:00.153439 (8.33%)
Wall time computing propagators: 0:00:01.320443 (71.65%)
Wall time computing forward propagation: 0:00:00.014841 (0.81%)
Wall time computing onward propagation: 0:00:00.014233 (0.77%)
Wall time computing gradient: 0:00:00.305658 (16.59%)

**** Iterations and function calls ****
Number of iterations: 43
Number of fidelity function calls: 95
Number of times fidelity is computed: 95
Number of gradient function calls: 95
Number of times gradients are computed: 95
Number of times timeslot evolution is recomputed: 95

**** Control amplitudes ****
Number of control amplitude updates: 94
Mean number of updates per iteration: 2.186046511627907
Number of timeslot values changed: 4700
Mean number of timeslot changes per update: 50.0
Number of amplitude values changed: 14718
Mean number of amplitude changes per update: 156.5744680851064
------------------------------------
Final evolution
Quantum object: dims = [[4], [4]], shape = [4, 4], type = oper, isherm = False
Qobj data =
[[-0.01716076-0.49172084j  0.20798976-0.25787853j  0.64235091-0.29300043j
   0.09595541-0.37481151j]
 [ 0.06651301-0.30493605j -0.02967478+0.05227976j -0.08265828+0.4614182j
  -0.74664327-0.34894259j]
 [ 0.56993572-0.41788874j  0.09618771+0.36054022j -0.21182798-0.46036628j
  -0.15161276+0.2854931j ]
 [ 0.18348661-0.35690487j -0.77936553-0.37412992j -0.01229437+0.15759938j
   0.17861792+0.1862114j ]]

********* Summary *****************
Final fidelity error 0.19185009310042855
Final gradient normal 1.2757724940135868
Terminated due to function converged
Number of iterations 43
Completed in 0:00:01.842925 HH:MM:SS.US
Initial amplitudes output to file: ctrl_amps_initial_QFT_n_ts50_ptypeLIN.txt
***********************************
Starting pulse optimisation for T=3
Final amplitudes output to file: ctrl_amps_final_QFT_n_ts50_ptypeLIN.txt
INFO:qutip.control.optimizer:Optimising pulse using L-BFGS-B

------------------------------------
---- Control optimisation stats ----
**** Timings (HH:MM:SS.US) ****
Total wall time elapsed during optimisation: 0:00:03.000021
Wall time computing Hamiltonians: 0:00:00.251382 (8.38%)
Wall time computing propagators: 0:00:02.177674 (72.59%)
Wall time computing forward propagation: 0:00:00.023806 (0.79%)
Wall time computing onward propagation: 0:00:00.023186 (0.77%)
Wall time computing gradient: 0:00:00.497174 (16.57%)

**** Iterations and function calls ****
Number of iterations: 39
Number of fidelity function calls: 52
Number of times fidelity is computed: 52
Number of gradient function calls: 52
Number of times gradients are computed: 52
Number of times timeslot evolution is recomputed: 52

**** Control amplitudes ****
Number of control amplitude updates: 51
Mean number of updates per iteration: 1.3076923076923077
Number of timeslot values changed: 7650
Mean number of timeslot changes per update: 150.0
Number of amplitude values changed: 30599
Mean number of amplitude changes per update: 599.9803921568628
------------------------------------
Final evolution
Quantum object: dims = [[4], [4]], shape = [4, 4], type = oper, isherm = False
Qobj data =
[[-0.19135424+0.46194414j -0.19134807+0.46192768j -0.19136193+0.46190902j
  -0.19134390+0.46196113j]
 [-0.19136028+0.46192437j -0.46193505-0.1913618j   0.19138278-0.46194621j
   0.46192211+0.19133765j]
 [-0.19134775+0.46195334j  0.19131498-0.46192232j -0.19133422+0.46193961j
   0.19132468-0.46196253j]
 [-0.19132518+0.46192869j  0.46197881+0.19133046j  0.19133309-0.46194552j
  -0.46193573-0.19130645j]]

********* Summary *****************
Final fidelity error 1.3450616176413632e-09
Final gradient normal 0.0009066544926169316
Terminated due to function converged
Number of iterations 39
Completed in 0:00:03.000021 HH:MM:SS.US
Initial amplitudes output to file: ctrl_amps_initial_QFT_n_ts50_ptypeLIN.txt
***********************************
Starting pulse optimisation for T=6
Final amplitudes output to file: ctrl_amps_final_QFT_n_ts50_ptypeLIN.txt

------------------------------------
---- Control optimisation stats ----
**** Timings (HH:MM:SS.US) ****
Total wall time elapsed during optimisation: 0:00:04.720934
Wall time computing Hamiltonians: 0:00:00.397070 (8.41%)
Wall time computing propagators: 0:00:03.428679 (72.63%)
Wall time computing forward propagation: 0:00:00.037914 (0.80%)
Wall time computing onward propagation: 0:00:00.037105 (0.79%)
Wall time computing gradient: 0:00:00.791080 (16.76%)

**** Iterations and function calls ****
Number of iterations: 30
Number of fidelity function calls: 40
Number of times fidelity is computed: 40
Number of gradient function calls: 40
Number of times gradients are computed: 40
Number of times timeslot evolution is recomputed: 40

**** Control amplitudes ****
Number of control amplitude updates: 39
Mean number of updates per iteration: 1.3
Number of timeslot values changed: 11700
Mean number of timeslot changes per update: 300.0
Number of amplitude values changed: 46800
Mean number of amplitude changes per update: 1200.0
------------------------------------
Final evolution
Quantum object: dims = [[4], [4]], shape = [4, 4], type = oper, isherm = False
Qobj data =
[[ 0.46196616+0.19133261j  0.46189141+0.19135946j  0.46196059+0.19133184j
   0.46195347+0.19131259j]
 [ 0.46189646+0.19134815j -0.19134122+0.46195251j -0.46197017-0.19135997j
   0.19133342-0.46193334j]
 [ 0.46194096+0.19133875j -0.46195837-0.19135204j  0.46190385+0.19133346j
  -0.46193897-0.19138344j]
 [ 0.46195960+0.1913374j   0.19134245-0.46194505j -0.46192036-0.19135148j
  -0.19130915+0.46194499j]]

********* Summary *****************
Final fidelity error 1.7512941097308499e-09
Final gradient normal 0.0020285252826135014
Terminated due to function converged
Number of iterations 30
Completed in 0:00:04.720934 HH:MM:SS.US

Plot the initial and final amplitudes

In [10]:
fig1 = plt.figure(figsize=(12,8))
for i in range(n_evo_times):
    #Initial amps
    ax1 = fig1.add_subplot(2, n_evo_times, i+1)
    ax1.set_title("Init amps T={}".format(evo_times[i]))
    # ax1.set_xlabel("Time")
    ax1.get_xaxis().set_visible(False)
    if i == 0:
        ax1.set_ylabel("Control amplitude")
    t = results[i].time[:-1]
    for j in range(n_ctrls):
        amps = results[i].initial_amps[:, j]
        ax1.plot(t, amps)
        
    ax2 = fig1.add_subplot(2, n_evo_times, i+n_evo_times+1)
    ax2.set_title("Final amps T={}".format(evo_times[i]))
    ax2.set_xlabel("Time")
    #Optimised amps
    if i == 0:
        ax2.set_ylabel("Control amplitude")
    for j in range(n_ctrls):
        amps = results[i].final_amps[:, j]
        ax2.plot(t, amps)

plt.show()

Versions

In [11]:
from qutip.ipynbtools import version_table

version_table()
Out[11]:
SoftwareVersion
QuTiP3.1.0
SciPy0.14.1
OSposix [linux]
Numpy1.9.1
Python3.4.0 (default, Apr 11 2014, 13:05:11) [GCC 4.8.2]
matplotlib1.4.2
Cython0.21.2
IPython2.3.1
Tue Jan 13 13:33:51 2015 JST