Quantum Tomography¶
TomographyBasis
¶
-
class
qinfer.tomography.
TomographyBasis
(data, dims, labels=None, superrep=None, name=None)[source]¶ Bases:
object
A basis of Hermitian operators used for representing tomographic objects (states and channels) as vectors of real elements. By assumption, a tomographic basis is taken to have an initial (0th) element proportional to the identity, and all other elements are taken to be traceless. For example, the Pauli matrices form a tomographic basis for qubits.
Instances of TomographyBasis convert between representations of tomographic objects as real vectors of model parameters and QuTiP
Qobj
instances. The latter is convienent for working with other libraries, and for reasoning about fidelities and other metrics, while model parameter representations are useful for defining prior distributions and tomographic models.Parameters: - data (np.ndarray) – Dense array of shape
(dim ** 2, dim, dim)
containing all elements of the new tomographic basis.data[alpha, i, j]
is the(i, j)
-th element of thealpha
-th matrix of the new basis. - dims (list) – Dimensions specification used in converting to QuTiP
representations. The product of all elements of
dims
must equal the dimension of axes 1 and 2 ofdata
. For instance,[2, 3]
specifies that the basis is over the tensor product of a qubit and a qutrit space. - labels (
str
orlist
ofstr
) – LaTeX-formatted labels for each basis element. If a singlestr
, a subscript is added to each basis element’s label. - superrep (str) – Superoperator representation to pass to QuTiP when reconstructing states.
-
data
= None¶ Dense matrix... TODO: document indices!
-
labels
= None¶ Labels for each basis element.
-
flat
()[source]¶ Returns a NumPy array that represents this operator basis in a flattened manner, such that
basis.flat()[i, j]
is the \(j\text{th}\) element of the flattened \(i\text{th}\) basis operator.
-
state_to_modelparams
(state)[source]¶ Converts a QuTiP-represented state into a model parameter vector.
Parameters: state (qutip.Qobj) – State to be converted. Return type: np.ndarray
Returns: The representation of the given state in this basis, as a vector of real parameters.
-
modelparams_to_state
(modelparams)[source]¶ Converts one or more vectors of model parameters into QuTiP-represented states.
Parameters: modelparams (np.ndarray) – Array of shape (basis.dim ** 2, )
or(n_states, basis.dim ** 2)
containing states represented as model parameter vectors in this basis.Return type: Qobj
orlist
ofQobj
instances.Returns: The given states represented as Qobj
instances.
- data (np.ndarray) – Dense array of shape
Built-in bases¶
-
qinfer.tomography.
gell_mann_basis
(dim)[source]¶ Returns a
TomographyBasis
on dim dimensions using the generalized Gell-Mann matrices.This implementation is based on a MATLAB-language implementation provided by Carlos Riofrío, Seth Merkel and Andrew Silberfarb. Used with permission.
Parameters: dim (int) – Dimension of the individual matrices making up the returned basis. Return type: TomographyBasis
Returns: A basis of dim * dim
Gell-Mann matrices.
DensityOperatorDistribution
¶
-
class
qinfer.tomography.
DensityOperatorDistribution
(basis)[source]¶ Bases:
qinfer.distributions.SingleSampleMixin
,qinfer.distributions.Distribution
Distribution over density operators parameterized in a given basis.
Parameters: basis ( int
orTomographyBasis
) – Basis to use in representing sampled density operators. If anint
, assumes a default (Gell-Mann) basis of that dimension.-
basis
¶ Basis used to represent sampled density operators as model parameter vectors.
-
Specific Distributions¶
See also
-
class
qinfer.tomography.
TensorProductDistribution
(factors)[source]¶ Bases:
qinfer.tomography.distributions.DensityOperatorDistribution
This class is implemented using QuTiP (v3.1.0 or later), and thus will not work unless QuTiP is installed.
Parameters: factors ( list
ofDensityOperatorDistribution
instances) – Distributions representing each factor of the tensor product used to generate samples.
-
class
qinfer.tomography.
GinibreDistribution
(basis, rank=None)[source]¶ Bases:
qinfer.tomography.distributions.DensityOperatorDistribution
Distribution over all trace-1 positive semidefinite operators of a given rank. Generalizes the Hilbert-Schmidt (full-rank) and Haar (rank-1) distributions.
Parameters: - basis (TomographyBasis) – Basis to use in generating samples.
- rank (int) – Rank of each sampled state. If
None
, defaults to full-rank.
-
class
qinfer.tomography.
GinibreReditDistribution
(basis, rank=None)[source]¶ Bases:
qinfer.tomography.distributions.DensityOperatorDistribution
Distribution over all real-valued trace-1 positive semidefinite operators of a given rank. Generalizes the Hilbert-Schmidt (full-rank) and Haar (rank-1) distributions. Useful for plotting.
Parameters: - basis (TomographyBasis) – Basis to use in generating samples.
- rank (int) – Rank of each sampled state. If
None
, defaults to full-rank.
-
class
qinfer.tomography.
BCSZChoiDistribution
(basis, rank=None, enforce_tp=True)[source]¶ Bases:
qinfer.tomography.distributions.DensityOperatorDistribution
Samples Choi states for completely-positive (CP) or CP and trace-preserving (CPTP) maps, as generated by the BCSZ prior [BCSZ09]. The sampled states are normalized as states (trace 1).
-
class
qinfer.tomography.
GADFLIDistribution
(fiducial_distribution, mean)[source]¶ Bases:
qinfer.tomography.distributions.DensityOperatorDistribution
Samples operators from the generalized amplitude damping prior for liklihood-based inference [GCC16], given a fiducial distribution and the desired mean for the prior.
Parameters: - fiducial_distribution (DensityOperatorDistribution) – Distribution from which samples are initially drawn before transformation under generalized amplitude damping.
- mean (qutip.Qobj) – State which will be the mean of the GAD-transformed samples.
Models¶
-
class
qinfer.tomography.
TomographyModel
(basis, allow_subnormalized=False)[source]¶ Bases:
qinfer.abstract_model.FiniteOutcomeModel
Model for tomographically learning a quantum state using two-outcome positive-operator valued measures (POVMs).
Parameters: - basis (TomographyBasis) – Basis used in representing states as model parameter vectors.
- allow_subnormalized (bool) – If
False
, states \(\rho\) are constrained during resampling such that \(\Tr(\rho) = 1\).
-
basis
¶ Basis used in converting between
Qobj
and model parameter vector representations of states.Type: TomographyBasis
-
n_modelparams
¶
-
modelparam_names
¶
-
is_n_outcomes_constant
¶
-
expparams_dtype
¶
-
canonicalize
(modelparams)[source]¶ Truncates negative eigenvalues and from each state represented by a tensor of model parameter vectors, and renormalizes as appropriate.
Parameters: modelparams (np.ndarray) – Array of shape (n_states, dim**2)
containing model parameter representations of each ofn_states
different states.Returns: The same model parameter tensor with all states truncated to be positive operators. If allow_subnormalized
isFalse
, all states are also renormalized to trace one.
-
trunc_neg_eigs
(particle)[source]¶ Given a state represented as a model parameter vector, returns a model parameter vector representing the same state with any negative eigenvalues set to zero.
Parameters: particle (np.ndarray) – Vector of length (dim ** 2, )
representing a state.Returns: The same state with any negative eigenvalues set to zero.
-
renormalize
(modelparams)[source]¶ Renormalizes one or more states represented as model parameter vectors, such that each state has trace 1.
Parameters: modelparams (np.ndarray) – Array of shape (n_states, dim ** 2)
representing one or more states as model parameter vectors.Returns: The same state, normalized to trace one.
Heuristics¶
Abstract Heuristics¶
-
class
qinfer.tomography.
StateTomographyHeuristic
(updater, basis=None, other_fields=None)[source]¶ Bases:
qinfer.expdesign.Heuristic
-
class
qinfer.tomography.
ProcessTomographyHeuristic
(updater, basis, other_fields=None)[source]¶ Bases:
qinfer.expdesign.Heuristic
-
class
qinfer.tomography.
BestOfKMetaheuristic
(updater, base_heuristic, k=3, other_fields=None)[source]¶ Bases:
qinfer.expdesign.Heuristic
Draws \(k\) different state or tomography measurements, then selects the one that has the largest expected value under the action of the covariance superoperator for the current posterior.
Specific Heuristics¶
-
class
qinfer.tomography.
RandomStabilizerStateHeuristic
(updater, basis=None, other_fields=None)[source]¶ Bases:
qinfer.tomography.expdesign.StateTomographyHeuristic
Randomly chooses rank-1 projectors onto a stabilizer state.
-
class
qinfer.tomography.
RandomPauliHeuristic
(updater, basis=None, other_fields=None)[source]¶ Bases:
qinfer.tomography.expdesign.StateTomographyHeuristic
Randomly chooses a Pauli measurement. Defined for qubits only.
-
class
qinfer.tomography.
ProductHeuristic
(updater, basis, prep_heuristic_class, meas_heuristic_class, other_fields=None)[source]¶ Bases:
qinfer.tomography.expdesign.ProcessTomographyHeuristic
Takes two heuristic classes, one for preparations and one for measurements, then returns a sample from each. The preparation heuristic is assumed to return only trace-1 Hermitian operators.
Plotting Functions¶
-
qinfer.tomography.
plot_cov_ellipse
(cov, pos, nstd=2, ax=None, **kwargs)[source]¶ Plots an
nstd
sigma error ellipse based on the specified covariance matrix (cov
). Additional keyword arguments are passed on to the ellipse patch artist.Parameters: - cov – The 2x2 covariance matrix to base the ellipse on.
- pos – The location of the center of the ellipse. Expects a 2-element
sequence of
[x0, y0]
. - nstd – The radius of the ellipse in numbers of standard deviations. Defaults to 2 standard deviations.
- ax – The axis that the ellipse will be plotted on. Defaults to the current axis.
Returns: A matplotlib ellipse artist.
-
qinfer.tomography.
plot_rebit_prior
(prior, rebit_axes=[1, 2], n_samples=2000, true_state=None, true_size=250, force_mean=None, legend=True, mean_color_index=2)[source]¶ Plots rebit states drawn from a given prior.
Parameters: - prior (qinfer.tomography.DensityOperatorDistribution) – Distribution over rebit states to plot.
- rebit_axes (list) – List containing indices for the \(x\) and \(z\) axes.
- n_samples (int) – Number of samples to draw from the prior.
- true_state (np.ndarray) – State to be plotted as a “true” state for comparison.
-
qinfer.tomography.
plot_rebit_posterior
(updater, prior=None, true_state=None, n_std=3, rebit_axes=[1, 2], true_size=250, legend=True, level=0.95, region_est_method='cov')[source]¶ Plots posterior distributions over rebits, including covariance ellipsoids
Parameters: - updater (qinfer.smc.SMCUpdater) – Posterior distribution over rebits.
- qinfer.tomography.DensityOperatorDistribution – Prior distribution over rebit states.
- true_state (np.ndarray) – Model parameters for “true” state to plot as comparison.
- n_std (float) – Number of standard deviations out from the mean
at which to draw the covariance ellipse. Only used if
region_est_method is
'cov'
. - level (float) – Credibility level to use for computing region estimators from convex hulls.
- rebit_axes (list) – List containing indices for the \(x\) and \(z\) axes.
- region_est_method (str) – Method to use to draw region estimation.
Must be one of None,
'cov'
or'hull'
.