quimb.tensor.tensor_1d

Classes and algorithms related to 1D tensor networks.

Attributes

Classes

TensorNetwork1D

Base class for tensor networks with a one-dimensional structure.

TensorNetwork1DVector

1D Tensor network which overall is like a vector with a single type of

TensorNetwork1DOperator

Base class for tensor networks with a one-dimensional structure.

TensorNetwork1DFlat

1D Tensor network which has a flat structure.

MatrixProductState

Initialise a matrix product state, with auto labelling and tagging.

MatrixProductOperator

Initialise a matrix product operator, with auto labelling and tagging.

Dense1D

Mimics other 1D tensor network structures, but really just keeps the

SuperOperator1D

A 1D tensor network super-operator class:

TNLinearOperator1D

A 1D tensor network linear operator like:

Functions

expec_TN_1D(*tns[, compress, eps])

Compute the expectation of several 1D TNs, using transfer matrix

maybe_factor_gate_into_tensor(G, phys_dim, nsites, where)

gate_TN_1D(tn, G, where[, contract, tags, ...])

Act with the gate G on sites where, maintaining the outer

superop_TN_1D(tn_super, tn_op[, upper_ind_id, ...])

Take a tensor network superoperator and act with it on a

parse_cur_orthog([cur_orthog, info])

convert_cur_orthog(fn)

set_default_compress_mode(opts[, cyclic])

Module Contents

quimb.tensor.tensor_1d.align_TN_1D[source]
quimb.tensor.tensor_1d.expec_TN_1D(*tns, compress=None, eps=1e-15)[source]

Compute the expectation of several 1D TNs, using transfer matrix compression if any are periodic.

Parameters:
  • tns (sequence of TensorNetwork1D) – The MPS and MPO to find expectation of. Should start and begin with an MPS e.g. (MPS, MPO, ...,  MPS).

  • compress ({None, False, True}, optional) – Whether to perform transfer matrix compression on cyclic systems. If set to None (the default), decide heuristically.

  • eps (float, optional) – The accuracy of the transfer matrix compression.

Returns:

x – The expectation value.

Return type:

float

quimb.tensor.tensor_1d.maybe_factor_gate_into_tensor(G, phys_dim, nsites, where)[source]
quimb.tensor.tensor_1d.gate_TN_1D(tn, G, where, contract=False, tags=None, propagate_tags='sites', info=None, inplace=False, cur_orthog=None, **compress_opts)[source]

Act with the gate G on sites where, maintaining the outer indices of the 1D tensor network:

contract=False       contract=True
    . .                    . .             <- where
o-o-o-o-o-o-o        o-o-o-GGG-o-o-o
| | | | | | |        | | | / \ | | |
    GGG
    | |


contract='split-gate'        contract='swap-split-gate'
      . .                          . .                      <- where
  o-o-o-o-o-o-o                o-o-o-o-o-o-o
  | | | | | | |                | | | | | | |
      G~G                          G~G
      | |                          \ /
                                    X
                                   / \

contract='swap+split'
        . .            <- where
  o-o-o-G=G-o-o-o
  | | | | | | | |

Note that the sites in where do not have to be contiguous. By default, site tags will be propagated to the gate tensors, identifying a ‘light cone’.

Parameters:
  • tn (TensorNetwork1DVector) – The 1D vector-like tensor network, for example, and MPS.

  • G (array) – A square array to act with on sites where. It should have twice the number of dimensions as the number of sites. The second half of these will be contracted with the MPS, and the first half indexed with the correct site_ind_id. Sites are read left to right from the shape. A two-dimensional array is permissible if each dimension factorizes correctly.

  • where (int or sequence of int) – Where the gate should act.

  • contract ({False, 'split-gate', 'swap-split-gate',) –

    ‘auto-split-gate’, True, ‘swap+split’}, optional Whether to contract the gate into the 1D tensor network. If,

    • False: leave the gate uncontracted, the default

    • ’split-gate’: like False, but split the gate if it is two-site.

    • ’swap-split-gate’: like ‘split-gate’, but decompose the gate as if a swap had first been applied

    • ’auto-split-gate’: automatically select between the above three options, based on the rank of the gate.

    • True: contract the gate into the tensor network, if the gate acts on more than one site, this will produce an ever larger tensor.

    • ’swap+split’: Swap sites until they are adjacent, then contract the gate and split the resulting tensor, then swap the sites back to their original position. In this way an MPS structure can be explicitly maintained at the cost of rising bond-dimension.

  • tags (str or sequence of str, optional) – Tag the new gate tensor with these tags.

  • propagate_tags ({'sites', 'register', False, True}, optional) –

    Add any tags from the sites to the new gate tensor (only matters if contract=False else tags are merged anyway):

    • If 'sites', then only propagate tags matching e.g. ‘I{}’ and ignore all others. I.e. just propagate the lightcone.

    • If 'register', then only propagate tags matching the sites of where this gate was actually applied. I.e. ignore the lightcone, just keep track of which ‘registers’ the gate was applied to.

    • If False, propagate nothing.

    • If True, propagate all tags.

  • inplace – Perform the gate in place.

  • bool – Perform the gate in place.

  • optional – Perform the gate in place.

  • compress_opts – Supplied to split() if contract='swap+split' or gate_with_auto_swap() if contract='swap+split'.

Return type:

TensorNetwork1DVector

Examples

>>> p = MPS_rand_state(3, 7)
>>> p.gate_(spin_operator('X'), where=1, tags=['GX'])
>>> p
<MatrixProductState(tensors=4, L=3, max_bond=7)>
>>> p.outer_inds()
('k0', 'k1', 'k2')
quimb.tensor.tensor_1d.superop_TN_1D(tn_super, tn_op, upper_ind_id='k{}', lower_ind_id='b{}', so_outer_upper_ind_id=None, so_inner_upper_ind_id=None, so_inner_lower_ind_id=None, so_outer_lower_ind_id=None)[source]

Take a tensor network superoperator and act with it on a tensor network operator, maintaining the original upper and lower indices of the operator:

outer_upper_ind_id                           upper_ind_id
   | | | ... |                               | | | ... |
   +----------+                              +----------+
   | tn_super +---+                          | tn_super +---+
   +----------+   |     upper_ind_id         +----------+   |
   | | | ... |    |      | | | ... |         | | | ... |    |
inner_upper_ind_id|     +-----------+       +-----------+   |
                  |  +  |   tn_op   |   =   |   tn_op   |   |
inner_lower_ind_id|     +-----------+       +-----------+   |
   | | | ... |    |      | | | ... |         | | | ... |    |
   +----------+   |      lower_ind_id        +----------+   |
   | tn_super +---+                          | tn_super +---+
   +----------+                              +----------+
   | | | ... | <--                           | | | ... |
outer_lower_ind_id                           lower_ind_id
Parameters:
  • tn_super (TensorNetwork) – The superoperator in the form of a 1D-like tensor network.

  • tn_op (TensorNetwork) – The operator to be acted on in the form of a 1D-like tensor network.

  • upper_ind_id (str, optional) – Current id of the upper operator indices, e.g. usually 'k{}'.

  • lower_ind_id (str, optional) – Current id of the lower operator indices, e.g. usually 'b{}'.

  • so_outer_upper_ind_id (str, optional) – Current id of the superoperator’s upper outer indices, these will be reindexed to form the new effective operators upper indices.

  • so_inner_upper_ind_id (str, optional) – Current id of the superoperator’s upper inner indices, these will be joined with those described by upper_ind_id.

  • so_inner_lower_ind_id (str, optional) – Current id of the superoperator’s lower inner indices, these will be joined with those described by lower_ind_id.

  • so_outer_lower_ind_id (str, optional) – Current id of the superoperator’s lower outer indices, these will be reindexed to form the new effective operators lower indices.

Returns:

KAK – The tensornetwork of the superoperator acting on the operator.

Return type:

TensorNetwork

quimb.tensor.tensor_1d.parse_cur_orthog(cur_orthog='calc', info=None)[source]
quimb.tensor.tensor_1d.convert_cur_orthog(fn)[source]
class quimb.tensor.tensor_1d.TensorNetwork1D(ts=(), *, virtual=False, check_collisions=True)[source]

Bases: quimb.tensor.tensor_arbgeom.TensorNetworkGen

Base class for tensor networks with a one-dimensional structure.

_NDIMS = 1
_EXTRA_PROPS = ('_site_tag_id', '_L')
_CONTRACT_STRUCTURED = True
_compatible_1d(other)[source]

Check whether self and other are compatible 2D tensor networks such that they can remain a 2D tensor network when combined.

combine(other, *, virtual=False, check_collisions=True)[source]

Combine this tensor network with another, returning a new tensor network. If the two are compatible, cast the resulting tensor network to a TensorNetwork1D instance.

Parameters:
  • other (TensorNetwork1D or TensorNetwork) – The other tensor network to combine with.

  • virtual (bool, optional) – Whether the new tensor network should copy all the incoming tensors (False, the default), or view them as virtual (True).

  • check_collisions (bool, optional) – Whether to check for index collisions between the two tensor networks before combining them. If True (the default), any inner indices that clash will be mangled.

Return type:

TensorNetwork1D or TensorNetwork

property L
The number of sites, i.e. length.
property nsites
The number of sites.
gen_site_coos()[source]

Generate the coordinates of all possible sites.

site_tag(i)[source]

The name of the tag specifiying the tensor at site i.

slice2sites(tag_slice)[source]

Take a slice object, and work out its implied start, stop and step, taking into account cyclic boundary conditions.

Examples

Normal slicing:

>>> p = MPS_rand_state(10, bond_dim=7)
>>> p.slice2sites(slice(5))
(0, 1, 2, 3, 4)
>>> p.slice2sites(slice(4, 8))
(4, 5, 6, 7)

Slicing from end backwards:

>>> p.slice2sites(slice(..., -3, -1))
(9, 8)

Slicing round the end:

>>> p.slice2sites(slice(7, 12))
(7, 8, 9, 0, 1)
>>> p.slice2sites(slice(-3, 2))
(7, 8, 9, 0, 1)

If the start point is > end point (before modulo n), then step needs to be negative to return anything.

maybe_convert_coo(x)[source]

Check if x is an integer and convert to the corresponding site tag if so.

contract_structured(tag_slice, structure_bsz=5, inplace=False, **opts)[source]

Perform a structured contraction, translating tag_slice from a slice or to a cumulative sequence of tags.

Parameters:
  • tag_slice (slice or ...) – The range of sites, or for all.

  • inplace (bool, optional) – Whether to perform the contraction inplace.

Returns:

The result of the contraction, still a TensorNetwork if the contraction was only partial.

Return type:

TensorNetwork, Tensor or scalar

See also

contract, contract_tags, contract_cumulative

_repr_info()[source]

General info to show in various reprs. Sublasses can add more relevant info to this dict.

class quimb.tensor.tensor_1d.TensorNetwork1DVector(ts=(), *, virtual=False, check_collisions=True)[source]

Bases: TensorNetwork1D, quimb.tensor.tensor_arbgeom.TensorNetworkGenVector

1D Tensor network which overall is like a vector with a single type of site ind.

_EXTRA_PROPS = ('_site_tag_id', '_site_ind_id', '_L')
reindex_sites(new_id, where=None, inplace=False)[source]

Update the physical site index labels to a new string specifier. Note that this doesn’t change the stored id string with the TN.

Parameters:
  • new_id (str) – A string with a format placeholder to accept an int, e.g. “ket{}”.

  • where (None or slice) – Which sites to update the index labels on. If None (default) all sites.

  • inplace (bool) – Whether to reindex in place.

reindex_sites_[source]
site_ind(i)[source]

Get the physical index name of site i.

gate(*args, inplace=False, **kwargs)[source]

Apply a gate to this vector tensor network at sites where. This is essentially a wrapper around gate_inds() apart from where can be specified as a list of sites, and tags can be optionally, intelligently propagated to the new gate tensor.

\[| \psi \rangle \rightarrow G_\mathrm{where} | \psi \rangle\]
Parameters:
  • G (array_ike) – The gate array to apply, should match or be factorable into the shape (*phys_dims, *phys_dims).

  • where (node or sequence[node]) – The sites to apply the gate to.

  • contract ({False, True, 'split', 'reduce-split', 'split-gate',) – ‘swap-split-gate’, ‘auto-split-gate’}, optional How to apply the gate, see gate_inds().

  • tags (str or sequence of str, optional) – Tags to add to the new gate tensor.

  • propagate_tags ({False, True, 'register', 'sites'}, optional) –

    Whether to propagate tags to the new gate tensor:

    - False: no tags are propagated
    - True: all tags are propagated
    - 'register': only site tags corresponding to ``where`` are
      added.
    - 'sites': all site tags on the current sites are propgated,
      resulting in a lightcone like tagging.
    

  • info (None or dict, optional) – Used to store extra optional information such as the singular values if not absorbed.

  • inplace (bool, optional) – Whether to perform the gate operation inplace on the tensor network or not.

  • compress_opts – Supplied to tensor_split() for any contract methods that involve splitting. Ignored otherwise.

Return type:

TensorNetworkGenVector

See also

TensorNetwork.gate_inds

gate_[source]
expec(*args, **kwargs)[source]
correlation(A, i, j, B=None, **expec_opts)[source]

Correlation of operator A between i and j.

Parameters:
  • A (array) – The operator to act with, can be multi site.

  • i (int or sequence of int) – The first site(s).

  • j (int or sequence of int) – The second site(s).

  • expec_opts – Supplied to expec_TN_1D().

Returns:

C – The correlation <A(i)> + <A(j)> - <A(ij)>.

Return type:

float

Examples

>>> ghz = (MPS_computational_state('0000') +
...        MPS_computational_state('1111')) / 2**0.5
>>> ghz.correlation(pauli('Z'), 0, 1)
1.0
>>> ghz.correlation(pauli('Z'), 0, 1, B=pauli('X'))
0.0
class quimb.tensor.tensor_1d.TensorNetwork1DOperator(ts=(), *, virtual=False, check_collisions=True)[source]

Bases: TensorNetwork1D, quimb.tensor.tensor_arbgeom.TensorNetworkGenOperator

Base class for tensor networks with a one-dimensional structure.

_EXTRA_PROPS = ('_site_tag_id', '_upper_ind_id', '_lower_ind_id', '_L')
reindex_lower_sites(new_id, where=None, inplace=False)[source]

Update the lower site index labels to a new string specifier.

Parameters:
  • new_id (str) – A string with a format placeholder to accept an int, e.g. "ket{}".

  • where (None or slice) – Which sites to update the index labels on. If None (default) all sites.

  • inplace (bool) – Whether to reindex in place.

reindex_lower_sites_[source]
reindex_upper_sites(new_id, where=None, inplace=False)[source]

Update the upper site index labels to a new string specifier.

Parameters:
  • new_id (str) – A string with a format placeholder to accept an int, e.g. “ket{}”.

  • where (None or slice) – Which sites to update the index labels on. If None (default) all sites.

  • inplace (bool) – Whether to reindex in place.

reindex_upper_sites_[source]
quimb.tensor.tensor_1d.set_default_compress_mode(opts, cyclic=False)[source]
class quimb.tensor.tensor_1d.TensorNetwork1DFlat(ts=(), *, virtual=False, check_collisions=True)[source]

Bases: TensorNetwork1D

1D Tensor network which has a flat structure.

_EXTRA_PROPS = ('_site_tag_id', '_L')
left_canonize_site(i, bra=None)[source]

Left canonize this TN’s ith site, inplace:

    i                i
   -o-o-            ->-s-
... | | ...  ==> ... | | ...
Parameters:
  • i (int) – Which site to canonize. The site at i + 1 also absorbs the non-isometric part of the decomposition of site i.

  • bra (None or matching TensorNetwork to self, optional) – If set, also update this TN’s data with the conjugate canonization.

right_canonize_site(i, bra=None)[source]

Right canonize this TN’s ith site, inplace:

      i                i
   -o-o-            -s-<-
... | | ...  ==> ... | | ...
Parameters:

i (int) –

Which site to canonize. The site at i - 1 also absorbs the

non-isometric part of the decomposition of site i.

braNone or matching TensorNetwork to self, optional

If set, also update this TN’s data with the conjugate canonization.

left_canonicalize(stop=None, start=None, normalize=False, bra=None, inplace=False)[source]

Left canonicalize all or a portion of this TN (i.e. sweep the orthogonality center to the right). If this is a MPS, this implies that:

              i              i
>->->->->->->-o-o-         +-o-o-
| | | | | | | | | ...  =>  | | | ...
>->->->->->->-o-o-         +-o-o-
Parameters:
  • start (int, optional) – If given, the site to start left canonicalizing at.

  • stop (int, optional) – If given, the site to stop left canonicalizing at.

  • normalize (bool, optional) – Whether to normalize the state, only works for OBC.

  • bra (MatrixProductState, optional) – If supplied, simultaneously left canonicalize this MPS too, assuming it to be the conjugate state.

  • inplace (bool, optional) – Whether to perform the operation inplace. If bra is supplied then it is always modifed inplace.

Return type:

TensorNetwork1DFlat

left_canonicalize_[source]
left_canonize[source]
right_canonicalize(stop=None, start=None, normalize=False, bra=None, inplace=False)[source]

Right canonicalize all or a portion of this TN (i.e. sweep the orthogonality center to the left). If this is a MPS,

i i

-o-o-<-<-<-<-<-<-< -o-o-+

… | | | | | | | | | -> … | | |

-o-o-<-<-<-<-<-<-< -o-o-+

Parameters:
  • start (int, optional) – If given, the site to start right canonizing at.

  • stop (int, optional) – If given, the site to stop right canonizing at.

  • normalize (bool, optional) – Whether to normalize the state.

  • bra (MatrixProductState, optional) – If supplied, simultaneously right canonicalize this MPS too, assuming it to be the conjugate state.

  • inplace (bool, optional) – Whether to perform the operation inplace. If bra is supplied then it is always modifed inplace.

Return type:

TensorNetwork1DFlat

right_canonicalize_[source]
right_canonize[source]
canonize_cyclic(i, bra=None, method='isvd', inv_tol=1e-10)[source]

Bring this MatrixProductState into (possibly only approximate) canonical form at site(s) i.

Parameters:
  • i (int or slice) – The site or range of sites to make canonical.

  • bra (MatrixProductState, optional) – Simultaneously canonize this state as well, assuming it to be the co-vector.

  • method ({'isvd', 'svds', ...}, optional) – How to perform the lateral compression.

  • inv_tol (float, optional) – Tolerance with which to invert the gauge.

shift_orthogonality_center(current, new, bra=None)[source]

Move the orthogonality center of this MPS.

Parameters:
  • current (int) – The current orthogonality center.

  • new (int) – The target orthogonality center.

  • bra (MatrixProductState, optional) – If supplied, simultaneously move the orthogonality center of this MPS too, assuming it to be the conjugate state.

canonicalize(where, cur_orthog='calc', info=None, bra=None, inplace=False)[source]

Gauge this MPS into mixed canonical form, implying:

              i                      i
>->->->->- ->-o-<- -<-<-<-<-<      +-o-+
| | | | |...| | |...| | | | |  ->  | | |
>->->->->- ->-o-<- -<-<-<-<-<      +-o-+

You can also supply a min/max of sites to orthogonalize around, and a current location of the orthogonality center for efficiency:

      current                              where
      .......                              .....
>->->-c-c-c-c-<-<-<-<-<-<      >->->->->->-w-<-<-<-<-<-<
| | | | | | | | | | | | |  ->  | | | | | | | | | | | | |
>->->-c-c-c-c-<-<-<-<-<-<      >->->->->->-w-<-<-<-<-<-<
   cmin     cmax                           i   j

This would only move cmin to i and cmax to j if necessary.

Parameters:
  • where (int or sequence of int) – Which site(s) to orthogonalize around. If a sequence of int then make sure that section from min(where) to max(where) is orthog.

  • info (dict, optional) – If supplied, will be used to infer and store various extra information. Currently, the key “cur_orthog” is used to store the current orthogonality center. Its input value can be "calc", a single site, or a pair of sites representing the min/max range, inclusive. It will be updated to the actual range after.

  • bra (MatrixProductState, optional) – If supplied, simultaneously mixed canonicalize this MPS too, assuming it to be the conjugate state.

  • inplace (bool, optional) – Whether to perform the operation inplace. If bra is supplied then it is always modifed inplace.

Returns:

The mixed canonical form MPS.

Return type:

MatrixProductState

canonicalize_[source]
canonize[source]
left_compress_site(i, bra=None, **compress_opts)[source]

Left compress this 1D TN’s ith site, such that the site is then left unitary with its right bond (possibly) reduced in dimension.

Parameters:
  • i (int) – Which site to compress.

  • bra (None or matching TensorNetwork to self, optional) – If set, also update this TN’s data with the conjugate compression.

  • compress_opts – Supplied to Tensor.split().

right_compress_site(i, bra=None, **compress_opts)[source]

Right compress this 1D TN’s ith site, such that the site is then right unitary with its left bond (possibly) reduced in dimension.

Parameters:
  • i (int) – Which site to compress.

  • bra (None or matching TensorNetwork to self, optional) – If set, update this TN’s data with the conjugate compression.

  • compress_opts – Supplied to Tensor.split().

left_compress(start=None, stop=None, bra=None, **compress_opts)[source]

Compress this 1D TN, from left to right, such that it becomes left-canonical (unless absorb != 'right').

Parameters:
  • start (int, optional) – Site to begin compressing on.

  • stop (int, optional) – Site to stop compressing at (won’t itself be an isometry).

  • bra (None or TensorNetwork like this one, optional) – If given, update this TN as well, assuming it to be the conjugate.

  • compress_opts – Supplied to Tensor.split().

right_compress(start=None, stop=None, bra=None, **compress_opts)[source]

Compress this 1D TN, from right to left, such that it becomes right-canonical (unless absorb != 'left').

Parameters:
  • start (int, optional) – Site to begin compressing on.

  • stop (int, optional) – Site to stop compressing at (won’t itself be an isometry).

  • bra (None or TensorNetwork like this one, optional) – If given, update this TN as well, assuming it to be the conjugate.

  • compress_opts – Supplied to Tensor.split().

compress(form=None, **compress_opts)[source]

Compress this 1D Tensor Network, possibly into canonical form.

Parameters:
  • form ({None, 'flat', 'left', 'right'} or int) – Output form of the TN. None left canonizes the state first for stability reasons, then right_compresses (default). 'flat' tries to distribute the singular values evenly – state will not be canonical. 'left' and 'right' put the state into left and right canonical form respectively with a prior opposite sweep, or an int will put the state into mixed canonical form at that site.

  • compress_opts – Supplied to Tensor.split().

compress_site(i, canonize=True, info=None, bra=None, **compress_opts)[source]

Compress the bonds adjacent to site i, by default first setting the orthogonality center to that site:

     i                     i
-o-o-o-o-o-    -->    ->->~o~<-<-
 | | | | |             | | | | |
Parameters:
  • i (int) – Which site to compress around

  • canonize (bool, optional) – Whether to first set the orthogonality center to site i.

  • info (dict, optional) – If supplied, will be used to infer and store various extra information. Currently, the key “cur_orthog” is used to store the current orthogonality center. Its input value can be "calc", a single site, or a pair of sites representing the min/max range, inclusive. It will be updated to the actual range after.

  • bra (MatrixProductState, optional) – The conjugate state to also apply the compression to.

  • compress_opts – Supplied to tensor_split().

bond(i, j)[source]

Get the name of the index defining the bond between sites i and j.

bond_size(i, j)[source]

Return the size of the bond between site i and j.

bond_sizes()[source]
amplitude(b)[source]

Compute the amplitude of configuration b.

Parameters:

b (sequence of int) – The configuration to compute the amplitude of.

Returns:

c_b

Return type:

scalar

singular_values(i, info=None, method='svd')[source]

Find the singular values associated with the ith bond:

....L....   i
o-o-o-o-o-l-o-o-o-o-o-o-o-o-o-o-o
| | | | |   | | | | | | | | | | |
       i-1  ..........R..........

Leaves the 1D TN in mixed canoncial form at bond i.

Parameters:
  • i (int) – Which bond, or equivalently, the number of sites in the left partition.

  • info (dict, optional) – If supplied, will be used to infer and store various extra information. Currently, the key “cur_orthog” is used to store the current orthogonality center. Its input value can be "calc", a single site, or a pair of sites representing the min/max range, inclusive. It will be updated to the actual range after.

Returns:

svals – The singular values.

Return type:

1d-array

expand_bond_dimension(new_bond_dim, rand_strength=0.0, bra=None, inplace=True)[source]

Expand the bond dimensions of this 1D tensor network to at least new_bond_dim.

Parameters:
  • new_bond_dim (int) – Minimum bond dimension to expand to.

  • inplace (bool, optional) – Whether to perform the expansion in place.

  • bra (MatrixProductState, optional) – Mirror the changes to bra inplace, treating it as the conjugate state.

  • rand_strength (float, optional) – If rand_strength > 0, fill the new tensor entries with gaussian noise of strength rand_strength.

Return type:

MatrixProductState

count_canonized()[source]
calc_current_orthog_center()[source]

Calculate the site(s) of the current orthogonality center.

Returns:

The min/max of sites around which the TN is currently orthogonal.

Return type:

(int, int)

as_cyclic(inplace=False)[source]

Convert this flat, 1D, TN into cyclic form by adding a dummy bond between the first and last sites.

show(max_width=None)[source]
class quimb.tensor.tensor_1d.MatrixProductState(arrays, *, sites=None, L=None, shape='lrp', tags=None, site_ind_id='k{}', site_tag_id='I{}', **tn_opts)[source]

Bases: TensorNetwork1DVector, TensorNetwork1DFlat

Initialise a matrix product state, with auto labelling and tagging.

Parameters:
  • arrays (sequence of arrays) – The tensor arrays to form into a MPS.

  • sites (sequence of int, optional) – Construct the MPO on these sites only. If not given, enumerate from zero. Should be monotonically increasing and match arrays.

  • L (int, optional) – The number of sites the MPO should be defined on. If not given, this is taken as the max sites value plus one (i.e.g the number of arrays if sites is not given).

  • shape (str, optional) – String specifying layout of the tensors. E.g. ‘lrp’ (the default) indicates the shape corresponds left-bond, right-bond, physical index. End tensors have either ‘l’ or ‘r’ dropped from the string.

  • tags (str or sequence of str, optional) – Global tags to attach to all tensors.

  • site_ind_id (str) – A string specifiying how to label the physical site indices. Should contain a '{}' placeholder. It is used to generate the actual indices like: map(site_ind_id.format, range(len(arrays))).

  • site_tag_id (str) – A string specifiying how to tag the tensors at each site. Should contain a '{}' placeholder. It is used to generate the actual tags like: map(site_tag_id.format, range(len(arrays))).

_EXTRA_PROPS = ('_site_tag_id', '_site_ind_id', 'cyclic', '_L')
arrays

Get the tuple of raw arrays containing all the tensor network data.

_L
_site_ind_id
_site_tag_id
cyclic
lrp_ord
tensors = []

Get the tuple of tensors in this tensor network.

tags
bonds
classmethod from_fill_fn(fill_fn, L, bond_dim, phys_dim=2, sites=None, cyclic=False, shape='lrp', site_ind_id='k{}', site_tag_id='I{}', tags=None)[source]

Create an MPS by supplying a ‘filling’ function to generate the data for each site.

Parameters:
  • fill_fn (callable) – A function with signature fill_fn(shape : tuple[int]) -> array_like.

  • L (int) – The number of sites.

  • bond_dim (int) – The bond dimension.

  • phys_dim (int or Sequence[int], optional) – The physical dimension(s) of each site, if a sequence it will be cycled over.

  • sites (None or sequence of int, optional) – Construct the MPS on these sites only. If not given, enumerate from zero.

  • cyclic (bool, optional) – Whether the MPS should be cyclic (periodic).

  • shape (str, optional) – What specific order to layout the indices in, should be a sequence of 'l', 'r', and 'p', corresponding to left, right, and physical indices respectively.

  • site_ind_id (str, optional) – How to label the physical site indices.

  • site_tag_id (str, optional) – How to tag the physical sites.

  • tags (str or sequence of str, optional) – Global tags to attach to all tensors.

Return type:

MatrixProductState

classmethod from_dense(psi, dims=2, tags=None, site_ind_id='k{}', site_tag_id='I{}', **split_opts)[source]

Create a MatrixProductState directly from a dense vector

Parameters:
  • psi (array_like) – The dense state to convert to MPS from.

  • dims (int or sequence of int) – Physical subsystem dimensions of each site. If a single int, all sites have this same dimension, by default, 2.

  • tags (str or sequence of str, optional) – Global tags to attach to all tensors.

  • site_ind_id (str, optional) – How to index the physical sites, see MatrixProductState.

  • site_tag_id (str, optional) – How to tag the physical sites, see MatrixProductState.

  • split_opts – Supplied to tensor_split() to in order to partition the dense vector into tensors. absorb='left' is set by default, to ensure the compression is canonical / optimal.

Return type:

MatrixProductState

Examples

>>> dims = [2, 2, 2, 2, 2, 2]
>>> psi = rand_ket(prod(dims))
>>> mps = MatrixProductState.from_dense(psi, dims)
>>> mps.show()
 2 4 8 4 2
o-o-o-o-o-o
| | | | | |
add_MPS(other, inplace=False, **kwargs)[source]

Add another MatrixProductState to this one.

add_MPS_[source]
permute_arrays(shape='lrp')[source]

Permute the indices of each tensor in this MPS to match shape. This doesn’t change how the overall object interacts with other tensor networks but may be useful for extracting the underlying arrays consistently. This is an inplace operation.

Parameters:

shape (str, optional) – A permutation of 'lrp' specifying the desired order of the left, right, and physical indices respectively.

normalize(bra=None, eps=1e-15, insert=None)[source]

Normalize this MPS, optional with co-vector bra. For periodic MPS this uses transfer matrix SVD approximation with precision eps in order to be efficient. Inplace.

Parameters:
  • bra (MatrixProductState, optional) – If given, normalize this MPS with the same factor.

  • eps (float, optional) – If cyclic, precision to approximation transfer matrix with. Default: 1e-14.

  • insert (int, optional) – Insert the corrective normalization on this site, random if not given.

Returns:

old_norm – The old norm self.H @ self.

Return type:

float

gate_split(G, where, inplace=False, **compress_opts)[source]

Apply a two-site gate and then split resulting tensor to retrieve a MPS form:

-o-o-A-B-o-o-
 | | | | | |            -o-o-GGG-o-o-           -o-o-X~Y-o-o-
 | | GGG | |     ==>     | | | | | |     ==>     | | | | | |
 | | | | | |                 i j                     i j
     i j

As might be found in TEBD.

Parameters:
  • G (array) – The gate, with shape (d**2, d**2) for physical dimension d.

  • where ((int, int)) – Indices of the sites to apply the gate to.

  • compress_opts – Supplied to tensor_split().

See also

gate, gate_with_auto_swap

gate_split_[source]
swap_sites_with_compress(i, j, info=None, inplace=False, **compress_opts)[source]

Swap sites i and j by contracting, then splitting with the physical indices swapped. If the sites are not adjacent, this will happen multiple times.

Parameters:
  • i (int) – The first site to swap.

  • j (int) – The second site to swap.

  • cur_orthog (int, sequence of int, or 'calc') – If known, the current orthogonality center.

  • info (dict, optional) – If supplied, will be used to infer and store various extra information. Currently, the key “cur_orthog” is used to store the current orthogonality center. Its input value can be "calc", a single site, or a pair of sites representing the min/max range, inclusive. It will be updated to the actual range after.

  • inplace (bond, optional) – Perform the swaps inplace.

  • compress_opts – Supplied to tensor_split().

swap_sites_with_compress_[source]
swap_site_to(i, f, info=None, inplace=False, **compress_opts)[source]

Swap site i to site f, compressing the bond after each swap:

      i       f
0 1 2 3 4 5 6 7 8 9      0 1 2 4 5 6 7 3 8 9
o-o-o-x-o-o-o-o-o-o      >->->->->->->-x-<-<
| | | | | | | | | |  ->  | | | | | | | | | |
Parameters:
  • i (int) – The site to move.

  • f (int) – The new location for site i.

  • info (dict, optional) – If supplied, will be used to infer and store various extra information. Currently, the key “cur_orthog” is used to store the current orthogonality center. Its input value can be "calc", a single site, or a pair of sites representing the min/max range, inclusive. It will be updated to the actual range after.

  • inplace (bond, optional) – Perform the swaps inplace.

  • compress_opts – Supplied to tensor_split().

swap_site_to_[source]
gate_with_auto_swap(G, where, info=None, swap_back=True, inplace=False, **compress_opts)[source]

Perform a two site gate on this MPS by, if necessary, swapping and compressing the sites until they are adjacent, using gate_split, then unswapping the sites back to their original position.

Parameters:
  • G (array) – The gate, with shape (d**2, d**2) for physical dimension d.

  • where ((int, int)) – Indices of the sites to apply the gate to.

  • info (dict, optional) – If supplied, will be used to infer and store various extra information. Currently, the key “cur_orthog” is used to store the current orthogonality center. Its input value can be "calc", a single site, or a pair of sites representing the min/max range, inclusive. It will be updated to the actual range after.

  • swap_back (bool, optional) – Whether to swap the sites back to their original position after applying the gate. If not, for sites i < j, the site j will remain swapped to i + 1, and sites between i + 1 and j will be shifted one place up.

  • inplace (bond, optional) – Perform the swaps inplace.

  • compress_opts – Supplied to tensor_split().

See also

gate, gate_split

gate_with_auto_swap_[source]
gate_with_submpo(submpo, where=None, method='direct', transpose=False, info=None, inplace=False, inplace_mpo=False, **compress_opts)[source]

Apply an MPO, which only acts on a subset of sites, to this MPS, compressing the MPS with the MPO only on the minimal set of sites covering where, keeping the MPS form:

    │   │ │
    A───A─A
    │   │ │         ->    │ │ │ │ │ │ │ │
                          >─>─O━O━O━O─<─<
│ │ │ │ │ │ │ │
o─o─o─o─o─o─o─o
Parameters:
  • submpo (MatrixProductOperator) – The MPO to apply.

  • where (sequence of int, optional) – The range of sites the MPO acts on, will be inferred from the support of the MPO if not given.

  • method ({'direct", 'dm', 'zipup', 'zipup-first', 'fit'}, optional) – The compression method to use.

  • transpose (bool, optional) – Whether to transpose the MPO before applying it. By default the lower inds of the MPO are contracted with the MPS, if transposed the upper inds are contracted.

  • info (dict, optional) – If supplied, will be used to infer and store various extra information. Currently, the key “cur_orthog” is used to store the current orthogonality center. Its input value can be "calc", a single site, or a pair of sites representing the min/max range, inclusive. It will be updated to the actual range after.

  • inplace (bool, optional) – Whether to perform the application and compression inplace.

  • compress_opts – Supplied to tensor_network_1d_compress().

Return type:

MatrixProductState

gate_with_submpo_[source]
gate_with_mpo(mpo, method='direct', transpose=False, inplace=False, inplace_mpo=False, **compress_opts)[source]

Gate this MPS with an MPO and compress the result with one of various methods back to MPS form:

│ │ │ │ │ │ │ │
A─A─A─A─A─A─A─A
│ │ │ │ │ │ │ │     ->    │ │ │ │ │ │ │ │
                          O━O━O━O━O━O━O━O
│ │ │ │ │ │ │ │
o─o─o─o─o─o─o─o
Parameters:
  • mpo (MatrixProductOperator) – The MPO to apply.

  • max_bond (int, optional) – A maximum bond dimension to keep when compressing.

  • cutoff (float, optional) – A singular value cutoff to use when compressing.

  • method ({'direct", 'dm', 'zipup', 'zipup-first', 'fit', ...}, optional) – The compression method to use.

  • transpose (bool, optional) – Whether to transpose the MPO before applying it. By default the lower inds of the MPO are contracted with the MPS, if transposed the upper inds are contracted.

  • inplace (bool, optional) – Whether to perform the compression inplace.

  • inplace_mpo (bool, optional) – Whether the modify the MPO in place, a minor performance gain.

  • compress_opts – Other options supplied to tensor_network_1d_compress().

Return type:

MatrixProductState

gate_with_mpo_[source]
gate_nonlocal(G, where, dims=None, method='direct', info=None, inplace=False, **compress_opts)[source]

Apply a potentially non-local gate to this MPS by first decomposing it into an MPO, then compressing the MPS with MPO only on the minimal set of sites covering where.

Parameters:
  • G (array_like) – The gate to apply.

  • where (sequence of int) – The sites to apply the gate to.

  • max_bond (int, optional) – A maximum bond dimension to keep when compressing.

  • cutoff (float, optional) – A singular value cutoff to use when compressing.

  • dims (sequence of int, optional) – The factorized dimensions of the gate G, which should match the physical dimensions of the sites it acts on. Calculated if not supplied. If a single int, all sites are assumed to have this same dimension.

  • method ({'direct", 'dm', 'zipup', 'zipup-first', 'fit', ...}, optional) – The compression method to use.

  • info (dict, optional) – If supplied, will be used to infer and store various extra information. Currently, the key “cur_orthog” is used to store the current orthogonality center. Its input value can be "calc", a single site, or a pair of sites representing the min/max range, inclusive. It will be updated to the actual range after.

  • inplace (bool, optional) – Whether to perform the compression inplace.

  • compress_opts – Supplied to tensor_network_1d_compress().

Return type:

MatrixProductState

gate_nonlocal_[source]
flip(inplace=False)[source]

Reverse the order of the sites in the MPS, such that site i is now at site L - i - 1.

magnetization(i, direction='Z', info=None)[source]

Compute the magnetization at site i.

schmidt_values(i, info=None, method='svd')[source]

Find the schmidt values associated with the bipartition of this MPS between sites on either site of i. In other words, i is the number of sites in the left hand partition:

....L....   i
o-o-o-o-o-S-o-o-o-o-o-o-o-o-o-o-o
| | | | |   | | | | | | | | | | |
       i-1  ..........R..........

The schmidt values, S, are the singular values associated with the (i - 1, i) bond, squared, provided the MPS is mixed canonized at one of those sites.

Parameters:
  • i (int) – The number of sites in the left partition.

  • info (dict, optional) – If given, will be used to infer and store various extra information. Currently the key “cur_orthog” is used to store the current orthogonality center.

Returns:

S – The schmidt values.

Return type:

1d-array

entropy(i, info=None, method='svd')[source]

The entropy of bipartition between the left block of i sites and the rest.

Parameters:
  • i (int) – The number of sites in the left partition.

  • info (dict, optional) – If given, will be used to infer and store various extra information. Currently the key “cur_orthog” is used to store the current orthogonality center.

Return type:

float

schmidt_gap(i, info=None, method='svd')[source]

The schmidt gap of bipartition between the left block of i sites and the rest.

Parameters:
  • i (int) – The number of sites in the left partition.

  • info (dict, optional) – If given, will be used to infer and store various extra information. Currently the key “cur_orthog” is used to store the current orthogonality center.

Return type:

float

partial_trace(keep, upper_ind_id='b{}', rescale_sites=True)[source]

Partially trace this matrix product state, producing a matrix product operator.

Parameters:
  • keep (sequence of int or slice) – Indicies of the sites to keep.

  • upper_ind_id (str, optional) – The ind id of the (new) ‘upper’ inds, i.e. the ‘bra’ inds.

  • rescale_sites (bool, optional) – If True (the default), then the kept sites will be rescaled to (0, 1, 2, ...) etc. rather than keeping their original site numbers.

Returns:

rho – The density operator in MPO form.

Return type:

MatrixProductOperator

ptr(keep, upper_ind_id='b{}', rescale_sites=True)[source]

Alias of partial_trace().

bipartite_schmidt_state(sz_a, get='ket', info=None)[source]

Compute the reduced state for a bipartition of an OBC MPS, in terms of the minimal left/right schmidt basis:

     A            B
 .........     ...........
 >->->->->--s--<-<-<-<-<-<    ->   +-s-+
 | | | | |     | | | | | |         |   |
k0 k1...                          kA   kB
Parameters:
  • sz_a (int) – The number of sites in subsystem A, must be 0 < sz_a < N.

  • get ({'ket', 'rho', 'ket-dense', 'rho-dense'}, optional) –

    Get the:

    • ’ket’: vector form as tensor.

    • ’rho’: density operator form, i.e. vector outer product

    • ’ket-dense’: like ‘ket’ but return qarray.

    • ’rho-dense’: like ‘rho’ but return qarray.

  • info (dict, optional) – If given, will be used to infer and store various extra information. Currently the key “cur_orthog” is used to store the current orthogonality center.

static _do_lateral_compress(mps, kb, section, leave_short, ul, ll, heps, hmethod, hmax_bond, verbosity, compressed, **compress_opts)[source]
static _do_vertical_decomp(mps, kb, section, sysa, sysb, compressed, ul, ur, ll, lr, vmethod, vmax_bond, veps, verbosity, **compress_opts)[source]
partial_trace_compress(sysa, sysb, eps=1e-08, method=('isvd', None), max_bond=(None, 1024), leave_short=True, renorm=True, lower_ind_id='b{}', verbosity=0, **compress_opts)[source]

Perform a compressed partial trace using singular value lateral then vertical decompositions of transfer matrix products:

        .....sysa......     ...sysb....
o-o-o-o-A-A-A-A-A-A-A-A-o-o-B-B-B-B-B-B-o-o-o-o-o-o-o-o-o
| | | | | | | | | | | | | | | | | | | | | | | | | | | | |

                          ==> form inner product

        ...............     ...........
o-o-o-o-A-A-A-A-A-A-A-A-o-o-B-B-B-B-B-B-o-o-o-o-o-o-o-o-o
| | | | | | | | | | | | | | | | | | | | | | | | | | | | |
o-o-o-o-A-A-A-A-A-A-A-A-o-o-B-B-B-B-B-B-o-o-o-o-o-o-o-o-o

                          ==> lateral SVD on each section

          .....sysa......     ...sysb....
          /\             /\   /\         /\
  ... ~~~E  A~~~~~~~~~~~A  E~E  B~~~~~~~B  E~~~ ...
          \/             \/   \/         \/

                          ==> vertical SVD and unfold on A & B

                  |                 |
          /-------A-------\   /-----B-----\
  ... ~~~E                 E~E             E~~~ ...
          \-------A-------/   \-----B-----/
                  |                 |

With various special cases including OBC or end spins included in subsytems.

Parameters:
  • sysa (sequence of int) – The sites, which should be contiguous, defining subsystem A.

  • sysb (sequence of int) – The sites, which should be contiguous, defining subsystem B.

  • eps (float or (float, float), optional) – Tolerance(s) to use when compressing the subsystem transfer matrices and vertically decomposing.

  • method (str or (str, str), optional) – Method(s) to use for laterally compressing the state then vertially compressing subsytems.

  • max_bond (int or (int, int), optional) – The maximum bond to keep for laterally compressing the state then vertially compressing subsytems.

  • leave_short (bool, optional) – If True (the default), don’t try to compress short sections.

  • renorm (bool, optional) – If True (the default), renomalize the state so that tr(rho)==1.

  • lower_ind_id (str, optional) – The index id to create for the new density matrix, the upper_ind_id is automatically taken as the current site_ind_id.

  • compress_opts (dict, optional) – If given, supplied to partial_trace_compress to govern how singular values are treated. See tensor_split.

  • verbosity ({0, 1}, optional) – How much information to print while performing the compressed partial trace.

Returns:

rho_ab – Density matrix tensor network with outer_inds = ('k0', 'k1', 'b0', 'b1') for example.

Return type:

TensorNetwork

logneg_subsys(sysa, sysb, compress_opts=None, approx_spectral_opts=None, verbosity=0, approx_thresh=2**12)[source]

Compute the logarithmic negativity between subsytem blocks, e.g.:

                   sysa         sysb
                 .........       .....
... -o-o-o-o-o-o-A-A-A-A-A-o-o-o-B-B-B-o-o-o-o-o-o-o- ...
     | | | | | | | | | | | | | | | | | | | | | | | |
Parameters:
  • sysa (sequence of int) – The sites, which should be contiguous, defining subsystem A.

  • sysb (sequence of int) – The sites, which should be contiguous, defining subsystem B.

  • eps (float, optional) – Tolerance to use when compressing the subsystem transfer matrices.

  • method (str or (str, str), optional) – Method(s) to use for laterally compressing the state then vertially compressing subsytems.

  • compress_opts (dict, optional) – If given, supplied to partial_trace_compress to govern how singular values are treated. See tensor_split.

  • approx_spectral_opts – Supplied to approx_spectral_function().

Returns:

ln – The logarithmic negativity.

Return type:

float

See also

MatrixProductState.partial_trace_compress, approx_spectral_function

measure(site, remove=False, outcome=None, renorm=True, info=None, get=None, inplace=False)[source]

Measure this MPS at site, including projecting the state. Optionally remove the site afterwards, yielding an MPS with one less site. In either case the orthogonality center of the returned MPS is min(site, new_L - 1).

Parameters:
  • site (int) – The site to measure.

  • remove (bool, optional) –

    Whether to remove the site completely after projecting the measurement. If True, sites greater than site will be retagged and reindex one down, and the MPS will have one less site. E.g:

    0-1-2-3-4-5-6
           / / /  - measure and remove site 3
    0-1-2-4-5-6
                  - reindex sites (4, 5, 6) to (3, 4, 5)
    0-1-2-3-4-5
    

  • outcome (None or int, optional) – Specify the desired outcome of the measurement. If None, it will be randomly sampled according to the local density matrix.

  • renorm (bool, optional) – Whether to renormalize the state post measurement.

  • info (dict, optional) – If given, will be used to infer and store various extra information. Currently the key “cur_orthog” is used to store the current orthogonality center.

  • get ({None, 'outcome'}, optional) – If 'outcome', simply return the outcome, and don’t perform any projection.

  • inplace (bool, optional) – Whether to perform the measurement in place or not.

Returns:

  • outcome (int) – The measurement outcome, drawn from range(phys_dim).

  • psi (MatrixProductState) – The measured state, if get != 'outcome'.

measure_[source]
sample_configuration(seed=None, info=None)[source]

Sample a configuration from this MPS.

Parameters:
  • seed (None, int, or np.random.Generator, optional) – A random seed or generator to use.

  • info (dict, optional) – If given, will be used to infer and store various extra information. Currently the key “cur_orthog” is used to store the current orthogonality center.

sample(C, seed=None, info=None)[source]

Generate C samples rom this MPS, along with their probabilities.

Parameters:
  • C (int) – The number of samples to generate.

  • seed (None, int, or np.random.Generator, optional) – A random seed or generator to use.

  • info (dict, optional) – If given, will be used to infer and store various extra information. Currently the key “cur_orthog” is used to store the current orthogonality center.

Yields:
  • config (sequence of int) – The sample configuration.

  • omega (float) – The probability of this configuration.

class quimb.tensor.tensor_1d.MatrixProductOperator(arrays, *, sites=None, L=None, shape='lrud', tags=None, upper_ind_id='k{}', lower_ind_id='b{}', site_tag_id='I{}', **tn_opts)[source]

Bases: TensorNetwork1DOperator, TensorNetwork1DFlat

Initialise a matrix product operator, with auto labelling and tagging.

Parameters:
  • arrays (sequence of arrays) – The tensor arrays to form into a MPO.

  • sites (sequence of int, optional) – Construct the MPO on these sites only. If not given, enumerate from zero. Should be monotonically increasing and match arrays.

  • L (int, optional) – The number of sites the MPO should be defined on. If not given, this is taken as the max sites value plus one (i.e.g the number of arrays if sites is not given).

  • shape (str, optional) – String specifying layout of the tensors. E.g. ‘lrud’ (the default) indicates the shape corresponds left-bond, right-bond, ‘up’ physical index, ‘down’ physical index. End tensors have either ‘l’ or ‘r’ dropped from the string.

  • tags (str or sequence of str, optional) – Global tags to attach to all tensors.

  • upper_ind_id (str) – A string specifiying how to label the upper physical site indices. Should contain a '{}' placeholder. It is used to generate the actual indices like: map(upper_ind_id.format, range(len(arrays))).

  • lower_ind_id (str) – A string specifiying how to label the lower physical site indices. Should contain a '{}' placeholder. It is used to generate the actual indices like: map(lower_ind_id.format, range(len(arrays))).

  • site_tag_id (str) – A string specifiying how to tag the tensors at each site. Should contain a '{}' placeholder. It is used to generate the actual tags like: map(site_tag_id.format, range(len(arrays))).

_EXTRA_PROPS = ('_site_tag_id', '_upper_ind_id', '_lower_ind_id', 'cyclic', '_L')
arrays

Get the tuple of raw arrays containing all the tensor network data.

_L
_upper_ind_id
_lower_ind_id
_site_tag_id
cyclic
lrud_order
tensors = []

Get the tuple of tensors in this tensor network.

tags
bonds
classmethod from_fill_fn(fill_fn, L, bond_dim, phys_dim=2, sites=None, cyclic=False, shape='lrud', tags=None, upper_ind_id='k{}', lower_ind_id='b{}', site_tag_id='I{}')[source]

Create an MPO by supplying a ‘filling’ function to generate the data for each site.

Parameters:
  • fill_fn (callable) – A function with signature fill_fn(shape : tuple[int]) -> array_like.

  • L (int) – The number of sites.

  • bond_dim (int) – The bond dimension.

  • phys_dim (int or Sequence[int], optional) – The physical dimension(s) of each site, if a sequence it will be cycled over.

  • sites (None or sequence of int, optional) – Construct the MPO on these sites only. If not given, enumerate from zero.

  • cyclic (bool, optional) – Whether the MPO should be cyclic (periodic).

  • shape (str, optional) – String specifying layout of the tensors. E.g. ‘lrud’ (the default) indicates the shape corresponds left-bond, right-bond, ‘up’ physical index, ‘down’ physical index. End tensors have either ‘l’ or ‘r’ dropped from the string.

  • tags (str or sequence of str, optional) – Global tags to attach to all tensors.

  • upper_ind_id (str) – A string specifiying how to label the upper physical site indices. Should contain a '{}' placeholder.

  • lower_ind_id (str) – A string specifiying how to label the lower physical site indices. Should contain a '{}' placeholder.

  • site_tag_id (str, optional) – How to tag the physical sites. Should contain a '{}' placeholder.

Return type:

MatrixProductState

classmethod from_dense(A, dims=2, sites=None, L=None, tags=None, site_tag_id='I{}', upper_ind_id='k{}', lower_ind_id='b{}', **split_opts)[source]

Build an MPO from a raw dense matrix.

Parameters:
  • A (array) – The dense operator, it should be reshapeable to (*dims, *dims).

  • dims (int, sequence of int, optional) – The physical subdimensions of the operator. If any integer, assume all sites have the same dimension. If a sequence, the dimension of each site. Default is 2.

  • sites (sequence of int, optional) – The sites to place the operator on. If None, will place it on first len(dims) sites.

  • L (int, optional) – The total number of sites in the MPO, if the operator represents only a subset.

  • tags (str or sequence of str, optional) – Global tags to attach to all tensors.

  • site_tag_id (str, optional) – The string to use to label the site tags.

  • upper_ind_id (str, optional) – The string to use to label the upper physical indices.

  • lower_ind_id (str, optional) – The string to use to label the lower physical indices.

  • split_opts – Supplied to tensor_split().

Return type:

MatrixProductOperator

fill_empty_sites(mode='full', phys_dim=None, fill_array=None, inplace=False)[source]

Fill any empty sites of this MPO with identity tensors, adding size 1 bonds or draping existing bonds where necessary such that the resulting tensor has nearest neighbor bonds only.

Parameters:
  • mode ({'full', 'minimal'}, optional) – Whether to fill in all sites, including at either end, or simply the minimal range covering the min to max current sites present.

  • phys_dim (int, optional) – The physical dimension of the identity tensors to add. If not specified, will use the upper physical dimension of the first present site.

  • fill_array (array, optional) – The array to use for the identity tensors. If not specified, will use the identity array of the same dtype as the first present site.

  • inplace (bool, optional) – Whether to perform the operation inplace.

Returns:

The modified MPO.

Return type:

MatrixProductOperator

fill_empty_sites_[source]
add_MPO(other, inplace=False, **kwargs)[source]
add_MPO_[source]
_apply_mps(other, compress=False, contract=True, **compress_opts)[source]
_apply_mpo(other, compress=False, contract=True, **compress_opts)[source]
apply(other, compress=False, **compress_opts)[source]

Act with this MPO on another MPO or MPS, such that the resulting object has the same tensor network structure/indices as other.

For an MPS:

       | | | | | | | | | | | | | | | | | |
 self: A-A-A-A-A-A-A-A-A-A-A-A-A-A-A-A-A-A
       | | | | | | | | | | | | | | | | | |
other: x-x-x-x-x-x-x-x-x-x-x-x-x-x-x-x-x-x

                       -->

       | | | | | | | | | | | | | | | | | |   <- other.site_ind_id
  out: y=y=y=y=y=y=y=y=y=y=y=y=y=y=y=y=y=y

For an MPO:

       | | | | | | | | | | | | | | | | | |
 self: A-A-A-A-A-A-A-A-A-A-A-A-A-A-A-A-A-A
       | | | | | | | | | | | | | | | | | |
other: B-B-B-B-B-B-B-B-B-B-B-B-B-B-B-B-B-B
       | | | | | | | | | | | | | | | | | |

                       -->

       | | | | | | | | | | | | | | | | | |   <- other.upper_ind_id
  out: C=C=C=C=C=C=C=C=C=C=C=C=C=C=C=C=C=C
       | | | | | | | | | | | | | | | | | |   <- other.lower_ind_id

The resulting TN will have the same structure/indices as other, but probably with larger bonds (depending on compression).

Parameters:
Return type:

MatrixProductOperator or MatrixProductState

dot[source]
permute_arrays(shape='lrud')[source]

Permute the indices of each tensor in this MPO to match shape. This doesn’t change how the overall object interacts with other tensor networks but may be useful for extracting the underlying arrays consistently. This is an inplace operation.

Parameters:

shape (str, optional) – A permutation of 'lrud' specifying the desired order of the left, right, upper and lower (down) indices respectively.

trace(left_inds=None, right_inds=None)[source]

Take the trace of this MPO.

partial_transpose(sysa, inplace=False)[source]

Perform the partial transpose on this MPO by swapping the bra and ket indices on sites in sysa.

Parameters:
  • sysa (sequence of int or int) – The sites to transpose indices on.

  • inplace (bool, optional) – Whether to perform the partial transposition inplace.

Return type:

MatrixProductOperator

rand_state(bond_dim, **mps_opts)[source]

Get a random vector matching this MPO.

identity(**mpo_opts)[source]

Get a identity matching this MPO.

show(max_width=None)[source]
class quimb.tensor.tensor_1d.Dense1D(array, phys_dim=2, tags=None, site_ind_id='k{}', site_tag_id='I{}', **tn_opts)[source]

Bases: TensorNetwork1DVector

Mimics other 1D tensor network structures, but really just keeps the full state in a single tensor. This allows e.g. applying gates in the same way for quantum circuit simulation as lazily represented hilbert spaces.

Parameters:
  • array (array_like) – The full hilbert space vector - assumed to be made of equal hilbert spaces each of size phys_dim and will be reshaped as such.

  • phys_dim (int, optional) – The hilbert space size of each site, default: 2.

  • tags (sequence of str, optional) – Extra tags to add to the tensor network.

  • site_ind_id (str, optional) – String formatter describing how to label the site indices.

  • site_tag_id (str, optional) – String formatter describing how to label the site tags.

  • tn_opts – Supplied to TensorNetwork.

_EXTRA_PROPS = ('_site_ind_id', '_site_tag_id', '_L')
_L
dims
data
_site_ind_id
site_inds

Return a tuple of all site indices.

_site_tag_id
site_tags

All of the site tags.

T
classmethod rand(n, phys_dim=2, dtype=float, **dense1d_opts)[source]

Create a random dense vector ‘tensor network’.

class quimb.tensor.tensor_1d.SuperOperator1D(arrays, shape='lrkud', site_tag_id='I{}', outer_upper_ind_id='kn{}', inner_upper_ind_id='k{}', inner_lower_ind_id='b{}', outer_lower_ind_id='bn{}', tags=None, tags_upper=None, tags_lower=None, **tn_opts)[source]

Bases: TensorNetwork1D

A 1D tensor network super-operator class:

0   1   2       n-1
|   |   |        |     <-- outer_upper_ind_id
O===O===O==     =O
|\  |\  |\       |\     <-- inner_upper_ind_id
  )   )   ) ...    )   <-- K (size of local Kraus sum)
|/  |/  |/       |/     <-- inner_lower_ind_id
O===O===O==     =O
|   | : |        |     <-- outer_lower_ind_id
      :
     chi (size of entangling bond dim)
Parameters:

arrays (sequence of arrays) – The data arrays defining the superoperator, this should be a sequence of 2n arrays, such that the first two correspond to the upper and lower operators acting on site 0 etc. The arrays should be 5 dimensional unless OBC conditions are desired, in which case the first two and last two should be 4-dimensional. The dimensions of array can be should match the shape option.

_EXTRA_PROPS = ('_site_tag_id', '_outer_upper_ind_id', '_inner_upper_ind_id', '_inner_lower_ind_id',...
arrays

Get the tuple of raw arrays containing all the tensor network data.

_L
_outer_upper_ind_id
_inner_upper_ind_id
_inner_lower_ind_id
_outer_lower_ind_id
sites_present
outer_upper_inds
inner_upper_inds
inner_lower_inds
outer_lower_inds
_site_tag_id
tags
tags_upper
tags_lower
cyclic
classmethod rand(n, K, chi, phys_dim=2, herm=True, cyclic=False, dtype=complex, **superop_opts)[source]
property outer_upper_ind_id
property inner_upper_ind_id
property inner_lower_ind_id
property outer_lower_ind_id
class quimb.tensor.tensor_1d.TNLinearOperator1D(tn, left_inds, right_inds, start, stop, ldims=None, rdims=None, is_conj=False, is_trans=False)[source]

Bases: scipy.sparse.linalg.LinearOperator

A 1D tensor network linear operator like:

         start                 stop - 1
           .                     .
         :-O-O-O-O-O-O-O-O-O-O-O-O-:                 --+
         : | | | | | | | | | | | | :                   |
         :-H-H-H-H-H-H-H-H-H-H-H-H-:    acting on    --V
         : | | | | | | | | | | | | :                   |
         :-O-O-O-O-O-O-O-O-O-O-O-O-:                 --+
left_inds^                         ^right_inds

Like TNLinearOperator, but performs a structured contract from one end to the other than can handle very long chains possibly more efficiently by contracting in blocks from one end.

Parameters:
  • tn (TensorNetwork) – The tensor network to turn into a LinearOperator.

  • left_inds (sequence of str) – The left indicies.

  • right_inds (sequence of str) – The right indicies.

  • start (int) – Index of starting site.

  • stop (int) – Index of stopping site (does not include this site).

  • ldims (tuple of int, optional) – If known, the dimensions corresponding to left_inds.

  • rdims (tuple of int, optional) – If known, the dimensions corresponding to right_inds.

See also

TNLinearOperator

tn
tags
is_conj
is_trans
_conj_linop = None
_adjoint_linop = None
_transpose_linop = None
_matvec(vec)[source]

Default matrix-vector multiplication handler.

If self is a linear operator of shape (M, N), then this method will be called on a shape (N,) or (N, 1) ndarray, and should return a shape (M,) or (M, 1) ndarray.

This default implementation falls back on _matmat, so defining that will define matrix-vector multiplication as well.

_matmat(mat)[source]

Default matrix-matrix multiplication handler.

Falls back on the user-defined _matvec method, so defining that will define matrix multiplication (though in a very suboptimal way).

copy(conj=False, transpose=False)[source]
conj()[source]
_transpose()[source]

Default implementation of _transpose; defers to rmatvec + conj

_adjoint()[source]

Hermitian conjugate of this TNLO.

to_dense()[source]
toarray()[source]
property A