tangelo.toolboxes.operators package

Subpackages

Submodules

tangelo.toolboxes.operators.multiformoperator module

Module defining the MultiformOperator class. It stores a qubit operator in many forms: tuples, numpy array of int and stabilizer formalism. Most of the internal methods use the numpy array of int. The main application for this class is to identify commutation relation faster with the stabilizer notation.

class tangelo.toolboxes.operators.multiformoperator.ConvertPauli(pauli_id)

Bases: object

Helper class to help convert from/to string <-> int <-> stabilizer. It aims to replace multiple dictionaries. It performs the conversion of

I <-> 0 <-> (0, 0) Z <-> 1 <-> (0, 1) X <-> 2 <-> (1, 0) Y <-> 3 <-> (1, 1)

The input can be anything found in the previous conversion table.

Parameters:

pauli_id (string, int or tuple of bool) – Single Pauli.

char

Pauli letter.

Type:

string

integer

Pauli integer.

Type:

int

tuple

Pauli binary representation.

Type:

tuple of bool

class tangelo.toolboxes.operators.multiformoperator.MultiformOperator(terms, n_qubits, factors, integer, binary)

Bases: QubitOperator

Construct integer and binary representations of the operator, as based on the user-specified input. Internal operations are mostly done with an integer array. The conversion is defined as: Letter | Integer | Binary - I | 0 | (0,0) - Z | 1 | (0,1) - X | 2 | (1,0) - Y | 3 | (1,1) Most of the algorithms use the integer representation to detect symmetries and performed operation (ex: multiplication). At the end, the operator can be output as a QubitOperator.

terms

Qubit terms.

Type:

dict

n_qubits

Number of qubits for this operator.

Type:

int

factors

Factors for each term. Without that, information would be lost when converting from a QubitOperator.

Type:

array-like of complex

integer

Array of 0s, 1s, 2s and 3s. Each line represents a term, and each term is a Pauli word.

Type:

array-like of np.int8

binary

Array of 0s (false) and 1s (true). Binary representation of the integer array.

Type:

array-like of bool

binary_swap

Copy of self.binary but with swapped columns.

Type:

array-like of bool

kernel

Null space of the binary representation + identity.

Type:

array-like

Properties:

n_terms (int): Number of terms in this qubit Hamiltonian. qubit_operator (QubitOperator): Self-explanatory.

static collapse(operator, factors)

Function to identify and sum over duplicate terms in an operator, to collapse a set of Pauli words to their minimal representation.

Parameters:
  • operator (array of int) – Operator in integer notation.

  • factors (array of complex) – Self-explanatory.

Returns:

(array of int, arrays of float) – Array of unique integer-encoded Pauli words, their factors in the operator.

compress(abs_tol=None, n_qubits=None)
Overloads the QubitOperator.compress method. It adds an update for

the attributes.

Parameters:
  • abs_tol (float) – Tolerance for the coefficients to be discarded or not. By default, it is EQ_TOLERANCE=1e-8 (in openfermion).

  • n_qubits (int) – Define n_qubits when compressing operator

classmethod from_binaryop(bin_op, factors)

Initialize MultiformOperator from a binary operator.

classmethod from_integerop(int_op, factors)

Initialize MultiformOperator from an integer operator.

classmethod from_qubitop(qubit_op, n_qubits=None)

Initialize MultiformOperator from a qubit operator.

get_kernel()

Get the kernel for a matrix of binary integers. The identity matrix is appended, and the extended matrix is then reduced to column-echelon form. The null space is then extracted.

Returns:

Array-like – 2d numpy array of binary integers.

property n_terms
property qubitoperator
remove_terms(indices)

Remove term(s) from operator. The associated entry is deleted from factors, integer, binary and binary_swap.

Parameters:

indices (int or list of int) – Remove all the terms corresponding to the indices.

tangelo.toolboxes.operators.multiformoperator.do_commute(hybrid_op_a, hybrid_op_b, term_resolved=False)

Check if two operators commute, using the stabilizer notation. With the term_resolved flag, the user can identify which terms in the parent operator do or do not commute with the target operator.

In stabilizer notation, (ax|az)(bz|bx) = -1^(ax.bz + az.bx) (bz|bx)(ax|az) define the commutation relations. By evaluating the sign of (ax.bz + az.bx), we recover commutation relation. For a given pair of terms in the operators, ax.bz + az.bx = sum(XOR (AND (a’ . b))) where a’ is the z-x swapped a-vector. We then apply an OR reduction over all of the b-terms to identify whether term a_i commutes with all of b_j. Applying once again, an OR reduction over all a_i, we get the commutator for the entirety of a,b.

Parameters:
  • hybrid_op_a – First MultiformOperator.

  • hybrid_op_b – Second MultiformOperator.

  • term_resolved (bool) – If True, get commutator for each term.

Returns:

bool or array of bool – If the operators commute or array of bool describing which terms are commuting.

tangelo.toolboxes.operators.multiformoperator.integer_to_binary(integer_op)

Perform conversion from integer array to binary (stabilizer) array.

Parameters:

interger_op (array-like of int) – 2-D numpy array of 0s (I), 1s (Z), 2s (X) and 3s (Y) representing a qubit operator (shape[0] = number of terms, shape[1] = number of qubits.

Returns:

array-like of bool – Array of 0s and 1s representing the operator with the stabilizer notation.

tangelo.toolboxes.operators.multiformoperator.integer_to_qubit_terms(integer_op, factors)

Perform conversion from integer array to qubit terms.

Parameters:
  • interger_op (array-like of int) – 2-D numpy array of 0s (I), 1s (Z), 2s (X) and 3s (Y) representing a qubit operator (shape[0] = number of terms, shape[1] = number of qubits.

  • factors (array-like of float) – Coefficients for the qubit terms.

Returns:

dict – Pauli terms and coefficients.

tangelo.toolboxes.operators.multiformoperator.qubit_to_integer(qubit_op, n_qubits=None)

Perform conversion from qubit operator to integer array. The integer attribute is instantiated, and populated from the QubitOperator attribute.

Parameters:

QubitOperator – self-explanatory.

Returns:

array-like of int – 2-D numpy array of 0s (I), 1s (Z), 2s (X) and 3s (Y) representing a qubit operator (shape[0] = number of terms, shape[1] = number of qubits.

tangelo.toolboxes.operators.operators module

This module defines various kinds of operators used in vqe. It can later be broken down in several modules if needed.

class tangelo.toolboxes.operators.operators.BosonOperator(term=None, coefficient=1.0)

Bases: BosonOperator

Currently, this class is coming from openfermion. Can be later on be replaced by our own implementation.

class tangelo.toolboxes.operators.operators.FermionOperator(term=None, coefficient=1.0, n_spinorbitals=None, n_electrons=None, spin=None)

Bases: FermionOperator

Custom FermionOperator class. Based on openfermion’s, with additional functionalities.

get_coeffs(coeff_threshold=1e-08, spatial=False)

Method to get the coefficient tensors from a fermion operator.

Parameters:
  • coeff_threshold (float) – Ignore coefficient below the threshold. Default value is 1e-8.

  • spatial (bool) – Spatial orbital or spin orbital.

Returns:

(float, array float, array of float) – Core constant, one- (N*N) and two-body coefficient matrices (N*N*N*N), where N is the number of spinorbitals or spatial orbitals.

to_openfermion()

Converts Tangelo FermionOperator to openfermion

class tangelo.toolboxes.operators.operators.QubitHamiltonian(term=None, coefficient=1.0, mapping=None, up_then_down=None)

Bases: QubitOperator

QubitHamiltonian objects are essentially openfermion.QubitOperator objects, with extra attributes. The mapping procedure (mapping) and the qubit ordering (up_then_down) are incorporated into the class. In addition to QubitOperator, several checks are done when performing arithmetic operations on QubitHamiltonians.

term

Same as openfermion term formats.

Type:

openfermion-like

coefficient

Coefficient for this term.

Type:

complex

mapping

Mapping procedure for fermionic to qubit encoding (ex: “JW”, “BK”, etc.).

Type:

string

up_then_down

Whether spin ordering is all up then all down.

Type:

bool

Properties:

n_terms (int): Number of terms in this qubit Hamiltonian.

property n_terms
to_qubitoperator()
class tangelo.toolboxes.operators.operators.QubitOperator(term=None, coefficient=1.0)

Bases: QubitOperator

Currently, this class is coming from openfermion. Can be later on be replaced by our own implementation.

frobenius_norm_compression(epsilon, n_qubits)

Reduces the number of operator terms based on its Frobenius norm and a user-defined threshold, epsilon. The eigenspectrum of the compressed operator will not deviate more than epsilon. For more details, see J. Chem. Theory Comput. 2020, 16, 2, 1055-1063.

Parameters:
  • epsilon (float) – Parameter controlling the degree of compression and resulting accuracy.

  • n_qubits (int) – Number of qubits in the register.

Returns:

QubitOperator – The compressed qubit operator.

classmethod from_openfermion(of_qop)

Enable instantiation of a QubitOperator from an openfermion QubitOperator object.

Parameters:

of_qop (openfermion QubitOperator) – an existing qubit operator defined with Openfermion

Returns:

corresponding QubitOperator object.

get_max_number_hamiltonian_terms(n_qubits)

Compute the possible number of terms for a qubit Hamiltonian. In the absence of an external magnetic field, each Hamiltonian term must have an even number of Pauli Y operators to preserve time-reversal symmetry. See J. Chem. Theory Comput. 2020, 16, 2, 1055-1063 for more details.

Parameters:

n_qubits (int) – Number of qubits in the register.

Returns:

int – The maximum number of possible qubit Hamiltonian terms.

property qubit_indices

Return a set of integers corresponding to qubit indices the qubit operator acts on.

Returns:

set – Set of qubit indices.

to_openfermion()

Converts Tangelo QubitOperator to openfermion

tangelo.toolboxes.operators.operators.count_qubits(qb_op)

Return the number of qubits used by the qubit operator based on the highest index found in the terms.

tangelo.toolboxes.operators.operators.list_to_fermionoperator(all_terms)

Input: a list of terms to generate FermionOperator

Returns:

FermionOperator – Single merged operator.

tangelo.toolboxes.operators.operators.normal_ordered(fe_op)

Input: a Fermionic operator of class toolboxes.operators.FermionicOperator or openfermion.FermionicOperator for reordering.

Returns:

FermionicOperator – Normal ordered operator.

tangelo.toolboxes.operators.operators.qubitop_to_qubitham(qubit_op, mapping, up_then_down)

Function to convert a QubitOperator into a QubitHamiltonian.

Parameters:
  • qubit_op (QubitOperator) – Self-explanatory.

  • mapping (string) – Qubit mapping procedure.

  • up_then_down (bool) – Whether or not spin ordering is all up then all down.

Returns:

QubitHamiltonian – Self-explanatory.

tangelo.toolboxes.operators.operators.squared_normal_ordered(all_terms)

Input: a list of terms to generate toolboxes.operators.FermionOperator or openfermion.FermionOperator

Returns:

FermionOperator – squared (i.e. fe_op*fe_op) and normal ordered.

tangelo.toolboxes.operators.taper_qubits module

Module defining a class to store objects related to qubit tapering, i.e. resource reduction through symmetries.

class tangelo.toolboxes.operators.taper_qubits.QubitTapering(qubit_operator, n_qubits, n_electrons, spin=0, mapping='JW', up_then_down=False)

Bases: object

Class keeping track of symmetries and operation to taper qubits in qubit operators. The taper function is kept as an attribute to taper other Pauli word generators (Hamiltonian dressing, ansatz construction, etc). Z2 tapering is implemented, but other tapering methods could be added in the core of this class.

initial_op

Qubit operator to be analyzed for symmetries.

Type:

MultiformOperator

initial_n_qubits

Number of qubits before tapering.

Type:

int

z2_taper

Function handle for tapering a MultiformOperator.

Type:

func

z2_tapered_op

Tapered operator with z2 symmetries.

Type:

MultiformOperator

z2_properties

Relevant quantities used to define the z2 taper function. Needed to back track a z2 tapered operator to the full operator.

Type:

dict

z2_tapering(qubit_operator, n_qubits=None)

Function to taper a qubit operator from symmetries found in self.initial_op.

Parameters:
  • qubit_operator (QubitOperator) – Self-explanatory.

  • n_qubits (int) – Self-explanatory.

Returns:

QubitOperator – The tapered qubit operator.

tangelo.toolboxes.operators.trim_trivial_qubits module

tangelo.toolboxes.operators.trim_trivial_qubits.is_bitflip_gate(gate: Gate, atol: float = 1e-05) bool

Check if a gate is a bitflip gate.

A gate is a bitflip gate if it satisfies one of the following conditions: 1. The gate name is either X, Y. 2. The gate name is RX or RY, and has a parameter that is an odd multiple of pi.

Parameters:
  • gate (Gate) – The gate to check.

  • atol (float) – Optional, the absolute tolerance for gate parameter

Returns:

bool – True if the gate is a single qubit bitflip gate, False otherwise.

tangelo.toolboxes.operators.trim_trivial_qubits.trim_trivial_circuit(circuit: Circuit) Tuple[Circuit, Dict[int, int]]

Split Circuit into entangled and unentangled components. Returns entangled Circuit, and the indices and states of unentangled qubits

Parameters:

circuit (Circuit) – circuit to be trimmed

Returns:
  • Circuit – Trimmed, entangled circuit

  • dict – dictionary mapping trimmed qubit indices to their states (0 or 1)

tangelo.toolboxes.operators.trim_trivial_qubits.trim_trivial_operator(qu_op: QubitOperator, trim_states: Dict[int, int], n_qubits: None | int = None, reindex: bool = True) QubitOperator

Calculate expectation values of a QubitOperator acting on qubits in a trivial |0> or |1> state. Return a trimmed QubitOperator with updated coefficients

Parameters:
  • qu_op (QubitOperator) – Operator to trim

  • trim_states (dict) – Dictionary mapping qubit indices to states to trim, e.g. {1: 0, 3: 1}

  • n_qubits (int) – Optional, number of qubits in full system

  • reindex (bool) – Optional, if True, remaining qubits will be reindexed

Returns:

QubitOperator – trimmed QubitOperator with updated coefficients

tangelo.toolboxes.operators.trim_trivial_qubits.trim_trivial_qubits(operator: QubitOperator, circuit: Circuit) Tuple[QubitOperator, Circuit]

Trim circuit and operator based on expectation values calculated from trivial components of the circuit.

Parameters:
Returns:
  • QubitOperator – Trimmed qubit operator

  • Circuit – Trimmed circuit

tangelo.toolboxes.operators.z2_tapering module

Module that defines helper functions to taper qubits in a molecular problem with Z2 tapering.

For all chemical Hamiltonians, there are at least two symmetries (electron number and spin conservation) that reduce the qubits count by two. Furthermore, molecular symmetries can lead to a reduction of the number of qubits required to encode a problem. Those symmetries can be interpreted as degenerated eigenvalues. In the real space, symmetry operations lead to the same total energy.

Ref:

Tapering off qubits to simulate fermionic Hamiltonians Sergey Bravyi, Jay M. Gambetta, Antonio Mezzacapo, Kristan Temme. arXiv:1701.08213

tangelo.toolboxes.operators.z2_tapering.get_clifford_operators(kernel)

Function to identify, with a kernel, suitable single-Pauli gates and the related unitary Clifford gates.

Parameters:

kernel (array-like of bool) – Array of M x 2N booleans, where M is the number of terms and N is the number of qubits. Refers to a qubit operator in the stabilizer notation.

Returns:

(list of MultiformOperator, list of int) – Encoded binary-encoded Pauli strings and symmetry indices.

tangelo.toolboxes.operators.z2_tapering.get_eigenvalues(symmetries, n_qubits, n_electrons, spin, mapping, up_then_down)

Get the initial state eigenvalues, as operated on by each of the symmetry operators. These are used to capture the action of each Pauli string in the Hamiltonian on the tapered qubits.

Parameters:
  • symmetries (array-like of bool) – Symmetries in binary encoding.

  • n_qubits (int) – Self-explanatory.

  • n_electrons (int) – Self-explanatory.

  • mapping (str) – Qubit mapping.

  • up_then_down (bool) – Whether or not spin ordering is all up then all down.

Returns:

array of +/-1 – Eigenvalues of an operator with symmetries.

tangelo.toolboxes.operators.z2_tapering.get_unitary(cliffords)

Function generating the product over multiple Clifford operators as a single unitary. The result is a MultiformOperator.

Parameters:

cliffords (list of MultiformOperator) – Encoded Clifford operators.

Returns:

MultiformOperator – Multiplication reflecting the composite operator.

tangelo.toolboxes.operators.z2_tapering.get_z2_taper_function(unitary, kernel, q_indices, n_qubits, n_symmetries, eigenvalues=None)

Defines a method for applying taper to an arbitrary operator in this space. The operator is first conditioned against the tapering-operators: terms which do not commute with the Hamiltonian’s symmetries are culled, and the remaining operator is rotated against the tapering unitary. After applying the eigenvalues to account for tapered elements, the operator is then finally tapered, with the relevant qubits removed.

Parameters:
  • unitary (array of float) – Unitary matrix to perform U*HU.

  • kernel (array of bool) – Kernel representing the NULL space for the operator.

  • q_indices (array of int) – Indices for the relevant columns in the operator array representation.

  • n_qubits (int) – Self-explanatory.

  • n_symmetries (int) – Number of qubits to remove.

  • eigenvalues (array of int) – Initial state eigenvalues for the qubits to be removed.

Returns:

function – function for tapering operator.

Module contents