Render on nbviewer Launch Binder

Tutorial

[1]:
# NBVAL_IGNORE_OUTPUT
%load_ext watermark
import numpy as np
import qutip
import matplotlib
import matplotlib.pylab as plt
import weylchamber
from weylchamber.visualize import WeylChamber
%watermark -v --iversions
Python implementation: CPython
Python version       : 3.9.16
IPython version      : 8.18.1

numpy      : 1.26.4
qutip      : 5.0.2
matplotlib : 3.9.0
weylchamber: 0.5.0+dev

\(\newcommand{Re}[0]{\operatorname{Re}} \newcommand{Im}[0]{\operatorname{Im}} \newcommand{dd}[0]{\,\text{d}} \newcommand{abs}[0]{\operatorname{abs}}\)

Every two-qubit gate is associated with a point in the “Weyl-chamber” that may be visualized in three dimensions as the following polyhedron:

[2]:
WeylChamber().plot()
_images/tutorial_4_0.png

Note: if you run this interactively, and switch to an interactive matplotlib backend, e.g.

%matplotlib tk

you will be able to rotate the 3D plot to get a better intuition.

Consider the following common two-qubit gates:

[3]:
IDENTITY = qutip.identity([2,2])
IDENTITY
[3]:
Quantum object: dims=[[2, 2], [2, 2]], shape=(4, 4), type='oper', dtype=Dia, isherm=True$$\left(\begin{array}{cc}1 & 0 & 0 & 0\\0 & 1 & 0 & 0\\0 & 0 & 1 & 0\\0 & 0 & 0 & 1\end{array}\right)$$
[4]:
CNOT = qutip.core.gates.cnot()
CNOT
[4]:
Quantum object: dims=[[2, 2], [2, 2]], shape=(4, 4), type='oper', dtype=CSR, isherm=True$$\left(\begin{array}{cc}1 & 0 & 0 & 0\\0 & 1 & 0 & 0\\0 & 0 & 0 & 1\\0 & 0 & 1 & 0\end{array}\right)$$
[5]:
CPHASE = qutip.core.gates.cphase(np.pi)
CPHASE
[5]:
Quantum object: dims=[[2, 2], [2, 2]], shape=(4, 4), type='oper', dtype=CSR, isherm=True$$\left(\begin{array}{cc}1 & 0 & 0 & 0\\0 & 1 & 0 & 0\\0 & 0 & 1 & 0\\0 & 0 & 0 & -1\end{array}\right)$$
[6]:
BGATE = qutip.core.gates.berkeley()
BGATE
[6]:
Quantum object: dims=[[2, 2], [2, 2]], shape=(4, 4), type='oper', dtype=Dense, isherm=False$$\left(\begin{array}{cc}0.924 & 0 & 0 & 0.383j\\0 & 0.383 & 0.924j & 0\\0 & 0.924j & 0.383 & 0\\0.383j & 0 & 0 & 0.924\end{array}\right)$$
[7]:
iSWAP = qutip.core.gates.iswap()
iSWAP
[7]:
Quantum object: dims=[[2, 2], [2, 2]], shape=(4, 4), type='oper', dtype=CSR, isherm=False$$\left(\begin{array}{cc}1 & 0 & 0 & 0\\0 & 0 & 1j & 0\\0 & 1j & 0 & 0\\0 & 0 & 0 & 1\end{array}\right)$$
[8]:
sqrtISWAP = qutip.core.gates.sqrtiswap()
sqrtISWAP
[8]:
Quantum object: dims=[[2, 2], [2, 2]], shape=(4, 4), type='oper', dtype=CSR, isherm=False$$\left(\begin{array}{cc}1 & 0 & 0 & 0\\0 & 0.707 & 0.707j & 0\\0 & 0.707j & 0.707 & 0\\0 & 0 & 0 & 1\end{array}\right)$$
[9]:
sqrtSWAP = qutip.core.gates.sqrtswap()
sqrtSWAP
[9]:
Quantum object: dims=[[2, 2], [2, 2]], shape=(4, 4), type='oper', dtype=CSR, isherm=False$$\left(\begin{array}{cc}1 & 0 & 0 & 0\\0 & (0.500+0.500j) & (0.500-0.500j) & 0\\0 & (0.500-0.500j) & (0.500+0.500j) & 0\\0 & 0 & 0 & 1\end{array}\right)$$
[10]:
MGATE = weylchamber.canonical_gate(3/4, 1/4, 0)
MGATE
[10]:
Quantum object: dims=[[2, 2], [2, 2]], shape=(4, 4), type='oper', dtype=Dense, isherm=False$$\left(\begin{array}{cc}0.707 & 0 & 0 & 0.707j\\0 & 0 & 1j & 0\\0 & 1j & 0 & 0\\0.707j & 0 & 0 & 0.707\end{array}\right)$$

All of these gates are situatated at special points in the Weyl chamber. We can print their Weyl chamber coordinates and add a point in the graphical representation

[11]:
w = WeylChamber();
list_of_gates = [
    ('Identity', IDENTITY),
    ('CNOT', CNOT), ('CPHASE', CPHASE), ('BGATE', BGATE),
    ('iSWAP', iSWAP), ('sqrtISWAP', sqrtISWAP),
    ('sqrtSWAP', sqrtSWAP), ('MGATE', MGATE)]
print("Weyl Chamber Coordinates")
print("----------------------------------")
for (name, gate) in list_of_gates:
    c1, c2, c3 = weylchamber.c1c2c3(gate)
    print("%10s: \t%.2fπ %.2fπ %.2fπ" % (name, c1, c2, c3))
    w.add_point(c1, c2, c3)
w.plot()
Weyl Chamber Coordinates
----------------------------------
  Identity:     0.00π 0.00π 0.00π
      CNOT:     0.50π 0.00π 0.00π
    CPHASE:     0.50π 0.00π 0.00π
     BGATE:     0.50π 0.25π 0.00π
     iSWAP:     0.50π 0.50π 0.00π
 sqrtISWAP:     0.25π 0.25π 0.00π
  sqrtSWAP:     0.75π 0.25π 0.25π
     MGATE:     0.75π 0.25π 0.00π
/Users/goerz/Documents/Programming/weylchamber/.venv/py39/lib/python3.9/site-packages/numpy/linalg/linalg.py:2180: RuntimeWarning: divide by zero encountered in det
  r = _umath_linalg.det(a, signature=signature)
/Users/goerz/Documents/Programming/weylchamber/.venv/py39/lib/python3.9/site-packages/numpy/linalg/linalg.py:2180: RuntimeWarning: invalid value encountered in det
  r = _umath_linalg.det(a, signature=signature)
_images/tutorial_16_2.png

The gates locally equivalent to the controlled-phase gates are on an the axis 0 - A1 in the Weyl chamber:

[12]:
w.scatter(*zip(*[
    weylchamber.c1c2c3(qutip.core.gates.cphase(phase))
    for phase in np.linspace(0, 2*np.pi, 20)]))
/Users/goerz/Documents/Programming/weylchamber/.venv/py39/lib/python3.9/site-packages/numpy/linalg/linalg.py:2180: RuntimeWarning: divide by zero encountered in det
  r = _umath_linalg.det(a, signature=signature)
/Users/goerz/Documents/Programming/weylchamber/.venv/py39/lib/python3.9/site-packages/numpy/linalg/linalg.py:2180: RuntimeWarning: invalid value encountered in det
  r = _umath_linalg.det(a, signature=signature)
[13]:
w.plot()
_images/tutorial_19_0.png

The Weyl chamber coordinates \((c_1, c_2, c_3)\) are closely associated with the local invariants \((g_1, g_2, g_3)\)

[14]:
print("Local Invariants")
print("----------------------------------")
for (name, gate) in list_of_gates:
    g1, g2, g3 = weylchamber.g1g2g3(gate)
    print("%10s: \t%5.2f %5.2f %5.2f" % (name, g1, g2, g3))
Local Invariants
----------------------------------
  Identity:      1.00  0.00  3.00
      CNOT:      0.00  0.00  1.00
    CPHASE:      0.00  0.00  1.00
     BGATE:      0.00  0.00  0.00
     iSWAP:      0.00  0.00 -1.00
 sqrtISWAP:      0.25  0.00  1.00
  sqrtSWAP:      0.00 -0.25  0.00
     MGATE:      0.25 -0.00  1.00
/Users/goerz/Documents/Programming/weylchamber/.venv/py39/lib/python3.9/site-packages/numpy/linalg/linalg.py:2180: RuntimeWarning: divide by zero encountered in det
  r = _umath_linalg.det(a, signature=signature)
/Users/goerz/Documents/Programming/weylchamber/.venv/py39/lib/python3.9/site-packages/numpy/linalg/linalg.py:2180: RuntimeWarning: invalid value encountered in det
  r = _umath_linalg.det(a, signature=signature)

This shows that the MGATE and \(\sqrt{\text{iSWAP}}\) are actually locally equivalent, despite being different Weyl chamber coordinates (M and Q, respectively)