Randomized Benchmarking#

Randomized benchmarking (RB) is a popular protocol for characterizing the error rate of quantum processors. An RB experiment consists of the generation of random Clifford circuits on the given qubits such that the unitary computed by the circuits is the identity. After running the circuits, the number of shots resulting in an error (i.e. an output different from the ground state) are counted, and from this data one can infer error estimates for the quantum device, by calculating the Error Per Clifford. See the Qiskit Textbook for an explanation on the RB method, which is based on Refs. [1] [2].

Note

This tutorial requires the qiskit-aer and qiskit-ibm-runtime packages to run simulations. You can install them with python -m pip install qiskit-aer qiskit-ibm-runtime.

import numpy as np
from qiskit_experiments.library import StandardRB, InterleavedRB
from qiskit_experiments.framework import ParallelExperiment, BatchExperiment
import qiskit.circuit.library as circuits

# For simulation
from qiskit_aer import AerSimulator
from qiskit_ibm_runtime.fake_provider import FakePerth

backend = AerSimulator.from_backend(FakePerth())

Standard RB experiment#

To run the RB experiment we need to provide the following RB parameters, in order to generate the RB circuits and run them on a backend:

  • qubits: The number of qubits or list of physical qubits for the experiment

  • lengths: A list of RB sequences lengths

  • num_samples: Number of samples to generate for each sequence length

  • seed: Seed or generator object for random number generation. If None then default_rng will be used

  • full_sampling: If True all Cliffords are independently sampled for all lengths. If False for sample of lengths longer sequences are constructed by appending additional Clifford samples to shorter sequences. The default is False

The analysis results of the RB Experiment may include:

  • EPC: The estimated Error Per Clifford

  • alpha: The depolarizing parameter. The fitting function is \(a \cdot \alpha^m + b\), where \(m\) is the Clifford length

  • EPG: The Error Per Gate calculated from the EPC, only for 1-qubit or 2-qubit quantum gates (see [3])

Running a 1-qubit RB experiment#

The standard RB experiment will provide you gate errors for every basis gate constituting an averaged Clifford gate. Note that you can only obtain a single EPC value \(\cal E\) from a single RB experiment. As such, computing the error values for multiple gates \(\{g_i\}\) requires some assumption of contribution of each gate to the total depolarizing error. This is provided by the gate_error_ratio analysis option.

Provided that we have \(n_i\) gates with independent error \(e_i\) per Clifford, the total EPC is estimated by the composition of error from every basis gate,

\[{\cal E} = 1 - \prod_{i} (1 - e_i)^{n_i} \sim \sum_{i} n_i e_i + O(e^2),\]

where \(e_i \ll 1\) and the higher order terms can be ignored.

We cannot distinguish \(e_i\) with a single EPC value \(\cal E\) as explained, however by defining an error ratio \(r_i\) with respect to some standard value \(e_0\), we can compute EPG \(e_i\) for each basis gate.

\[{\cal E} \sim e_0 \sum_{i} n_i r_i\]

The EPG of the \(i\) th basis gate will be

\[e_i \sim r_i e_0 = \dfrac{r_i{\cal E}}{\sum_{i} n_i r_i}.\]

Because EPGs are computed based on this simple assumption, this is not necessarily representing the true gate error on the hardware. If you have multiple kinds of basis gates with unclear error ratio \(r_i\), interleaved RB experiment will always give you accurate error value \(e_i\).

lengths = np.arange(1, 800, 200)
num_samples = 10
seed = 1010
qubits = [0]

# Run an RB experiment on qubit 0
exp1 = StandardRB(qubits, lengths, num_samples=num_samples, seed=seed)
expdata1 = exp1.run(backend).block_for_results()
results1 = expdata1.analysis_results()

# View result data
print("Gate error ratio: %s" % expdata1.experiment.analysis.options.gate_error_ratio)
display(expdata1.figure(0))
for result in results1:
    print(result)
Gate error ratio: {'x': 1.0, 'rz': 0.0, 'sx': 1.0}
../../_images/randomized_benchmarking_1_1.png
AnalysisResult
- name: @Parameters_RBAnalysis
- value: CurveFitResult:
 - fitting method: least_squares
 - number of sub-models: 1
  * F_rb_decay(x) = a * alpha ** x + b
 - success: True
 - number of function evals: 116
 - degree of freedom: 1
 - chi-square: 0.02868487067150859
 - reduced chi-square: 0.02868487067150859
 - Akaike info crit.: -13.750719242535872
 - Bayesian info crit.: -15.5918361591762
 - init params:
  * a = 0.47670899078379886
  * alpha = 0.9997406700260948
  * b = 0.5
 - fit params:
  * a = 0.12548979353786505 ± 0.03932169989755437
  * alpha = 0.9986695726771834 ± 0.0005852345418361386
  * b = 0.8512989596781368 ± 0.039845172383574726
 - correlations:
  * (a, b) = -0.9989613330649482
  * (alpha, b) = -0.9919116315680421
  * (a, alpha) = 0.9876301630518417
- quality: good
- extra: <2 items>
- device_components: ['Q0']
- verified: False
AnalysisResult
- name: alpha
- value: 0.9987+/-0.0006
- χ²: 0.02868487067150859
- quality: good
- extra: <2 items>
- device_components: ['Q0']
- verified: False
AnalysisResult
- name: EPC
- value: 0.00067+/-0.00029
- χ²: 0.02868487067150859
- quality: good
- extra: <2 items>
- device_components: ['Q0']
- verified: False
AnalysisResult
- name: EPG_x
- value: 0.00072+/-0.00032
- χ²: 0.02868487067150859
- quality: good
- extra: <2 items>
- device_components: ['Q0']
- verified: False
AnalysisResult
- name: EPG_rz
- value: 0.0+/-0
- χ²: 0.02868487067150859
- quality: good
- extra: <2 items>
- device_components: ['Q0']
- verified: False
AnalysisResult
- name: EPG_sx
- value: 0.00072+/-0.00032
- χ²: 0.02868487067150859
- quality: good
- extra: <2 items>
- device_components: ['Q0']
- verified: False

Running a 2-qubit RB experiment#

In the same way we can compute EPC for two-qubit RB experiment. However, the EPC value obtained by the experiment indicates a depolarization which is a composition of underlying error channels for 2Q gates and 1Q gates in each qubit. Usually 1Q gate contribution is small enough to ignore, but in case this contribution is significant comparing to the 2Q gate error, we can decompose the contribution of 1Q gates [3].

\[\alpha_{2Q,C} = \frac{1}{5} \left( \alpha_0^{N_1/2} + \alpha_1^{N_1/2} + 3 \alpha_0^{N_1/2} \alpha_1^{N_1/2} \right) \alpha_{01}^{N_2},\]

where \(\alpha_i\) is the single qubit depolarizing parameter of channel \(i\), and \(\alpha_{01}\) is the two qubit depolarizing parameter of interest. \(N_1\) and \(N_2\) are total count of single and two qubit gates, respectively.

Note that the single qubit gate sequence in the channel \(i\) may consist of multiple kinds of basis gates \(\{g_{ij}\}\) with different EPG \(e_{ij}\). Therefore the \(\alpha_i^{N_1/2}\) should be computed from EPGs, rather than directly using the \(\alpha_i\), which is usually a composition of depolarizing maps of every single qubit gate. As such, EPGs should be measured in the separate single-qubit RBs in advance.

\[\alpha_i^{N_1/2} = \alpha_{i0}^{n_{i0}} \cdot \alpha_{i1}^{n_{i1}} \cdot ...,\]

where \(\alpha_{ij}^{n_{ij}}\) indicates a depolarization due to a particular basis gate \(j\) in the channel \(i\). Here we assume EPG \(e_{ij}\) corresponds to the depolarizing probability of the map of \(g_{ij}\), and thus we can express \(\alpha_{ij}\) with EPG.

\[e_{ij} = \frac{2^n - 1}{2^n} (1 - \alpha_{ij}) = \frac{1 - \alpha_{ij}}{2},\]

for the single qubit channel \(n=1\). Accordingly,

\[\alpha_i^{N_1/2} = \prod_{j} (1 - 2 e_{ij})^{n_{ij}},\]

as a composition of depolarization from every primitive gates per qubit. This correction will give you two EPC values as a result of the two-qubit RB experiment. The corrected EPC must be closer to the outcome of interleaved RB. The EPGs of two-qubit RB are analyzed with the corrected EPC if available.

lengths_2_qubit = np.arange(1, 200, 30)
lengths_1_qubit = np.arange(1, 800, 200)
num_samples = 10
seed = 1010
qubits = (1, 2)

# Run a 1-qubit RB experiment on qubits 1, 2 to determine the error-per-gate of 1-qubit gates
single_exps = BatchExperiment(
    [
        StandardRB((qubit,), lengths_1_qubit, num_samples=num_samples, seed=seed)
        for qubit in qubits
    ],
    flatten_results=True,
)
expdata_1q = single_exps.run(backend).block_for_results()
# Run an RB experiment on qubits 1, 2
exp_2q = StandardRB(qubits, lengths_2_qubit, num_samples=num_samples, seed=seed)

# Use the EPG data of the 1-qubit runs to ensure correct 2-qubit EPG computation
exp_2q.analysis.set_options(epg_1_qubit=expdata_1q.analysis_results())

# Run the 2-qubit experiment
expdata_2q = exp_2q.run(backend).block_for_results()

# View result data
print("Gate error ratio: %s" % expdata_2q.experiment.analysis.options.gate_error_ratio)
display(expdata_2q.figure(0))
for result in expdata_2q.analysis_results():
    print(result)
Gate error ratio: {'cx': 1.0}
../../_images/randomized_benchmarking_3_1.png
AnalysisResult
- name: @Parameters_RBAnalysis
- value: CurveFitResult:
 - fitting method: least_squares
 - number of sub-models: 1
  * F_rb_decay(x) = a * alpha ** x + b
 - success: True
 - number of function evals: 20
 - degree of freedom: 4
 - chi-square: 0.5019395620662545
 - reduced chi-square: 0.1254848905165636
 - Akaike info crit.: -12.446299969188175
 - Bayesian info crit.: -12.608569522022236
 - init params:
  * a = 0.7148792344061657
  * alpha = 0.983485233440491
  * b = 0.25
 - fit params:
  * a = 0.7076352351330613 ± 0.006658200700388068
  * alpha = 0.9817628340683558 ± 0.000541832570933822
  * b = 0.25831058308903687 ± 0.005922782827216247
 - correlations:
  * (alpha, b) = -0.8701524047380615
  * (a, b) = -0.7540518607044796
  * (a, alpha) = 0.5351514768336156
- quality: good
- extra: <2 items>
- device_components: ['Q1', 'Q2']
- verified: False
AnalysisResult
- name: alpha
- value: 0.9818+/-0.0005
- χ²: 0.1254848905165636
- quality: good
- extra: <2 items>
- device_components: ['Q1', 'Q2']
- verified: False
AnalysisResult
- name: EPC
- value: 0.0137+/-0.0004
- χ²: 0.1254848905165636
- quality: good
- extra: <2 items>
- device_components: ['Q1', 'Q2']
- verified: False
AnalysisResult
- name: EPC_corrected
- value: 0.0121+/-0.0007
- χ²: 0.1254848905165636
- quality: good
- extra: <2 items>
- device_components: ['Q1', 'Q2']
- verified: False
AnalysisResult
- name: EPG_cx
- value: 0.0078+/-0.0004
- χ²: 0.1254848905165636
- quality: good
- extra: <2 items>
- device_components: ['Q1', 'Q2']
- verified: False

Note that EPC_corrected value is smaller than one of raw EPC, which indicates contribution of depolarization from single-qubit error channels.

Displaying the RB circuits#

The default RB circuit output shows Clifford blocks:

# Run an RB experiment on qubit 0
exp = StandardRB(physical_qubits=(0,), lengths=[2], num_samples=1, seed=seed)
c = exp.circuits()[0]
c.draw(output="mpl", style="iqp")
../../_images/randomized_benchmarking_4_0.png

You can decompose the circuit into underlying gates:

c.decompose().draw(output="mpl", style="iqp")
../../_images/randomized_benchmarking_5_0.png

And see the transpiled circuit using the basis gate set of the backend:

from qiskit import transpile
transpile(c, backend, **vars(exp.transpile_options)).draw(output="mpl", style="iqp", idle_wires=False)
../../_images/randomized_benchmarking_6_0.png

Note

In 0.5.0, the default value of optimization_level in transpile_options changed from 0 to 1 for RB experiments. Transpiled circuits may have less number of gates after the change.

Interleaved RB experiment#

The interleaved RB experiment is used to estimate the gate error of the interleaved gate (see [4]). In addition to the usual RB parameters, we also need to provide:

  • interleaved_element: the element to interleave, given either as a group element or as an instruction/circuit

The analysis results of the RB Experiment includes the following:

  • EPC: The estimated error of the interleaved gate

  • alpha and alpha_c: The depolarizing parameters of the original and interleaved RB sequences respectively

Extra analysis results include

  • EPC_systematic_err: The systematic error of the interleaved gate error [4]

  • EPC_systematic_bounds: The systematic error bounds of the interleaved gate error [4]

Let’s run an interleaved RB experiment on two qubits:

lengths = np.arange(1, 200, 30)
num_samples = 10
seed = 1010
qubits = (1, 2)

# The interleaved gate is the CX gate
int_exp2 = InterleavedRB(
    circuits.CXGate(), qubits, lengths, num_samples=num_samples, seed=seed)

int_expdata2 = int_exp2.run(backend).block_for_results()
int_results2 = int_expdata2.analysis_results()
# View result data
display(int_expdata2.figure(0))
for result in int_results2:
    print(result)
../../_images/randomized_benchmarking_8_0.png
AnalysisResult
- name: @Parameters_InterleavedRBAnalysis
- value: CurveFitResult:
 - fitting method: least_squares
 - number of sub-models: 2
  * F_standard(x) = a * alpha ** x + b
  * F_interleaved(x) = a * (alpha_c * alpha) ** x + b
 - success: True
 - number of function evals: 25
 - degree of freedom: 10
 - chi-square: 6.962276804258971
 - reduced chi-square: 0.6962276804258971
 - Akaike info crit.: -1.7797109440046057
 - Bayesian info crit.: 0.7765183744564279
 - init params:
  * a = 0.715868928301002
  * alpha = 0.9826502782588401
  * b = 0.25
  * alpha_c = 0.9926745948870679
 - fit params:
  * a = 0.7170318652192937 ± 0.003716834913456795
  * alpha = 0.9824792987889924 ± 0.0003264305528128247
  * b = 0.2493127162977085 ± 0.002808438033659369
  * alpha_c = 0.9901050271912178 ± 0.0004943959020440037
 - correlations:
  * (alpha, b) = -0.686999804971292
  * (a, b) = -0.6220651509122922
  * (b, alpha_c) = -0.35628815838422295
  * (alpha, alpha_c) = -0.06697360770949136
  * (a, alpha_c) = 0.1733242327482011
  * (a, alpha) = 0.26923021317291834
- quality: good
- extra: <4 items>
- device_components: ['Q1', 'Q2']
- verified: False
AnalysisResult
- name: alpha
- value: 0.98248+/-0.00033
- χ²: 0.6962276804258971
- quality: good
- extra: <4 items>
- device_components: ['Q1', 'Q2']
- verified: False
AnalysisResult
- name: alpha_c
- value: 0.9901+/-0.0005
- χ²: 0.6962276804258971
- quality: good
- extra: <4 items>
- device_components: ['Q1', 'Q2']
- verified: False
AnalysisResult
- name: EPC
- value: 0.0074+/-0.0004
- χ²: 0.6962276804258971
- quality: good
- extra: <4 items>
- device_components: ['Q1', 'Q2']
- verified: False

References#

See also#