mdopt.mps.canonical module#

This module contains the CanonicalMPS class. Hereafter, saying a MPS is in a canonical form will mean one of the following.

  1. Right-canonical: all tensors are right isometries, i.e.:

    ---( )---     ---
        |   |       |
        |   |  ==   |
        |   |       |
    ---(*)---     ---
    
  2. Left-canonical: all tensors are left isometries, i.e.:

    ---( )---     ---
    |   |         |
    |   |    ==   |
    |   |         |
    ---(*)---     ---
    

3) Mixed-canonical: all but one tensors are left or right isometries. This exceptional tensor will be hereafter called the orthogonality centre.

Note, that in the diagrams, a tensor with a star inside means that it is complex-conjugated.

A Matrix Product State is thus stored as a list of three-dimensional tensors as shown in the following diagram:

       i    i+1
...---( )---( )---...
       |     |
       |     |

Essentially, this corresponds to storing each A[i] or B[i] as shown in fig.4c in reference [1].

Note, that we enumerate the bonds from the right side of the tensors. For example, the bond 0 is a bond to the right of tensor 0.

class mdopt.mps.canonical.CanonicalMPS(tensors: List[ndarray], orth_centre: int | None = None, tolerance: float = 1e-12, chi_max: int = 10000)#

Class for finite-size canonical matrix product states with open boundary conditions.

tensors#

The tensors of the MPS, one per each physical site. Each tensor has legs (virtual left, physical, virtual right), in short (vL, i, vR).

Type:

List[np.ndarray]

orth_centre#

Position of the orthogonality centre, does not support negative indexing. As a convention, this attribute is taken 0 for a right-canonical form, len(tensors) - 1 for a left-canonical form, None for a product state.

Type:

Optional[int]

tolerance#

Numerical tolerance to zero out the singular values in Singular Value Decomposition.

Type:

float

chi_max#

The maximum bond dimension to keep in Singular Value Decompositions.

Type:

int

bond_dimensions#

The list of all bond dimensions of the MPS.

Type:

List[int]

phys_dimensions#

The list of all physical dimensions of the MPS.

Type:

List[int]

num_sites#

Number of sites.

Type:

int

num_bonds#

Number of bonds. Note, that the “ghost” bonds are not included.

Type:

int

Raises:
  • ValueError – If the orthogonality centre position does not correspond to the number of sites.

  • ValueError – If any of the tensors does not have three dimensions.

property all_dimensions: List[Tuple[int, ...]]#

Returns the list of all dimensions of the MPS.

property bond_dimensions: List[int]#

Returns the list of all bond dimensions of the MPS.

check_orth_centre() int | List[int] | None#

Checks the current position of the orthogonality centre by checking each tensor for the isometry conditions. Note, this method does not update the current instance’s orth_centre attribute.

compress(chi_max: int = 10000, cut: float = 1e-17, renormalise: bool = False, strategy: str = 'svd', return_truncation_errors: bool = False) Tuple[CanonicalMPS, List[float | None]]#

Compresses the MPS, i.e., runs the compress_bond method for each bond.

Parameters:
  • chi_max (int) – The maximum bond dimension to keep.

  • cut (float) – Singular values smaller than this will be discarded.

  • renormalise (bool) – Whether to renormalise the singular value spectrum after the cut.

  • strategy (str) – Which strategy to use for decomposition at the bond. Available options: svd, qr and svd_advanced.

  • return_truncation_errors (bool) – Whether to return the list of truncation errors (for each bond).

Returns:

  • compressed_mps (CanonicalMPS) – The new compressed Matrix Product State.

  • truncation_errors (Optional[List[float]]) – The truncation errors. Returned only if return_truncation_errors is set to True.

Raises:

ValueError – If the compression strategy is not supported.

compress_bond(bond: int, chi_max: int = 10000, cut: float = 1e-17, renormalise: bool = False, strategy: str = 'svd', return_truncation_error: bool = False) Tuple[CanonicalMPS, float | None]#

Compresses the bond at a given site, i.e., reduces its bond dimension.

Parameters:
  • bond (int) – The index of the bond to compress.

  • chi_max (int) – The maximum bond dimension to keep.

  • cut (float) – Singular values smaller than this will be discarded.

  • renormalise (bool) – Whether to renormalise the singular value spectrum after the cut.

  • strategy (str) – Which strategy to use for decomposition at the bond. Available options: svd, qr and svd_advanced.

  • return_truncation_error (bool) – Whether to return the truncation error.

Returns:

  • compressed_mps (CanonicalMPS) – The new compressed Matrix Product State.

  • truncation_error (Optional[float]) – The truncation error. Returned only if return_truncation_error is set to True.

Raises:
  • ValueError – If the bond index is out of range.

  • ValueError – If the compression strategy is not supported.

Notes

1) The bonds are being enumerated from the right side of the tensors. For example, the bond 0 is a bond to the right of tensor 0. 2) The compression scheme svd_advanced follows the scheme outlined in Fig.2 of https://arxiv.org/abs/1708.08932. This strategy can give speed ups for large bond dimensions by doing two SVDs on moderate-size matrices instead of one SVD on a large matrix.

conjugate() CanonicalMPS#

Returns a complex-conjugated version of the current MPS.

copy() CanonicalMPS#

Returns a copy of the current MPS.

dense(flatten: bool = True, renormalise: bool = False, norm: None | float | Literal['fro', 'nuc'] = 2) ndarray#

Returns a dense representation of the MPS.

Warning: this method can cause memory overload for number of sites > ~20!

Parameters:
  • flatten (bool) – Whether to merge all the physical indices to form a vector or not.

  • renormalise (bool) – Whether to renormalise the resulting tensor.

  • norm (Any) – The order of the norm to use while renormalising, see numpy.linalg.norm.

density_mpo() List[ndarray]#

Returns the MPO representation (as a list of tensors) of the density matrix defined by the MPS in a canonical form. Each tensor in the MPO list has legs (vL, vR, pU, pD), where v stands for “virtual”, p – for “physical”, and L, R, U, D stand for “left”, “right”, “up”, “down”.

Notes

This operation is depicted in the following diagram:

       i     j
 a     |     |    c                 i     j
...---(*)---(*)---...           ab  |     |  cd
                        -->  ...---[ ]---[ ]---...
...---( )---( )---...               |     |
 b     |     |    d                 k     l
       k     l

In the cartoon, {i,j,k,l} and {a,b,c,d} are indices. Here, the ( )’s represent the MPS tensors, the [ ]’s —the MPO tensors. The MPS with the physical legs up is complex-conjugated element-wise. The empty line between the MPS and its complex-conjugated version stands in fact for the tensor (kronecker) product.

Warning, this object can be memory-intensive for large bond dimensions!

entanglement_entropy(tolerance: float = 1e-12) ndarray#

Returns the entanglement entropy for bipartitions at each of the bonds.

explicit(tolerance: float = 1e-12, renormalise: bool = True) ExplicitMPS#

Transforms a CanonicalMPS instance into a ExplicitMPS instance. Essentially, retrieves each Γ[i] and Λ[i] from A[i] or B[i]. See fig.4b in [1] for reference.

Parameters:
  • tolerance (float) – Numerical tolerance for the singular values.

  • renormalise (bool) – Whether to renormalise singular values during each SVD.

left_canonical(renormalise: bool = True) CanonicalMPS#

Returns the current MPS in the left-canonical form. See eq.19 in [1] for reference.

marginal(sites_to_marginalise: List[int], renormalise: bool = False) CanonicalMPS#

Computes a marginal over a subset of sites of an MPS. Returns a new CanonicalMPS instance.

Parameters:
  • sites_to_marginalise (List[int]) – The sites to marginalise over.

  • renormalise (bool) – Whether to renormalise the resulting MPS.

Notes

Marginalizes (traces out) selected sites from the MPS by first contracting them with the trace tensor and then absorbing them. The contraction is done in the following way: if a traced tensor has a right neighbor, it is absorbed to the right. (If it is the last tensor, it is merged with its left neighbor.)

In the absorption step the physical leg is traced out: (vL, i, vR) –> (vL, vR) and then the traced tensor is contracted with the neighbor so that new_tensor = tensordot(traced, neighbor, axes=([1],[0])) yields the expected dimensions.

An example of marginalising is shown in the following diagram:

---(0)---(1)---(2)---(3)---    ---(2)---(3)---
    |     |     |     |     ->     |     |
    |     |     |     |            |     |
   (t)   (t)

Here, the (t) (trace) tensor is a tensor consisting af all 1’s and normalised by the square root of the physical dimension.

mixed_canonical(orth_centre: int, renormalise: bool = True) CanonicalMPS#

Returns the current MPS in the mixed-canonical form with the orthogonality centre being located at orth_centre.

Parameters:

orth_centre_index (int) – An integer which can take values 0, 1, ..., num_sites-1. Denotes the position of the orthogonality centre – the only non-isometry in the new canonical MPS.

move_orth_centre(final_pos: int, return_singular_values: bool = False, renormalise: bool = True) CanonicalMPS | Tuple[CanonicalMPS, List[list]]#

Moves the orthogonality centre from its current position to final_pos.

Returns a new version of the current CanonicalMPS instance with the orthogonality centre moved from self.orth_centre to final_pos, returns also the singular value tensors from every covered bond as well.

Parameters:
  • final_pos (int) – Final position of the orthogonality centre.

  • return_singular_values (bool) – Whether to return the singular values obtained at each involved bond.

  • renormalise (bool) – Whether to renormalise singular values during each SVD.

Raises:
  • ValueError – If final_pos does not match the MPS length.

  • ValueError – If self.orth_centre is still None after the search.

move_orth_centre_to_border(renormalise: bool = True) Tuple[CanonicalMPS, str]#

Moves the orthogonality centre from its current position to the closest border.

Parameters:

renormalise (bool) – Whether to renormalise singular values during each SVD.

Notes

Returns a new version of the current CanonicalMPS instance with the orthogonality centre moved to the closest (from the current position) border.

norm() float#

Computes the norm of the current MPS, that is, the modulus squared of its inner product with itself.

one_site_expectation_value(site: int, operator: ndarray) float | complex128#

Computes an expectation value of an arbitrary one-site operator (not necessarily unitary) on the given site.

Parameters:
  • site (int) – The site where the operator is applied.

  • operator (np.ndarray) – The one-site operator

Notes

An example of a one-site expectation value is shown in the following diagram:

( )---( )---( )---( )
 |     |     |     |
 |    (o)    |     |
 |     |     |     |
( )---( )---( )---( )
one_site_tensor(site: int) ndarray#

Returs a particular MPS tensor located at the corresponding site.

Parameters:

site (int) – The site index of the tensor.

one_site_tensor_iter() Iterable#

Returns an iterator over the one-site tensors for every site.

property phys_dimensions: List[int]#

Returns the list of all physical dimensions of the MPS.

reverse() CanonicalMPS#

Returns a reversed version of the current MPS.

right_canonical(renormalise: bool = True) CanonicalMPS#

Returns the current MPS in the right-canonical form. See eq.19 in [1] for reference.

two_site_expectation_value(site: int, operator: ndarray) float | complex128#

Computes an expectation value of an arbitrary two-site operator (not necessarily unitary) on the given site and its next neighbour. The operator has legs (UL, DL, UR, DR), where L, R, U, D stand for “left”, “right”, “up”, “down” accordingly.

Parameters:
  • site (int) – The first site where the operator is applied, the second site to be site + 1.

  • operator (np.ndarray) – The two-site operator.

Notes

An example of a two-site expectation value is shown in the following diagram:

( )---( )---( )---( )
 |     |     |     |
 |    (operator)   |
 |     |     |     |
( )---( )---( )---( )
two_site_tensor_next(site: int) ndarray#

Computes a two-site tensor on a given site and the next one.

Parameters:

site (int) – The site index of the tensor.

two_site_tensor_next_iter() Iterable#

Returns an iterator over the two-site tensors for every site and its right neighbour.

two_site_tensor_prev(site: int) ndarray#

Computes a two-site tensor on a given site and the previous one.

Parameters:

site (int) – The site index of the tensor.

two_site_tensor_prev_iter() Iterable#

Returns an iterator over the two-site tensors for every site and its left neighbour.