pygsti.modelmembers.operations.lindbladcoefficients

Module Contents

Classes

LindbladCoefficientBlock

SCRATCH:

Attributes

IMAG_TOL

pygsti.modelmembers.operations.lindbladcoefficients.IMAG_TOL = 1e-07
class pygsti.modelmembers.operations.lindbladcoefficients.LindbladCoefficientBlock(block_type, basis, basis_element_labels=None, initial_block_data=None, param_mode='static', truncate=False)

Bases: pygsti.baseobjs.nicelyserializable.NicelySerializable

SCRATCH: This routine computes the Hamiltonian and Non-Hamiltonian (“other”) superoperator generators which correspond to the terms of the Lindblad expression:

L(rho) = sum_i( h_i [A_i,rho] ) +
sum_ij( o_ij * (B_i rho B_j^dag -

0.5( rho B_j^dag B_i + B_j^dag B_i rho) ) )

where {A_i} and {B_i} are bases (possibly the same) for Hilbert Schmidt (density matrix) space with the identity element removed so that each A_i and B_i are traceless. If we write L(rho) in terms of superoperators H_i and O_ij,

L(rho) = sum_i( h_i H_i(rho) ) + sum_ij( o_ij O_ij(rho) )

then this function computes the matrices for H_i and O_ij using the given density matrix basis. Thus, if dmbasis is expressed in the standard basis (as it should be), the returned matrices are also in this basis.

If these elements are used as projectors it may be usedful to normalize them (by setting normalize=True). Note, however, that these projectors are not all orthogonal - in particular the O_ij’s are not orthogonal to one another.

Parameters

dmbasis_hamlist

A list of basis matrices {B_i} including the identity as the first element, for the returned Hamiltonian-type error generators. This argument is easily obtained by call to pp_matrices() or a similar function. The matrices are expected to be in the standard basis, and should be traceless except for the identity. Matrices should be NumPy arrays or SciPy CSR sparse matrices.

dmbasis_otherlist

A list of basis matrices {B_i} including the identity as the first element, for the returned Stochastic-type error generators. This argument is easily obtained by call to pp_matrices() or a similar function. The matrices are expected to be in the standard basis, and should be traceless except for the identity. Matrices should be NumPy arrays or SciPy CSR sparse matrices.

normalizebool

Whether or not generators should be normalized so that numpy.linalg.norm(generator.flat) == 1.0 Note that the generators will still, in general, be non-orthogonal.

other_mode{“diagonal”, “diag_affine”, “all”}

Which non-Hamiltonian Lindblad error generators to construct. Allowed values are: “diagonal” (only the diagonal Stochastic generators are returned; that is, the generators corresponding to the i==j terms in the Lindblad expression.), “diag_affine” (diagonal + affine generators), and “all” (all generators).

property basis_element_labels
property num_params
property elementary_errorgen_indices

TODO docstring - rewrite this docstring - especially return value! Constructs a dictionary mapping Lindblad term labels to projection coefficients.

This method is used for finding the index of a particular error generator coefficient in the 1D array formed by concatenating the Hamiltonian and flattened stochastic projection arrays.

Parameters
ham_basis{‘std’, ‘gm’, ‘pp’, ‘qt’}, list of matrices, or Basis object

The basis used to construct ham_projs. Allowed values are Matrix-unit (std), Gell-Mann (gm), Pauli-product (pp), and Qutrit (qt), list of numpy arrays, or a custom basis object.

other_basis{‘std’, ‘gm’, ‘pp’, ‘qt’}, list of matrices, or Basis object

The basis used to construct other_projs. Allowed values are Matrix-unit (std), Gell-Mann (gm), Pauli-product (pp), and Qutrit (qt), list of numpy arrays, or a custom basis object.

other_mode{“diagonal”, “diag_affine”, “all”}

Which non-Hamiltonian Lindblad error projections other_projs includes. Allowed values are: “diagonal” (only the diagonal Stochastic), “diag_affine” (diagonal + affine generators), and “all” (all generators).

Returns
Ltermdictdict

Keys are (termType, basisLabel1, <basisLabel2>) tuples, where termType is “H” (Hamiltonian), “S” (Stochastic), or “A” (Affine). Hamiltonian and Affine terms always have a single basis label (so key is a 2-tuple) whereas Stochastic tuples have 1 basis label to indicate a diagonal term and otherwise have 2 basis labels to specify off-diagonal non-Hamiltonian Lindblad terms. Basis labels are taken from ham_basis and other_basis. Values are integer indices.

property elementary_errorgens

Converts a set of coefficients for this block into a linear combination of elementary error generators.

This linear combination is given by a dictionary with keys equal to elementary error generator labels and values equal to their coefficients in the linear combination.

Parameters
block_datanumpy.ndarray

A 1- or 2-dimensional array with each dimension of size len(self.basis_element_labels), specifying the coefficients of this block. Array is 1-dimensional when this block is of type ‘ham’ or ‘other_diagonal’ and is 2-dimensional for type ‘other’.

Returns
elementary_errorgensdict

Specifies block_data as a linear combination of elementary error generators. Keys are LocalElementaryErrorgenLabel objects and values are floats.

property coefficient_labels

Labels for the elements of self.block_data (flattened if relevant)

property param_labels

Generate human-readable labels for the parameters of this block.

Returns
param_labelslist

A list of strings that describe each parameter.

create_lindblad_term_superoperators(mx_basis='pp', sparse='auto', include_1norms=False, flat=False)

Compute the superoperator-generators corresponding to the Lindblad coefficiens in this block. TODO: docstring update

Returns
generatorsnumpy.ndarray or list of SciPy CSR matrices

If parent holds a dense basis, this is an array of shape (n,d,d) for ‘ham’ and ‘other_diagonal’ blocks or (n,n,d,d) for ‘other’ blocks, where d is the dimension of the relevant basis and n is the number of basis elements included in this coefficient block. In the case of ‘ham’ and ‘other_diagonal’ type blocks, generators[i] gives the superoperator matrix corresponding to the block’s i-th coefficient (and i-th basis element). In the ‘other’ case, generators[i,j] corresponds to the (i,j)-th coefficient. If the parent holds a sparse basis, the dimensions of size n are lists of CSR matrices with shape (d,d).

create_lindblad_term_objects(parameter_index_offset, max_polynomial_vars, evotype, state_space)
set_elementary_errorgens(elementary_errorgens, on_missing='ignore', truncate=False)
set_from_errorgen_projections(errorgen, errorgen_basis='pp', return_projected_errorgen=False, truncate=False)
to_vector()

Compute parameter values for this coefficient block.

Returns
param_valsnumpy.ndarray

A 1D array of real parameter values. Length is len(self.basis_element_labels) in the case of ‘ham’ or ‘other_diagonal’ blocks, and len(self.basis_element_labels)**2 in the case of ‘other’ blocks.

from_vector(v)

Construct Lindblad coefficients (for this block) from a set of parameter values.

This function essentially performs the inverse of coefficients_to_paramvals().

Parameters
vnumpy.ndarray

A 1D array of real parameter values.

deriv_wrt_params(v=None)

Construct derivative of Lindblad coefficients (for this block) from a set of parameter values.

This function gives the Jacobian of what is returned by paramvals_to_coefficients() (as a function of the parameters).

Parameters
vnumpy.ndarray, optional

A 1D array of real parameter values. If not specified, then self.to_vector() is used.

Returns
block_data_derivnumpy.ndarray

A real array of shape (nBEL,nP) or (nBEL,nBEL,nP), depending on the block type, where nBEL is this block’s number of basis elements (see self.basis_element_labels) and nP is the number of parameters (the length of parameter_values).

elementary_errorgen_deriv_wrt_params(v=None)
superop_deriv_wrt_params(superops, v=None, superops_are_flat=False)

TODO: docstring

superopsnumpy.ndarray

Output of create_lindblad_term_superoperators (with flat=True if superops_are_flat==True), so that this is a 3- or 4-dimensional array indexed by (iSuperop, superop_row, superop_col) or (iSuperop1, iSuperop2, superop_row, superop_col).

Returns
numpy.ndarray

per-superop-element derivative, indexed by (superop_row, superop_col, parameter_index) or (superop_row, superop_col, parameter_index1, parameter_index2) where there are two parameter indices because parameters are indexed by an (i,j) pair rather than a single index.

superop_hessian_wrt_params(superops, v=None, superops_are_flat=False)

TODO: docstring

Returns
numpy.ndarray

Indexed by (superop_row, superop_col, param1, param2).

is_similar(other_coeff_block)

TODO: docstring

convert(param_mode)

TODO: docstring - return a new LindbladCoefficientBlock with the same block type and data, but with the given parameterization mode.