Tomorun, the main tomography routine (tomographer.tomorun)

Perform a random in the full state space of a quantum system according to our practical, reliable procedure, and collect a histogram of a specific figure of merit.

exception tomographer.tomorun.TomorunInvalidInputError

Bases: tomographer.TomographerCxxError

Exception which gets raised if invalid input is supplied to the tomographer.tomorun.tomorun() function.

tomographer.tomorun.tomorun(dim, ...)

Produce a histogram of a figure of merit during a random walk in quantum state space according to the distribution \(\mu_{B^n}(\cdot)\) defined in Ref. [1]. The likelihood function is specified with independent POVM effects (see below). This python function provides comparable functionality to the tomorun executable program, and allows for a better seamless interoperability with NumPy—all data matrices here are specified as NumPy arrays.

Parameters:
  • dim – The dimension of the quantum system
  • Exn – The observed POVM effects, specified as a matrix in which each row is the X-parameterization of a POVM effect. You may want to specify Emn instead, which may be simpler.
  • Emn – The observed POVM effects, specified as a list of \(\textit{dim}\times\textit{dim}\) matrices.
  • Nm – the list of observed frequency counts for each POVM effect in Emn or Exn.
  • fig_of_merit – The choice of the figure of merit to study. This is either a Python string or a Python callable. If it is a string, it must be one of ‘obs-value’, ‘fidelity’, ‘tr-dist’ or ‘purif-dist’ (see below for more info). If it is a callable, it should accept a single argument, the T-parameterization of the density matrix, and should calculate and return the figure of merit. The T-parameterization is a matrix \(T\) that \(\rho=TT^\dagger\).
  • ref_state – For figures of merit which compare to a reference state (‘fidelity’, ‘tr-dist’, and ‘purif-dist’), this is the reference state to calculate the figure of merit with, specified as a density matrix.
  • observable – For the ‘obs-value’ figure of merit, specify the observable here as a matrix.
  • hist_params – The requested range of values to look at when collecting a histogram of the figure of mert. This should be a tomographer.HistogramParams instance.
  • mhrw_params – The parameters of the random walk, including the step size, the sweep size, the number of thermalization sweeps, and the number of live sweeps. Specify a tomographer.MHRWParams instance here.
  • binning_num_levels – The number of levels in the binning analysis [2]. One should make sure that there are enough bins at the last level to estimate the standard deviation. This is done automatically by default (or if you specify the value -1), so in normal circumstances you won’t have to change the default value.
  • num_repeats – The number of independent random walks to run in parallel.
  • progress_fn – A python callback function to monitor progress. The function should accept a single argument of type tomographer.multiproc.FullStatusReport. Check out tomographer.jpyutil.RandWalkProgressBar if you are using a Jupyter notebook. See below for more information on status progress reporting.
  • progress_interval_ms – The approximate time interval in milliseconds between two progress reports.
  • jumps_method – Method to use for the jumps in the random walk. This may be either “full” or “light”. The “full” method is the one described in the paper, with a jump corresponding to moving the purified bipartite state vector uniformly on the hypersphere. The “light” method is an optimized version, where only an “elementary rotation”—a simple qubit rotation—is applied onto two randomly chosen computational basis elements on the purified bipartite state vector. In the end the random walk explores the same space with the same distribution. The “light” can go much faster especially for large dimensions, but may be slower to converge.
  • ctrl_step_size_params

    A python dict with parameters to set up the controller which dynamically adjusts the step size of the random walk during the thermalization runs. The possible keys are:

    • ’enabled’ (set to ‘True’ or ‘False’): Whether to enable the controller or not. If disabled, the step size is not automatically adjusted.
    • ’desired_accept_ratio_min’, ‘desired_accept_ratio_max’: The range in which we would like to keep the acceptance ratio, by adapting the step size.
    • ’acceptable_accept_ratio_min’, ‘acceptable_accept_ratio_max’: The range of values which the acceptance ratio is not to exceed.
    • ’ensure_n_therm_fixed_params_fraction’: Whenever the step size is adjusted, the controller guarantees that at least this fraction of the given number n_therm of thermalization sweeps is carried out before finishing the thermalization phase.
    • ’num_samples’: The acceptance ratio is computed over this number of previous iterations (deafult 2048). If you notice that the program has trouble tracking the acceptance rate, you might have to increase this value.

    See also this doc of the corresponding C++ controller class.

  • ctrl_converged_params

    A dictionary to set up the controller which dynamically keeps the random walk running while the error bars from binning haven’t converged as required. The possible keys are:

    • ’enabled’ (set to ‘True’ or ‘False’): Whether to enable the controller or not. If disabled, the error bars of the histogram bins will not be checked for convergence before terminating the random walk.
    • ’max_allowed_unknown’, ‘max_allowed_unknown_notisolated’, ‘max_allowed_not_converged’: The maximum allowed number of bins for which the error bars via binning analysis have the respective convergence status. Only after all these requirements are met will the random walk be allowed to finish (or until ‘max_add_run_iters’ faction of run sweeps is exceeded). Default: max_allowed_unknown = 1 + 2% of num_bins, max_allowed_unknown_notisolated = 1 + 1% of bins and max_allowed_not_converged = 1 + .5% of bins.
    • ’check_frequency_sweeps’: How often to check for the convergence of the binning analysis error bars (in number of sweeps).
    • ’max_add_run_iters’: End the random walk after a certain amount runs regardless of bins error bars convergence status. Specify the amount as a fraction of the set number of run sweeps, e.g. a value of 1.5 prolongs the random walk by at most 50% of the run sweeps. Set to a negative value to run as long as necessary to make error bars converge as requested.

    See also this doc of the corresponding C++ controller class.

  • rng_base_seed

    The base seed to use for the random number generator. If None (the default), a base seed is generated based on the current time. Each different run of the random walk is seeded with incremental seeds starting with the base seed (e.g., if the base seed is 910533, the different tasks are given the seeds 910533, 910534, 910535, …). Alternatively, you may specify rng_base_seed=’random-device’ to read a single base seed from a system random device (e.g. /dev/random), or rng_base_seed=’random-device-all’ to read a list of seeds from the random device to use for each individual random walk, instead of using incremental seeds.

    New in version 5.3: Added the rng_base_seed argument

Figures of merit

The value of the fig_of_merit argument may be a Python string, in which case it should be one of the following:

  • “obs-value”: the expectation value of an observable. You should specify the argument observable as a 2-D NumPy array specifying the observable you are interested in.
  • “tr-dist”: the trace distance to a reference state. You should specify the argument ref_state as a 2-D NumPy array specifying the density matrix of the state which should serve as reference state.
  • “fidelity”: the (root) fidelity to a reference state [3]. You should specify the argument ref_state as a 2-D NumPy array specifying the density matrix of the state which should serve as reference state. (For the squared fidelity to a pure reference state, see note below.)
  • “purif-dist”: the purified distance to a reference state [5]. You should specify the argument ref_state as a 2-D NumPy array specifying the density matrix of the state which should serve as reference state.

Note

For the squared fidelity to a pure state (usually preferred in experimental papers), you should use “obs-value” with the observable being the density matrix of the reference state [4].

The value of the fig_of_merit argument may also be a Python callable which directly calculates the figure of merit. It should accept a single argument, the T-parameterization of the density matrix given as a NumPy array (defined such that \(\rho=TT^\dagger\)), and should return the value of the figure of merit. For example, to calculate the purity of the state \(\operatorname{tr}(\rho^2)\):

import numpy as np
import numpy.linalg as npl
...
r = tomographer.tomorun.tomorun(...,
                                fig_of_merit=lambda T: npl.norm(np.dot(T,T.T.conj())),
                                ...)

Return value

The tomorun() function returns a Python dictionary with the following keys and values set:

  • final_histogram: a HistogramWithErrorBars instance with the final histogram data. The histogram has the parameters specified in the hist_params argument. The histogram is NOT normalized to a probabilty density; you should call its normalized() method if you need a normalized histogram.
  • simple_final_histogram: a HistogramWithErrorBars obtained from averaging the raw histograms from each task run, ignoring their error bars from the binning analysis. Under normal circumstances there is no reason you should ignore the binning analysis, so normally you should not be using this member. This member is only useful if you want to test the error bars from the binning analysis against “naive” error bars
  • elapsed_seconds: the total time elapsed while running the random walks, in seconds.
  • final_report_runs: a human-readable summary report of each task run. Allows the user to visually check that all error bars have converged in the binning analysis, and to get an approximate visual representation of what the histogram looks like for each run.
  • final_report: a human-readable summary of the whole procedure. This includes the final report of all the runs contained in final_report_runs, as well as a visual representation of the final averaged histogram.
  • runs_results: a list of all the raw results provided by each task run. Each item of the list is an instance of tomographer.mhrwtasks.MHRandomWalkTaskResult, with its stats_results member being a instance of tomographer.ValueHistogramWithBinningMHRWStatsCollectorResult.

Status reporting

You may receive periodic status reports via a custom Python callback, so that you can stay informed of the overall progress. The callback specified to progress_fn will be called approximately every progress_interval_ms milliseconds with information on the overall progress given as a tomographer.multiproc.FullStatusReport object. Total progress value may exceed 100% to indicate that as many run sweeps as prescribed have been completed, but that bins have not yet converged satisfactorily (see ctrl_converged_params argument).

The individual workers provide the following additional information, formatted within the data dictionary attribute of each WorkerStatusReport object:

  • data['mhrw_params'] – a MHRWParams instance with the current parameters of the random walk
  • data['acceptance_ratio'] – the current acceptance ratio of the Metropolis-Hastings random walk, as a real value between zero and one. You should try to keep this value around ~0.2-0.4. The acceptance ratio is not available during the thermalizing runs.
  • data['kstep'] – the current iteration step number (an iteration corresponds to creating a jump proposal, and to jump with a certain probability)
  • data['n_total_iters'] – the total number of iterations this random walk is going to complete. This is equal to n_sweep*(n_therm + n_run).

Footnotes and references

[1] Christandl and Renner, Phys. Rev. Lett. 12:120403 (2012), arXiv:1108.5329
[2] Ambegaokar and Troyer, Am. J. Phys., 78(2):150 (2010), arXiv:0906.0943
[3] The root fidelity is defined as \(F(\rho,\sigma)=\left\Vert\rho^{1/2}\sigma^{1/2}\right\Vert_1\), as in Nielsen and Chuang, “Quantum Computation and Quantum Information”.
[4] Indeed, for pure \(\rho_\mathrm{ref}\), \(F^2(\rho,\rho_\mathrm{ref}) = \mathrm{tr}(\rho\rho_\mathrm{ref})\).
[5] The purified distance, also called “infidelity” in the literature, is defined as \(P(\rho,\sigma) = \sqrt{1 - F^2(\rho,\sigma)}\).