quimb.tensor.tensor_arbgeom

Classes and algorithms related to arbitrary geometry tensor networks.

Module Contents

Classes

TensorNetworkGen

A tensor network which notionally has a single tensor per 'site',

TensorNetworkGenVector

A tensor network which notionally has a single tensor and outer index

TensorNetworkGenOperator

A tensor network which notionally has a single tensor and two outer

Functions

get_coordinate_formatter(ndims)

tensor_network_align(*tns[, ind_ids, trace, inplace])

Align an arbitrary number of tensor networks in a stack-like geometry:

tensor_network_apply_op_vec(A, x[, which_A, contract, ...])

Apply a general a general tensor network representing an operator (has

tensor_network_apply_op_op(A, B[, which_A, which_B, ...])

Apply the operator (has upper and lower site inds) represented by tensor

create_lazy_edge_map(tn[, site_tags])

Given a tensor network, where each tensor is in exactly one group or

tensor_network_ag_sum(tna, tnb[, site_tags, negate, ...])

Add two tensor networks with arbitrary, but matching, geometries. They

gauge_product_boundary_vector(tn, tags[, which, ...])

_compute_expecs_maybe_in_parallel(fn, tn, terms[, ...])

Unified helper function for the various methods that compute many

_tn_local_expectation(tn, *args, **kwargs)

Define as function for pickleability.

_tn_local_expectation_cluster(tn, *args, **kwargs)

Define as function for pickleability.

_tn_local_expectation_exact(tn, *args, **kwargs)

Define as function for pickleability.

Attributes

quimb.tensor.tensor_arbgeom.get_coordinate_formatter(ndims)[source]
quimb.tensor.tensor_arbgeom.tensor_network_align(*tns, ind_ids=None, trace=False, inplace=False)[source]

Align an arbitrary number of tensor networks in a stack-like geometry:

a-a-a-a-a-a-a-a-a-a-a-a-a-a-a-a-a-a
| | | | | | | | | | | | | | | | | | <- ind_ids[0] (defaults to 1st id)
b-b-b-b-b-b-b-b-b-b-b-b-b-b-b-b-b-b
| | | | | | | | | | | | | | | | | | <- ind_ids[1]
               ...
| | | | | | | | | | | | | | | | | | <- ind_ids[-2]
y-y-y-y-y-y-y-y-y-y-y-y-y-y-y-y-y-y
| | | | | | | | | | | | | | | | | | <- ind_ids[-1]
z-z-z-z-z-z-z-z-z-z-z-z-z-z-z-z-z-z
Parameters:
  • tns (sequence of TensorNetwork) – The TNs to align, should be structured and either effective ‘vectors’ (have a site_ind_id) or ‘operators’ (have a up_ind_id and lower_ind_id).

  • ind_ids (None, or sequence of str) – String with format specifiers to id each level of sites with. Will be automatically generated like (tns[0].site_ind_id, "__ind_a{}__", "__ind_b{}__", ...) if not given.

  • inplace (bool) – Whether to modify the input tensor networks inplace.

Returns:

tns_aligned

Return type:

sequence of TensorNetwork

quimb.tensor.tensor_arbgeom.tensor_network_apply_op_vec(A, x, which_A='lower', contract=False, fuse_multibonds=True, compress=False, inplace=False, inplace_A=False, **compress_opts)[source]

Apply a general a general tensor network representing an operator (has upper_ind_id and lower_ind_id) to a tensor network representing a vector (has site_ind_id), by contracting each pair of tensors at each site then compressing the resulting tensor network. How the compression takes place is determined by the type of tensor network passed in. The returned tensor network has the same site indices as x, and it is the lower_ind_id of A that is contracted.

This is like performing A.to_dense() @ x.to_dense(), or the transpose thereof, depending on the value of which_A.

Parameters:
  • A (TensorNetworkGenOperator) – The tensor network representing the operator.

  • x (TensorNetworkGenVector) – The tensor network representing the vector.

  • which_A ({"lower", "upper"}, optional) – Whether to contract the lower or upper indices of A with the site indices of x.

  • contract (bool) – Whether to contract the tensors at each site after applying the operator, yielding a single tensor at each site.

  • fuse_multibonds (bool) – If contract=True, whether to fuse any multibonds after contracting the tensors at each site.

  • compress (bool) – Whether to compress the resulting tensor network.

  • inplace (bool) – Whether to modify x, the input vector tensor network inplace.

  • inplace_A (bool) – Whether to modify A, the operator tensor network inplace.

  • compress_opts – Options to pass to tn.compress, where tn is the resulting tensor network, if compress=True.

Returns:

The same type as x.

Return type:

TensorNetworkGenVector

quimb.tensor.tensor_arbgeom.tensor_network_apply_op_op(A, B, which_A='lower', which_B='upper', contract=False, fuse_multibonds=True, compress=False, inplace=False, inplace_A=False, **compress_opts)[source]

Apply the operator (has upper and lower site inds) represented by tensor network A to the operator represented by tensor network B. The resulting tensor network has the same upper and lower indices as B. Optionally contract the tensors at each site, fuse any multibonds, and compress the resulting tensor network.

This is like performing A.to_dense() @ B.to_dense(), or various combinations of tranposes thereof, depending on the values of which_A and which_B.

Parameters:
  • A (TensorNetworkGenOperator) – The tensor network representing the operator to apply.

  • B (TensorNetworkGenOperator) – The tensor network representing the target operator.

  • which_A ({"lower", "upper"}, optional) – Whether to contract the lower or upper indices of A.

  • which_B ({"lower", "upper"}, optional) – Whether to contract the lower or upper indices of B.

  • contract (bool) – Whether to contract the tensors at each site after applying the operator, yielding a single tensor at each site.

  • fuse_multibonds (bool) – If contract=True, whether to fuse any multibonds after contracting the tensors at each site.

  • compress (bool) – Whether to compress the resulting tensor network.

  • inplace (bool) – Whether to modify B, the target tensor network inplace.

  • inplace_A (bool) – Whether to modify A, the applied operator tensor network inplace.

  • compress_opts – Options to pass to tn.compress, where tn is the resulting tensor network, if compress=True.

Returns:

The same type as B.

Return type:

TensorNetworkGenOperator

quimb.tensor.tensor_arbgeom.create_lazy_edge_map(tn, site_tags=None)[source]

Given a tensor network, where each tensor is in exactly one group or ‘site’, compute which sites are connected to each other, without checking each pair.

Parameters:
  • tn (TensorNetwork) – The tensor network to analyze.

  • site_tags (None or sequence of str, optional) – Which tags to consider as ‘sites’, by default uses tn.site_tags.

Returns:

  • edges (dict[tuple[str, str], list[str]]) – Each key is a sorted pair of tags, which are connected, and the value is a list of the indices connecting them.

  • neighbors (dict[str, list[str]]) – For each site tag, the other site tags it is connected to.

quimb.tensor.tensor_arbgeom.tensor_network_ag_sum(tna, tnb, site_tags=None, negate=False, compress=False, inplace=False, **compress_opts)[source]

Add two tensor networks with arbitrary, but matching, geometries. They should have the same site tags, with a single tensor per site and sites connected by a single index only (but the name of this index can differ in the two TNs).

Parameters:
  • tna (TensorNetworkGen) – The first tensor network to add.

  • tnb (TensorNetworkGen) – The second tensor network to add.

  • site_tags (None or sequence of str, optional) – Which tags to consider as ‘sites’, by default uses tna.site_tags.

  • negate (bool, optional) – Whether to negate the second tensor network before adding.

  • compress (bool, optional) – Whether to compress the resulting tensor network, by calling the compress method with the given options.

  • inplace (bool, optional) – Whether to modify the first tensor network inplace.

Returns:

The resulting tensor network.

Return type:

TensorNetworkGen

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

Bases: quimb.tensor.tensor_core.TensorNetwork

A tensor network which notionally has a single tensor per ‘site’, though these could be labelled arbitrarily could also be linked in an arbitrary geometry by bonds.

property nsites

The total number of sites.

property sites

Tuple of the possible sites in this tensor network.

property site_tag_id

The string specifier for tagging each site of this tensor network.

property site_tags

All of the site tags.

property site_tags_present

All of the site tags still present in this tensor network.

_NDIMS = 1
_EXTRA_PROPS = ('_sites', '_site_tag_id')
retag_all_[source]
align_[source]
_compatible_arbgeom(other)[source]

Check whether self and other represent the same set of sites and are tagged equivalently.

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 TensorNetworkGen instance.

Parameters:
  • other (TensorNetworkGen 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:

TensorNetworkGen or TensorNetwork

gen_site_coos()[source]

Generate the coordinates of all sites, same as self.sites.

_get_site_set()[source]

The set of all sites.

gen_sites_present()[source]

Generate the sites which are currently present (e.g. if a local view of a larger tensor network), based on whether their tags are present.

Examples

>>> tn = qtn.TN3D_rand(4, 4, 4, 2)
>>> tn_sub = tn.select_local('I1,2,3', max_distance=1)
>>> list(tn_sub.gen_sites_present())
[(0, 2, 3), (1, 1, 3), (1, 2, 2), (1, 2, 3), (1, 3, 3), (2, 2, 3)]
site_tag(site)[source]

The name of the tag specifiying the tensor at site.

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

Modify the site tags for all or some tensors in this tensor network (without changing the site_tag_id).

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

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

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

retag_all(new_id, inplace=False)[source]

Retag all sites and change the site_tag_id.

_get_site_tag_set()[source]

The oset of all site tags.

filter_valid_site_tags(tags)[source]

Get the valid site tags from tags.

maybe_convert_coo(x)[source]

Check if x is a valid site and convert to the corresponding site tag if so, else return x.

gen_tags_from_coos(coos)[source]

Generate the site tags corresponding to the given coordinates.

_get_tids_from_tags(tags, which='all')[source]

This is the function that lets coordinates such as site be used for many ‘tag’ based functions.

reset_cached_properties()[source]

Reset any cached properties, one should call this when changing the actual geometry of a TN inplace, for example.

align(*args, inplace=False, **kwargs)[source]
__add__(other)[source]
__sub__(other)[source]
__iadd__(other)[source]
__isub__(other)[source]
quimb.tensor.tensor_arbgeom.gauge_product_boundary_vector(tn, tags, which='all', max_bond=1, smudge=1e-06, canonize_distance=0, select_local_distance=None, select_local_opts=None, **contract_around_opts)[source]
quimb.tensor.tensor_arbgeom._VALID_GATE_PROPAGATE
quimb.tensor.tensor_arbgeom._LAZY_GATE_CONTRACT
class quimb.tensor.tensor_arbgeom.TensorNetworkGenVector(ts=(), *, virtual=False, check_collisions=True)[source]

Bases: TensorNetworkGen

A tensor network which notionally has a single tensor and outer index per ‘site’, though these could be labelled arbitrarily and could also be linked in an arbitrary geometry by bonds.

property site_ind_id

The string specifier for the physical indices.

property site_inds

Return a tuple of all site indices.

property site_inds_present

All of the site inds still present in this tensor network.

_EXTRA_PROPS = ('_sites', '_site_tag_id', '_site_ind_id')
reindex_sites_[source]
reindex_all_[source]
to_qarray[source]
gate_with_op_lazy_[source]
gate_[source]
local_expectation_simple[source]
compute_local_expectation_simple[source]
compute_local_expectation_rehearse[source]
compute_local_expectation_tn[source]
site_ind(site)[source]
reset_cached_properties()[source]

Reset any cached properties, one should call this when changing the actual geometry of a TN inplace, for example.

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

Modify the site indices for all or some tensors in this vector tensor network (without changing the site_ind_id).

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

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

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

reindex_all(new_id, inplace=False)[source]

Reindex all physical sites and change the site_ind_id.

gen_inds_from_coos(coos)[source]

Generate the site inds corresponding to the given coordinates.

phys_dim(site=None)[source]

Get the physical dimension of site, defaulting to the first site if not specified.

to_dense(*inds_seq, to_qarray=False, to_ket=None, **contract_opts)[source]

Contract this tensor network ‘vector’ into a dense array. By default, turn into a ‘ket’ qarray, i.e. column vector of shape (d, 1).

Parameters:
  • inds_seq (sequence of sequences of str) – How to group the site indices into the dense array. By default, use a single group ordered like sites, but only containing those sites which are still present.

  • to_qarray (bool) – Whether to turn the dense array into a qarray, if the backend would otherwise be 'numpy'.

  • to_ket (None or str) – Whether to reshape the dense array into a ket (shape (d, 1) array). If None (default), do this only if the inds_seq is not supplied.

  • contract_opts – Options to pass to contract().

Return type:

array

gate_with_op_lazy(A, transpose=False, inplace=False)[source]

Act lazily with the operator tensor network A, which should have matching structure, on this vector/state tensor network, like A @ x. The returned tensor network will have the same structure as this one, but with the operator gated in lazily, i.e. uncontracted.

\[| x \rangle \rightarrow A | x \rangle\]

or (if transpose=True):

\[| x \rangle \rightarrow A^T | x \rangle\]
Parameters:
  • A (TensorNetworkGenOperator) – The operator tensor network to gate with, or apply to this tensor network.

  • transpose (bool, optional) – Whether to contract the lower or upper indices of A with the site indices of x. If False (the default), the lower indices of A will be contracted with the site indices of x, if True the upper indices of A will be contracted with the site indices of x, which is like applying A.T @ x.

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

Return type:

TensorNetworkGenVector

gate(G, where, contract=False, tags=None, propagate_tags=False, info=None, inplace=False, **compress_opts)[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_simple_(G, where, gauges, renorm=True, **gate_opts)[source]

Apply a gate to this vector tensor network at sites where, using simple update style gauging of the tensors first, as supplied in gauges. The new singular values for the bond are reinserted into gauges.

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

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

  • gauges (dict[str, array_like]) – The store of gauge bonds, the keys being indices and the values being the vectors. Only bonds present in this dictionary will be used.

  • renorm (bool, optional) – Whether to renormalise the singular after the gate is applied, before reinserting them into gauges.

gate_fit_local_(G, where, max_distance=0, fillin=0, gauges=None, **fit_opts)[source]
local_expectation_cluster(G, where, normalized=True, max_distance=0, fillin=False, gauges=None, optimize='auto', max_bond=None, rehearse=False, **contract_opts)[source]

Approximately compute a single local expectation value of the gate G at sites where, either treating the environment beyond max_distance as the identity, or using simple update style bond gauges as supplied in gauges.

This selects a local neighbourhood of tensors up to distance max_distance away from where, then traces over dangling bonds after potentially inserting the bond gauges, to form an approximate version of the reduced density matrix.

\[\langle \psi | G | \psi \rangle \approx \frac{ \mathrm{Tr} [ G \tilde{\rho}_\mathrm{where} ] }{ \mathrm{Tr} [ \tilde{\rho}_\mathrm{where} ] }\]

assuming normalized==True.

Parameters:
  • G (array_like) – The gate to compute the expecation of.

  • where (node or sequence[node]) – The sites to compute the expectation at.

  • normalized (bool, optional) – Whether to locally normalize the result, i.e. divide by the expectation value of the identity.

  • max_distance (int, optional) – The maximum graph distance to include tensors neighboring where when computing the expectation. The default 0 means only the tensors at sites where are used.

  • fillin (bool or int, optional) – When selecting the local tensors, whether and how many times to ‘fill-in’ corner tensors attached multiple times to the local region. On a lattice this fills in the corners. See select_local().

  • gauges (dict[str, array_like], optional) – The store of gauge bonds, the keys being indices and the values being the vectors. Only bonds present in this dictionary will be used.

  • optimize (str or PathOptimizer, optional) – The contraction path optimizer to use, when exactly contracting the local tensors.

  • max_bond (None or int, optional) – If specified, use compressed contraction.

  • rehearse ({False, 'tn', 'tree', True}, optional) –

    Whether to perform the computations or not:

    - False: perform the computation.
    - 'tn': return the tensor networks of each local expectation,
      without running the path optimizer.
    - 'tree': run the path optimizer and return the
      ``cotengra.ContractonTree`` for each local expectation.
    - True: run the path optimizer and return the ``PathInfo`` for
      each local expectation.
    

Returns:

expectation

Return type:

float

compute_local_expectation_cluster(terms, *, max_distance=0, fillin=False, normalized=True, gauges=None, optimize='auto', max_bond=None, return_all=False, rehearse=False, executor=None, progbar=False, **contract_opts)[source]

Compute all local expectations of the given terms, either treating the environment beyond max_distance as the identity, or using simple update style bond gauges as supplied in gauges.

This selects a local neighbourhood of tensors up to distance max_distance away from each term’s sites, then traces over dangling bonds after potentially inserting the bond gauges, to form an approximate version of the reduced density matrix.

\[\sum_\mathrm{i} \langle \psi | G_\mathrm{i} | \psi \rangle \approx \sum_\mathrm{i} \frac{ \mathrm{Tr} [ G_\mathrm{i} \tilde{\rho}_\mathrm{i} ] }{ \mathrm{Tr} [ \tilde{\rho}_\mathrm{i} ] }\]

assuming normalized==True.

Parameters:
  • terms (dict[node or (node, node), array_like]) – The terms to compute the expectation of, with keys being the sites and values being the local operators.

  • max_distance (int, optional) – The maximum graph distance to include tensors neighboring each term’s sites when computing the expectation. The default 0 means only the tensors at sites of each term are used.

  • fillin (bool or int, optional) – When selecting the local tensors, whether and how many times to ‘fill-in’ corner tensors attached multiple times to the local region. On a lattice this fills in the corners. See select_local().

  • normalized (bool, optional) – Whether to locally normalize the result, i.e. divide by the expectation value of the identity. This implies that a different normalization factor is used for each term.

  • gauges (dict[str, array_like], optional) – The store of gauge bonds, the keys being indices and the values being the vectors. Only bonds present in this dictionary will be used.

  • optimize (str or PathOptimizer, optional) – The contraction path optimizer to use, when exactly contracting the local tensors.

  • max_bond (None or int, optional) – If specified, use compressed contraction.

  • return_all (bool, optional) – Whether to return all results, or just the summed expectation.

  • rehearse ({False, 'tn', 'tree', True}, optional) –

    Whether to perform the computations or not:

    - False: perform the computation.
    - 'tn': return the tensor networks of each local expectation,
      without running the path optimizer.
    - 'tree': run the path optimizer and return the
      ``cotengra.ContractonTree`` for each local expectation.
    - True: run the path optimizer and return the ``PathInfo`` for
      each local expectation.
    

  • executor (Executor, optional) – If supplied compute the terms in parallel using this executor.

  • progbar (bool, optional) – Whether to show a progress bar.

  • contract_opts – Supplied to contract().

Returns:

expecs – If return_all==False, return the summed expectation value of the given terms. Otherwise, return a dictionary mapping each term’s location to the expectation value.

Return type:

float or dict[node or (node, node), float]

local_expectation_exact(G, where, optimize='auto-hq', normalized=True, rehearse=False, **contract_opts)[source]

Compute the local expectation of operator G at site(s) where by exactly contracting the full overlap tensor network.

compute_local_expectation_exact(terms, optimize='auto-hq', *, normalized=True, return_all=False, rehearse=False, executor=None, progbar=False, **contract_opts)[source]

Compute the local expectations of many operators, by exactly contracting the full overlap tensor network.

Parameters:
  • terms (dict[node or (node, node), array_like]) – The terms to compute the expectation of, with keys being the sites and values being the local operators.

  • optimize (str or PathOptimizer, optional) – The contraction path optimizer to use, when exactly contracting the full tensor network.

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

  • return_all (bool, optional) – Whether to return all results, or just the summed expectation.

  • rehearse ({False, 'tn', 'tree', True}, optional) –

    Whether to perform the computations or not:

    - False: perform the computation.
    - 'tn': return the tensor networks of each local expectation,
      without running the path optimizer.
    - 'tree': run the path optimizer and return the
      ``cotengra.ContractonTree`` for each local expectation.
    - True: run the path optimizer and return the ``PathInfo`` for
      each local expectation.
    

  • executor (Executor, optional) – If supplied compute the terms in parallel using this executor.

  • progbar (bool, optional) – Whether to show a progress bar.

  • contract_opts – Supplied to contract().

Returns:

expecs – If return_all==False, return the summed expectation value of the given terms. Otherwise, return a dictionary mapping each term’s location to the expectation value.

Return type:

float or dict[node or (node, node), float]

partial_trace(keep, max_bond, optimize, flatten=True, reduce=False, normalized=True, symmetrized='auto', rehearse=False, method='contract_compressed', **contract_compressed_opts)[source]

Partially trace this tensor network state, keeping only the sites in keep, using compressed contraction.

Parameters:
  • keep (iterable of hashable) – The sites to keep.

  • max_bond (int) – The maximum bond dimensions to use while compressed contracting.

  • optimize (str or PathOptimizer, optional) – The contraction path optimizer to use, should specifically generate contractions paths designed for compressed contraction.

  • flatten ({False, True, 'all'}, optional) – Whether to force ‘flattening’ (contracting all physical indices) of the tensor network before contraction, whilst this makes the TN generally more complex to contract, the accuracy is usually improved. If 'all' also flatten the tensors in keep.

  • reduce (bool, optional) – Whether to first ‘pull’ the physical indices off their respective tensors using QR reduction. Experimental.

  • normalized (bool, optional) – Whether to normalize the reduced density matrix at the end.

  • symmetrized ({'auto', True, False}, optional) – Whether to symmetrize the reduced density matrix at the end. This should be unecessary if flatten is set to True.

  • rehearse ({False, 'tn', 'tree', True}, optional) –

    Whether to perform the computation or not:

    - False: perform the computation.
    - 'tn': return the tensor network without running the path
      optimizer.
    - 'tree': run the path optimizer and return the
      ``cotengra.ContractonTree``..
    - True: run the path optimizer and return the ``PathInfo``.
    

  • contract_compressed_opts (dict, optional) – Additional keyword arguments to pass to contract_compressed().

Returns:

rho – The reduce density matrix of sites in keep.

Return type:

array_like

local_expectation(G, where, max_bond, optimize, flatten=True, normalized=True, symmetrized='auto', reduce=False, rehearse=False, **contract_compressed_opts)[source]

Compute the local expectation of operator G at site(s) where by approximately contracting the full overlap tensor network.

Parameters:
  • G (array_like) – The local operator to compute the expectation of.

  • where (node or sequence of nodes) – The sites to compute the expectation for.

  • max_bond (int) – The maximum bond dimensions to use while compressed contracting.

  • optimize (str or PathOptimizer, optional) – The contraction path optimizer to use, should specifically generate contractions paths designed for compressed contraction.

  • method ({'rho', 'rho-reduced'}, optional) – The method to use to compute the expectation value.

  • flatten (bool, optional) – Whether to force ‘flattening’ (contracting all physical indices) of the tensor network before contraction, whilst this makes the TN generally more complex to contract, the accuracy is usually much improved.

  • normalized (bool, optional) – If computing via partial_trace, whether to normalize the reduced density matrix at the end.

  • symmetrized ({'auto', True, False}, optional) – If computing via partial_trace, whether to symmetrize the reduced density matrix at the end. This should be unecessary if flatten is set to True.

  • rehearse ({False, 'tn', 'tree', True}, optional) –

    Whether to perform the computation or not:

    - False: perform the computation.
    - 'tn': return the tensor network without running the path
      optimizer.
    - 'tree': run the path optimizer and return the
      ``cotengra.ContractonTree``..
    - True: run the path optimizer and return the ``PathInfo``.
    

  • contract_compressed_opts (dict, optional) – Additional keyword arguments to pass to contract_compressed().

Returns:

expec

Return type:

float

compute_local_expectation(terms, max_bond, optimize, *, flatten=True, normalized=True, symmetrized='auto', reduce=False, return_all=False, rehearse=False, executor=None, progbar=False, **contract_compressed_opts)[source]

Compute the local expectations of many local operators, by approximately contracting the full overlap tensor network.

Parameters:
  • terms (dict[node or (node, node), array_like]) – The terms to compute the expectation of, with keys being the sites and values being the local operators.

  • max_bond (int) – The maximum bond dimension to use during contraction.

  • optimize (str or PathOptimizer) – The compressed contraction path optimizer to use.

  • method ({'rho', 'rho-reduced'}, optional) –

    The method to use to compute the expectation value.

    • ’rho’: compute the expectation value via the reduced density matrix.

    • ’rho-reduced’: compute the expectation value via the reduced density matrix, having reduced the physical indices onto the bonds first.

  • flatten (bool, optional) – Whether to force ‘flattening’ (contracting all physical indices) of the tensor network before contraction, whilst this makes the TN generally more complex to contract, the accuracy can often be much improved.

  • normalized (bool, optional) – Whether to locally normalize the result.

  • symmetrized ({'auto', True, False}, optional) – Whether to symmetrize the reduced density matrix at the end. This should be unecessary if flatten is set to True.

  • return_all (bool, optional) – Whether to return all results, or just the summed expectation. If rehease is not False, this is ignored and a dict is always returned.

  • rehearse ({False, 'tn', 'tree', True}, optional) –

    Whether to perform the computations or not:

    - False: perform the computation.
    - 'tn': return the tensor networks of each local expectation,
      without running the path optimizer.
    - 'tree': run the path optimizer and return the
      ``cotengra.ContractonTree`` for each local expectation.
    - True: run the path optimizer and return the ``PathInfo`` for
      each local expectation.
    

  • executor (Executor, optional) – If supplied compute the terms in parallel using this executor.

  • progbar (bool, optional) – Whether to show a progress bar.

  • contract_compressed_opts – Supplied to contract_compressed().

Returns:

expecs – If return_all==False, return the summed expectation value of the given terms. Otherwise, return a dictionary mapping each term’s location to the expectation value.

Return type:

float or dict[node or (node, node), float]

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

Bases: TensorNetworkGen

A tensor network which notionally has a single tensor and two outer indices per ‘site’, though these could be labelled arbitrarily and could also be linked in an arbitrary geometry by bonds. By convention, if converted to a dense matrix, the ‘upper’ indices would be on the left and the ‘lower’ indices on the right.

property upper_ind_id

The string specifier for the upper phyiscal indices.

property upper_inds

Return a tuple of all upper indices.

property upper_inds_present

Return a tuple of all upper indices still present in the tensor network.

property lower_ind_id

The string specifier for the lower phyiscal indices.

property lower_inds

Return a tuple of all lower indices.

property lower_inds_present

Return a tuple of all lower indices still present in the tensor network.

_EXTRA_PROPS = ('_sites', '_site_tag_id', '_upper_ind_id', '_lower_ind_id')
reindex_upper_sites_[source]
reindex_lower_sites_[source]
to_qarray[source]
gate_upper_with_op_lazy_[source]
gate_lower_with_op_lazy_[source]
gate_sandwich_with_op_lazy_[source]
upper_ind(site)[source]

Get the upper physical index name of site.

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

Modify the upper site indices for all or some tensors in this operator tensor network (without changing the upper_ind_id).

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

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

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

lower_ind(site)[source]

Get the lower physical index name of site.

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

Modify the lower site indices for all or some tensors in this operator tensor network (without changing the lower_ind_id).

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

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

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

to_dense(*inds_seq, to_qarray=False, **contract_opts)[source]

Contract this tensor network ‘operator’ into a dense array.

Parameters:
  • inds_seq (sequence of sequences of str) – How to group the site indices into the dense array. By default, use a single group ordered like sites, but only containing those sites which are still present.

  • to_qarray (bool) – Whether to turn the dense array into a qarray, if the backend would otherwise be 'numpy'.

  • contract_opts – Options to pass to contract().

Return type:

array

phys_dim(site=None, which='upper')[source]

Get the physical dimension of site.

gate_upper_with_op_lazy(A, transpose=False, inplace=False)[source]

Act lazily with the operator tensor network A, which should have matching structure, on this operator tensor network (B), like A @ B. The returned tensor network will have the same structure as this one, but with the operator gated in lazily, i.e. uncontracted.

\[B \rightarrow A B\]

or (if transpose=True):

\[B \rightarrow A^T B\]
Parameters:
  • A (TensorNetworkGenOperator) – The operator tensor network to gate with, or apply to this tensor network.

  • transpose (bool, optional) – Whether to contract the lower or upper indices of A with the upper indices of B. If False (the default), the lower indices of A will be contracted with the upper indices of B, if True the upper indices of A will be contracted with the upper indices of B, which is like applying the transpose first.

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

Return type:

TensorNetworkGenOperator

gate_lower_with_op_lazy(A, transpose=False, inplace=False)[source]

Act lazily ‘from the right’ with the operator tensor network A, which should have matching structure, on this operator tensor network (B), like B @ A. The returned tensor network will have the same structure as this one, but with the operator gated in lazily, i.e. uncontracted.

\[B \rightarrow B A\]

or (if transpose=True):

\[B \rightarrow B A^T\]
Parameters:
  • A (TensorNetworkGenOperator) – The operator tensor network to gate with, or apply to this tensor network.

  • transpose (bool, optional) – Whether to contract the upper or lower indices of A with the lower indices of this TN. If False (the default), the upper indices of A will be contracted with the lower indices of B, if True the lower indices of A will be contracted with the lower indices of this TN, which is like applying the transpose first.

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

Return type:

TensorNetworkGenOperator

gate_sandwich_with_op_lazy(A, inplace=False)[source]

Act lazily with the operator tensor network A, which should have matching structure, on this operator tensor network (B), like \(B \rightarrow A B A^\dagger\). The returned tensor network will have the same structure as this one, but with the operator gated in lazily, i.e. uncontracted.

Parameters:
  • A (TensorNetworkGenOperator) – The operator tensor network to gate with, or apply to this tensor network.

  • inplace (bool, optional) – Whether to perform the gate operation inplace on this tensor

Return type:

TensorNetworkGenOperator

quimb.tensor.tensor_arbgeom._compute_expecs_maybe_in_parallel(fn, tn, terms, return_all=False, executor=None, progbar=False, **kwargs)[source]

Unified helper function for the various methods that compute many expectations, possibly in parallel, possibly with a progress bar.

quimb.tensor.tensor_arbgeom._tn_local_expectation(tn, *args, **kwargs)[source]

Define as function for pickleability.

quimb.tensor.tensor_arbgeom._tn_local_expectation_cluster(tn, *args, **kwargs)[source]

Define as function for pickleability.

quimb.tensor.tensor_arbgeom._tn_local_expectation_exact(tn, *args, **kwargs)[source]

Define as function for pickleability.