Helper functions¶
-
helper.
Fidelity
(psi, phi)[source]¶ Returns the fidelity bewteen two quantum states psi, phi defined as |<psi|phi>|^2
- Parameters
psi (complex np.array or quimb.core.qarray) -- Quantum state
phi (complex np.array or quimb.core.qarray) -- Quantum state
- Returns
Fidelity -- Fidelity of the two quantum states
- Return type
double real
-
helper.
get_figsize
(wf=0.5, hf=0.6180339887498949)[source]¶ - Parameters
[float] (columnwidth) --
[float] -- Set by default to golden ratio.
[float] -- using showthecolumnwidth , optional
- Returns
fig_size -- fig_size [width, height] that should be given to matplotlib
- Return type
list of float
-
helper.
left_contract
(states)[source]¶ Given the N left-most states of a MPS, computes their contraction with themselves, as in the Example.
- Parameters
states (list of tensors) -- N left-most states of a MPS
- Returns
Matrix -- Returns a order 2 tensor.
- Return type
order 2 tensor
Examples
>>> o - o - o - o - >>> | | | | >>> o.H - o.H - o.H - o.H -
-
helper.
print_state
(dense_state)[source]¶ Prints a dense_state with kets. Compatible with quimb states.
- Parameters
dense_state (array_like) -- Dense representation of a quantum state
- Returns
None
- Return type
None
-
helper.
qiskit_get_statevect
(qc)[source]¶ Returns the statevector of the qiskit quantum circuit qc in the reversed order
- Parameters
qc (Quantum circuit) -- Quantum circuit of which we want the statevector
- Returns
st -- Statevector of the quantum circuit after the application of the reverse operation on the qubit's ordering
- Return type
array_like
-
helper.
reverse_qiskit
(qiskit_dense, N)[source]¶ Reverse the ordering of a qiskit state, in which the 0-th qubit is the leftmost (for example in |1000> the 0-th qubit is in the state 1), to match the more usual notation in which the 0-th qubit is the rightmost.
- qiskit_dense: array_like
Dense qiskit state vector
- N: int
Number of qubits which forms qiskit_dense
- Returns
qiskit_reordered -- Reordered dense state_vector
- Return type
array_like
Examples
>>> qiskit_result = [0, 0, 1, 0] # |01>, but usually recognized as |10> >>> reverse_qiskit(qiskit_result, 2) [0, 1, 0, 0] # |01> with the usual representation
-
helper.
right_contract
(states)[source]¶ Given the N right-most states of a MPS, computes their contraction with themselves, as in the Example
- Parameters
states (list of tensors) -- N right-most states of a MPS
- Returns
Matrix -- Returns a order 2 tensor.
- Return type
order 2 tensor
Examples
>>> - o - o - o - o >>> | | | | >>> o.H - o.H - o.H - o.H
-
helper.
to_approx_MPS
(dense_state, N, d=2, chi=2)[source]¶ Converts a dense_state of a N-body system made by d-dimensional sites into a Matrix Product State in left-canonical form, with the size of links bounded by chi.
- Parameters
dense_state (ndarray of shape (d^N,)) -- Input dense state, such that the (i,j,k...) entry in dense_state.reshape([d]*N) is the (i,j,k...) coefficient of the state in the computational basis.
N (integer > 0) -- Number of particles/sites
d (integer > 0, optional) -- Local dimension of each particle/site. For a qubit, d=2.
chi (integer > 0, optional) -- Maximum bond dimension
- Returns
MPS -- List of N tensors containing the left-canonical MPS. The first and last tensors are of order 2 (matrices), while all the others are of order 3. The shapes are not fixed, but they are (a_i, d, a_{i+1}), with a_i, a_{i+1} <= chi for the order 3 tensors, and (d, a_1) or (a_{N-1}, d) for the order 2 tensors at the boundaries. We can see the representation in the first example. The index ordering convention is from left-to-right. For instance, the "left" index of U2 is the first, the "bottom" one is the second, and the "right" one is the third.
- Return type
List of N tensors containing the left-canonical MPS
Examples
>>> U1 - U2 - U3 - ... - UN >>> | | | | # For d=2, N=7 and chi=5, the tensor network is as follows: >>> U1 -2- U2 -4- U3 -5- U4 -5- U5 -4- U6 -2- U7 >>> | | | | | | | # where -x- denotes the bounds' dimension (all the "bottom-facing" indices are of dimension d=2). Thus, the shapes # of the returned tensors are as follows: >>> U1 U2 U3 U4 U5 U6 U7 >>> [(2, 2), (2, 2, 4), (4, 2, 5), (5, 2, 5), (5, 2, 4), (4, 2, 2), (2, 2)]
-
helper.
to_dense
(MPS)[source]¶ Given a list of N tensors MPS [U1, U2, ..., UN] , representing a Matrix Product State, perform the contraction in the Examples, leading to a single tensor of order N, representing a dense state.
The index ordering convention is from left-to-right. For instance, the "left" index of U2 is the first, the "bottom" one is the second, and the "right" one is the third.
- Parameters
MPS (list of ndarrays) -- List of tensors. First and last should be of order 2, all the others of order 3. The last dimension of MPS[i] should be the same of the first dimension of MPS[i+1], for all i.
- Returns
dense_state -- N-order tensor representing the dense state.
- Return type
ndarray of shape ([d] ^ N)
Examples
>>> U1 - U2 - ... - UN >>> | | |
-
helper.
to_full_MPS
(dense_state, N, d=2)[source]¶ Converts a dense_state of a N-body system made by d-dimensional sites into a Matrix Product State in left-canonical form, with sufficiently sized bonds so that exactness is maintained.
- Parameters
dense_state (ndarray of shape (d^N,)) -- Input dense state, such that the (i,j,k...) entry in dense_state.reshape([d]^N) is the (i,j,k...) coefficient of the state in the computational basis.
N (integer > 0) -- Number of particles/sites
d (integer > 0, optional) -- Local dimension of each particle/site. For a qubit, d=2.
- Returns
MPS -- List of N tensors containing the left-canonical MPS. The first and last tensors are of order 2 (matrices), while all the others are of order 3. We can see the sketch in the Example. The index ordering convention is from left-to-right. For instance, the "left" index of U2 is the first, the "bottom" one is the second, and the "right" one is the third.
- Return type
list of N tensors
Examples
>>> U1 - U2 - U3 - ... - UN >>> | | | |