pygsti.models.explicitmodel
Defines the ExplicitOpModel class and supporting functionality.
Module Contents
Classes
Encapsulates a set of gate, state preparation, and POVM effect operations. |
|
Rule: layer must appear explicitly as a "primitive op" |
Functions
|
Return a copy of mdl whose members have been gauge-transformed by s, |
- class pygsti.models.explicitmodel.ExplicitOpModel(state_space, basis='pp', default_gate_type='full', default_prep_type='auto', default_povm_type='auto', default_instrument_type='auto', prep_prefix='rho', effect_prefix='E', gate_prefix='G', povm_prefix='M', instrument_prefix='I', simulator='auto', evotype='default')
Bases:
pygsti.models.model.OpModelEncapsulates a set of gate, state preparation, and POVM effect operations.
An ExplictOpModel stores a set of labeled LinearOperator objects and provides dictionary-like access to their matrices. State preparation and POVM effect operations are represented as column vectors.
Parameters
- state_spaceStateSpace
The state space for this model.
- basis{“pp”,”gm”,”qt”,”std”,”sv”} or Basis, optional
The basis used for the state space by dense superoperator representations.
- default_param{“full”, “TP”, “CPTP”, etc.}, optional
Specifies the default gate and SPAM vector parameterization type. Can be any value allowed by
set_all_parameterizations(), which also gives a description of each parameterization type.- prep_prefix: string, optional
Key prefixe for state preparations, allowing the model to determing what type of object a key corresponds to.
- effect_prefixstring, optional
Key prefix for POVM effects, allowing the model to determing what type of object a key corresponds to.
- gate_prefixstring, optional
Key prefix for gates, allowing the model to determing what type of object a key corresponds to.
- povm_prefixstring, optional
Key prefix for POVMs, allowing the model to determing what type of object a key corresponds to.
- instrument_prefixstring, optional
Key prefix for instruments, allowing the model to determing what type of object a key corresponds to.
- simulatorForwardSimulator or {“auto”, “matrix”, “map”}
The circuit simulator used to compute any requested probabilities, e.g. from
probs()orbulk_probs(). The default value of “auto” automatically selects the simulation type, and is usually what you want. Other special allowed values are:“matrix” : op_matrix-op_matrix products are computed and cached to get composite gates which can then quickly simulate a circuit for any preparation and outcome. High memory demand; best for a small number of (1 or 2) qubits.
“map” : op_matrix-state_vector products are repeatedly computed to simulate circuits. Slower for a small number of qubits, but faster and more memory efficient for higher numbers of qubits (3+).
- evotypeEvotype or str, optional
The evolution type of this model, describing how states are represented. The special value “default” is equivalent to specifying the value of pygsti.evotypes.Evotype.default_evotype.
Creates a new OpModel. Rarely used except from derived classes __init__ functions. Parameters ———- state_space : StateSpace
The state space for this model.
- basisBasis
The basis used for the state space by dense operator representations.
- evotypeEvotype or str, optional
The evolution type of this model, describing how states are represented. The special value “default” is equivalent to specifying the value of pygsti.evotypes.Evotype.default_evotype.
- layer_rulesLayerRules
The “layer rules” used for constructing operators for circuit layers. This functionality is essential to using this model to simulate ciruits, and is typically supplied by derived classes.
- simulatorForwardSimulator or {“auto”, “matrix”, “map”}
The forward simulator (or typ) that this model should use. “auto” tries to determine the best type automatically.
- property default_gauge_group
Gets the default gauge group for performing gauge transformations on this Model.
Returns
GaugeGroup
- property prep
The unique state preparation in this model, if one exists.
If not, a ValueError is raised.
Returns
State
- property effects
The effect vectors from the unique POVM in this model, if one exists.
If not, a ValueError is raised.
Returns
list of POVMEffects
- property num_elements
Return the number of total operation matrix and spam vector elements in this model.
This is in general different from the number of parameters in the model, which are the number of free variables used to generate all of the matrix and vector elements.
Returns
- int
the number of model elements.
- property num_nongauge_params
Return the number of non-gauge parameters in this model.
Returns
- int
the number of non-gauge model parameters.
- property num_gauge_params
Return the number of gauge parameters in this model.
Returns
- int
the number of gauge model parameters.
- preps
- povms
- operations
- instruments
- factories
- effects_prefix = "'E'"
- convert_members_inplace(to_type, categories_to_convert='all', labels_to_convert='all', ideal_model=None, flatten_structure=False, set_default_gauge_group=False, cptp_truncation_tol=1e-07, spam_cp_penalty=1e-07)
Method for converting the parameterizations of modelmembers within this model to new ones in-place.
Parameters
- to_typestr
String specifier for the parameterization type to convert to.
- categories_to_convertstr or list of str, optional (default ‘all’)
Categories of modelmembers to perform conversion on. Allowed options are: ‘all’, ‘ops’ or ‘operations’ (these two are aliases for the same option), ‘instruments’, ‘povms’ or ‘preps’.
- labels_to_convertstr or list of Label, optional (default ‘all’)
A string specifier, or list of Label objects, specifying the set of objects (state preparations, operations, instruments, etc.) within the model to apply the conversion to.
- ideal_modelModel, optional (default None)
A model containing modelmembers instantiated such that they all correspond to the ideal actions of the given gate set elements. It is recommended that this be specified when converting to an error-generator-based parameterization.
- flatten_structurebool, optional (default False)
When False, the sub-members of composed and embedded operations are separately converted, leaving the original modelmember structure unchanged. When True, composed and embedded operations are “flattened” into a single modelmember parameterized according to the requested to_type.
- set_default_gauge_groupbool, optional (default False)
A flag specifying whether the default gauge group for the model should be updated to the default value associated with the specified value of to_type. See set_default_gauge_group_for_member_type for more on these default gauge groups.
- cptp_truncation_tolfloat, optional (default 1e-7)
Tolerance term used to enforce the CPTP constraint on gates when moving between different parameterizations.
- spam_cp_penaltyfloat, optional (default 1e-7)
Penalty term used to enforce the CP constraint on SPAM when moving between different parameterizations.
- set_default_gauge_group_for_member_type(member_type)
Updates the default gauge group to the default value for the specified modelmember type.
Parameters
- member_typestr
A string specifier for the modelmember type used to select the gauge group type. Mapping is the following:
‘full’ -> FullGaugeGroup
‘full TP’, ‘TP’, TPGaugeGroup
‘CPTP’ or Anything that is a valid lindblad type -> UnitaryGaugeGroup
Otherwise -> TrivialGaugeGroup
- set_all_parameterizations(gate_type, prep_type='auto', povm_type='auto', instrument_type='auto', ideal_model=None, cptp_truncation_tol=1e-06, spam_cp_penalty=1e-07)
Convert all gates, states, and POVMs to a specific parameterization type.
Parameters
- parameterization_typestring
The gate, state, and POVM parameterization type. Allowed values are (where ‘*’ means “ terms” and “ clifford terms” evolution-type suffixes are allowed):
“full” : each gate / state / POVM effect element is an independent parameter
“TP” : Trace-Preserving gates and state preps
“static” : no parameters
“static unitary” : no parameters; convert superops to unitaries
“clifford” : no parameters; convert unitaries to Clifford symplecitics.
“GLND*” : General unconstrained Lindbladian
“CPTP*” : Completely-Positive-Trace-Preserving
“H+S+A*” : Hamiltoian, Pauli-Stochastic, and Affine errors
“H+S*” : Hamiltonian and Pauli-Stochastic errors
“S+A*” : Pauli-Stochastic and Affine errors
“S*” : Pauli-Stochastic errors
“H+D+A*” : Hamiltoian, Depolarization, and Affine errors
“H+D*” : Hamiltonian and Depolarization errors
“D+A*” : Depolarization and Affine errors
“D*” : Depolarization errors
Any of the above with “S” replaced with “s” or “D” replaced with “d”. This removes the CPTP constraint on the gates and SPAM operations (and as such is seldom used).
- ideal_modelModel, optional
This may specify an ideal model of unitary gates and pure state vectors to be used as the ideal operation of each gate/SPAM operation, which is particularly useful as target for CPTP-based conversions.
- cptp_truncation_tolfloat, optional (default 1e-6)
Tolerance used for conversion to CPTP parameterizations. When converting to CPTP models negative eigenvalues of the choi matrix representation of a superoperator are truncated, which can result in a change in the PTM for that operator. This tolerance indicates the maximum amount of truncation induced deviation from the original operations (measured by frobenius distance) we’re willing to accept without marking the conversion as failed.
- spam_cp_penaltyfloat, optional (default 0.5)
Converting SPAM operations to an error generator representation may introduce trivial gauge degrees of freedom. These gauge degrees of freedom are called trivial because they quite literally do not change the dense representation (i.e. Hilbert-Schmidt vectors) at all. Despite being trivial, error generators along this trivial gauge orbit may be non-CP, so this cptp penalty is used to favor channels within this gauge orbit which are CPTP.
Returns
None
- deriv_wrt_params()
The element-wise derivative of all this models’ operations.
Constructs a matrix whose columns are the vectorized derivatives of all the model’s raw matrix and vector elements (placed in a vector) with respect to each single model parameter.
Thus, each column has length equal to the number of elements in the model, and there are num_params() columns. In the case of a “fully parameterized model” (i.e. all operation matrices and SPAM vectors are fully parameterized) then the resulting matrix will be the (square) identity matrix.
Returns
- numpy array
2D array of derivatives.
- compute_nongauge_and_gauge_spaces(item_weights=None, non_gauge_mix_mx=None)
TODO: docstring
- compute_nongauge_projector(item_weights=None, non_gauge_mix_mx=None)
Construct a projector onto the non-gauge parameter space.
Useful for isolating the gauge degrees of freedom from the non-gauge degrees of freedom.
Parameters
- item_weightsdict, optional
Dictionary of weighting factors for individual gates and spam operators. Keys can be gate, state preparation, POVM effect, spam labels, or the special strings “gates” or “spam” which represent the entire set of gate or SPAM operators, respectively. Values are floating point numbers. These weights define the metric used to compute the non-gauge space, orthogonal the gauge space, that is projected onto.
- non_gauge_mix_mxnumpy array, optional
An array of shape (n_non_gauge_params,n_gauge_params) specifying how to mix the non-gauge degrees of freedom into the gauge degrees of freedom that are projected out by the returned object. This argument essentially sets the off-diagonal block of the metric used for orthogonality in the “gauge + non-gauge” space. It is for advanced usage and typically left as None (the default).
Returns
- numpy array
The projection operator as a N x N matrix, where N is the number of parameters (obtained via num_params()). This projector acts on parameter-space, and has rank equal to the number of non-gauge degrees of freedom.
- transform_inplace(s)
Gauge transform this model.
Update each of the operation matrices G in this model with inv(s) * G * s, each rhoVec with inv(s) * rhoVec, and each EVec with EVec * s
Parameters
- sGaugeGroupElement
A gauge group element which specifies the “s” matrix (and it’s inverse) used in the above similarity transform.
Returns
None
- frobeniusdist(other_model, transform_mx=None, item_weights=None, normalize=True)
Compute the weighted frobenius norm of the difference between this model and other_model.
Differences in each corresponding gate matrix and spam vector element are squared, weighted (using item_weights as applicable), then summed. The value returned is the square root of this sum, or the square root of this sum divided by the number of summands if normalize == True.
Parameters
- other_modelModel
the other model to difference against.
- transform_mxnumpy array, optional
if not None, transform this model by G => inv(transform_mx) * G * transform_mx, for each operation matrix G (and similar for rho and E vectors) before taking the difference. This transformation is applied only for the difference and does not alter the values stored in this model.
- item_weightsdict, optional
Dictionary of weighting factors for individual gates and spam operators. Weights are applied multiplicatively to the squared differences, i.e., (before the final square root is taken). Keys can be gate, state preparation, POVM effect, or spam labels, as well as the two special labels “gates” and “spam” which apply to all of the gate or SPAM elements, respectively (but are overridden by specific element values). Values are floating point numbers. By default, all weights are 1.0.
- normalizebool, optional
if True (the default), the sum of weighted squared-differences is divided by the weighted number of differences before the final square root is taken. If False, the division is not performed.
Returns
float
- residuals(other_model, transform_mx=None, item_weights=None)
Compute the weighted residuals between two models.
Residuals are the differences in corresponding operation matrix and spam vector elements.
Parameters
- other_modelModel
the other model to difference against.
- transform_mxnumpy array, optional
if not None, transform this model by G => inv(transform_mx) * G * transform_mx, for each operation matrix G (and similar for rho and E vectors) before taking the difference. This transformation is applied only for the difference and does not alter the values stored in this model.
- item_weightsdict, optional
Dictionary of weighting factors for individual gates and spam operators. Weights applied such that they act multiplicatively on the squared differences, so that the residuals themselves are scaled by the square roots of these weights. Keys can be gate, state preparation, POVM effect, or spam labels, as well as the two special labels “gates” and “spam” which apply to all of the gate or SPAM elements, respectively (but are overridden by specific element values). Values are floating point numbers. By default, all weights are 1.0.
Returns
- residualsnumpy.ndarray
A 1D array of residuals (differences w.r.t. other)
- nSummandsint
The (weighted) number of elements accounted for by the residuals.
- jtracedist(other_model, transform_mx=None, include_spam=True)
Compute the Jamiolkowski trace distance between this model and other_model.
This is defined as the maximum of the trace distances between each corresponding gate, including spam gates.
Parameters
- other_modelModel
the other model to difference against.
- transform_mxnumpy array, optional
if not None, transform this model by G => inv(transform_mx) * G * transform_mx, for each operation matrix G (and similar for rho and E vectors) before taking the difference. This transformation is applied only for the difference and does not alter the values stored in this model.
- include_spambool, optional
Whether to add to the max-trace-distance the frobenius distances between corresponding SPAM operations.
Returns
float
- diamonddist(other_model, transform_mx=None, include_spam=True)
Compute the diamond-norm distance between this model and other_model.
This is defined as the maximum of the diamond-norm distances between each corresponding gate, including spam gates.
Parameters
- other_modelModel
the other model to difference against.
- transform_mxnumpy array, optional
if not None, transform this model by G => inv(transform_mx) * G * transform_mx, for each operation matrix G (and similar for rho and E vectors) before taking the difference. This transformation is applied only for the difference and does not alter the values stored in this model.
- include_spambool, optional
Whether to add to the max-diamond-distance the frobenius distances between corresponding SPAM operations.
Returns
float
- strdiff(other_model, metric='frobenius')
Return a string describing the distances between this model and other_model.
The returned string displays differences between each corresponding gate, state prep, and POVM effect.
Parameters
- other_modelModel
the other model to difference against.
- metric{‘frobenius’, ‘infidelity’, ‘diamond’}
Which distance metric to use.
Returns
str
- all_objects()
Iterate over all of the (label, operator object) entities in this model.
This iterator runs over all state preparations, POVMS, operations, and instruments.
- depolarize(op_noise=None, spam_noise=None, max_op_noise=None, max_spam_noise=None, seed=None)
Apply depolarization uniformly or randomly to this model’s gate and/or SPAM elements.
The result is returned without modifying the original (this) model. You must specify either op_noise or max_op_noise (for the amount of gate depolarization), and either spam_noise or max_spam_noise (for spam depolarization).
Parameters
- op_noisefloat, optional
- apply depolarizing noise of strength
1-op_noiseto all gates in the model. (Multiplies each assumed-Pauli-basis operation matrix by the diagonal matrix with
(1.0-op_noise)along all the diagonal elements except the first (the identity).
- apply depolarizing noise of strength
- spam_noisefloat, optional
apply depolarizing noise of strength
1-spam_noiseto all SPAM opeations (state and POVM effects) in the model. (Multiplies the non-identity part of each assumed-Pauli-basis state preparation vector and measurement vector by(1.0-spam_noise)).- max_op_noisefloat, optional
specified instead of op_noise; apply a random depolarization with maximum strength
1-max_op_noiseto each gate in the model.- max_spam_noisefloat, optional
specified instead of spam_noise; apply a random depolarization with maximum strength
1-max_spam_noiseto each state preparation and POVM in the model.- seedint, optional
if not
None, seed numpy’s random number generator with this value before generating random depolarizations.
Returns
- Model
the depolarized Model
- rotate(rotate=None, max_rotate=None, seed=None)
Apply a rotation uniformly or randomly to this model.
Uniformly means the same rotation applied to each gate and randomly means different random rotations are applied to each gate of this model. The result is returned without modifying the original (this) model.
You must specify either rotate or max_rotate. This method currently only works on n-qubit models.
Parameters
- rotatetuple of floats, optional
If you specify the rotate argument, then the same rotation operation is applied to each gate. That is, each gate’s matrix G is composed with a rotation operation R (so G -> dot(R, G) ) where R is the unitary superoperator corresponding to the unitary operator U = exp( sum_k( i * rotate[k] / 2.0 * Pauli_k ) ). Here Pauli_k ranges over all of the non-identity un-normalized Pauli operators (e.g. {X,Y,Z} for 1 qubit, {IX, IY, IZ, XI, XX, XY, XZ, YI, YX, YY, YZ, ZI, ZX, ZY, ZZ} for 2 qubits).
- max_rotatefloat, optional
If max_rotate is specified (instead of rotate), then pyGSTi randomly generates a different rotate tuple, and applies the corresponding rotation, to each gate in this Model. Each component of each tuple is drawn uniformly from [0, max_rotate).
- seedint, optional
if not None, seed numpy’s random number generator with this value before generating random depolarizations.
Returns
- Model
the rotated Model
- randomize_with_unitary(scale, seed=None, rand_state=None)
Create a new model with random unitary perturbations.
Apply a random unitary to each element of a model, and return the result, without modifying the original (this) model. This method works on Model as long as the dimension is a perfect square.
Parameters
- scalefloat
maximum element magnitude in the generator of each random unitary transform.
- seedint, optional
if not None, seed numpy’s random number generator with this value before generating random depolarizations.
- rand_statenumpy.random.RandomState
A RandomState object to generate samples from. Can be useful to set instead of seed if you want reproducible distribution samples across multiple random function calls but you don’t want to bother with manually incrementing seeds between those calls.
Returns
- Model
the randomized Model
- increase_dimension(new_dimension)
Enlarge the dimension of this model.
Enlarge the spam vectors and operation matrices of model to a specified dimension, and return the resulting inflated model. Spam vectors are zero-padded and operation matrices are padded with 1’s on the diagonal and zeros on the off-diagonal (effectively padded by identity operation).
Parameters
- new_dimensionint
the dimension of the returned model. That is, the returned model will have rho and E vectors that have shape (new_dimension,1) and operation matrices with shape (new_dimension,new_dimension)
Returns
- Model
the increased-dimension Model
- kick(absmag=1.0, bias=0, seed=None)
“Kick” this model by adding to each gate a random matrix.
The random matrices have values uniformly distributed in the interval [bias-absmag,bias+absmag].
Parameters
- absmagfloat, optional
The maximum magnitude of the entries in the “kick” matrix relative to bias.
- biasfloat, optional
The bias of the entries in the “kick” matrix.
- seedint, optional
if not None, seed numpy’s random number generator with this value before generating random depolarizations.
Returns
- Model
the kicked model.
- compute_clifford_symplectic_reps(oplabel_filter=None)
Constructs a dictionary of the symplectic representations for all the Clifford gates in this model.
Non-
StaticCliffordOpgates will be ignored and their entries omitted from the returned dictionary.Parameters
- oplabel_filteriterable, optional
A list, tuple, or set of operation labels whose symplectic representations should be returned (if they exist).
Returns
- dict
keys are operation labels and/or just the root names of gates (without any state space indices/labels). Values are (symplectic_matrix, phase_vector) tuples.
- print_info()
Print to stdout relevant information about this model.
This information includes the Choi matrices and their eigenvalues.
Returns
None
- create_processor_spec(qudit_labels='auto')
Create a processor specification from this model with the given qudit labels.
Currently this only works for models on qudits.
Parameters
- qudit_labelstuple or “auto”, optional
A tuple of qudit labels, e.g. (‘Q0’, ‘Q1’) or (0, 1). “auto” uses the labels in this model’s state space labels.
Returns
QuditProcessorSpec or QubitProcessorSpec
- create_modelmember_graph()
Generate a ModelMemberGraph for the model.
Returns
- ModelMemberGraph
A directed graph capturing dependencies among model members
- errorgen_coefficients(normalized_elem_gens=True)
TODO: docstring - returns a nested dict containing all the error generator coefficients for all the operations in this model.
- pygsti.models.explicitmodel.transform_composed_model(mdl: pygsti.models.ExplicitOpModel, s: pygsti.models.gaugegroup.GaugeGroupElement) pygsti.models.ExplicitOpModel
Return a copy of mdl whose members have been gauge-transformed by s, while retaining the parameterization of mdl.
This function’s implementation requires that mdl use ComposedState for stateprep and ComposedPOVM for measurements. It does NOT require that operations be represented with ComposedOp. It ignores any factories that might be present in mdl.
- class pygsti.models.explicitmodel.ExplicitLayerRules
Bases:
pygsti.models.layerrules.LayerRulesRule: layer must appear explicitly as a “primitive op”
- prep_layer_operator(model, layerlbl, caches)
Create the operator corresponding to layerlbl.
Parameters
- layerlblLabel
A circuit layer label.
Returns
State