quimb.tensor.tensor_3d#

Classes and algorithms related to 3D tensor networks.

Module Contents#

Classes#

Rotator3D

Object for rotating coordinates and various contraction functions so

TensorNetwork3D

Mixin class for tensor networks with a cubic lattice three-dimensional

TensorNetwork3DVector

Mixin class for a 3D square lattice vector TN, i.e. one with a single

TensorNetwork3DFlat

Mixin class for a 3D square lattice tensor network with a single tensor

PEPS3D

Projected Entangled Pair States object (3D).

Functions#

gen_3d_bonds(Lx, Ly, Lz[, steppers, coo_filter, cyclic])

Convenience function for tiling pairs of bond coordinates on a 3D

gen_3d_plaquette(coo0, steps)

Generate a plaquette at site coo0 by stepping first in steps and

gen_3d_plaquettes(Lx, Ly, Lz[, tiling])

Generate a tiling of plaquettes in a cubic 3D lattice.

gen_3d_strings(Lx, Ly, Lz)

Generate all length-wise strings in a cubic 3D lattice.

parse_boundary_sequence(sequence)

is_lone_coo(where)

Check if where has been specified as a single coordinate triplet.

calc_cell_sizes(coo_groups[, autogroup])

cell_to_sites(p)

Turn a cell ((i0, j0, k0), (di, dj, dk)) into the sites it contains.

sites_to_cell(sites)

Get the minimum covering cell for sites.

calc_cell_map(cells)

Attributes#

quimb.tensor.tensor_3d.gen_3d_bonds(Lx, Ly, Lz, steppers=None, coo_filter=None, cyclic=False)[source]#

Convenience function for tiling pairs of bond coordinates on a 3D lattice given a function like lambda i, j, k: (i + 1, j + 1, k + 1).

Parameters:
  • Lx (int) – The number of x-slices.

  • Ly (int) – The number of y-slices.

  • Lz (int) – The number of z-slices.

  • steppers (callable or sequence of callable) – Function(s) that take args (i, j, k) and generate another coordinate, thus defining a bond. Only valid steps are taken. If not given, defaults to nearest neighbor bonds.

  • coo_filter (callable) – Function that takes args (i, j, k) and only returns True if this is to be a valid starting coordinate.

Yields:

bond (tuple[tuple[int, int, int], tuple[int, int, int]]) – A pair of coordinates.

Examples

Generate nearest neighbor bonds:

>>> for bond in gen_3d_bonds(2, 2, 2, [lambda i, j, k: (i + 1, j, k),
...                                    lambda i, j, k: (i, j + 1, k),
...                                    lambda i, j, k: (i, j, k + 1)]):
...     print(bond)
((0, 0, 0), (1, 0, 0))
((0, 0, 0), (0, 1, 0))
((0, 0, 0), (0, 0, 1))
((0, 0, 1), (1, 0, 1))
((0, 0, 1), (0, 1, 1))
((0, 1, 0), (1, 1, 0))
((0, 1, 0), (0, 1, 1))
((0, 1, 1), (1, 1, 1))
((1, 0, 0), (1, 1, 0))
((1, 0, 0), (1, 0, 1))
((1, 0, 1), (1, 1, 1))
((1, 1, 0), (1, 1, 1))
quimb.tensor.tensor_3d.gen_3d_plaquette(coo0, steps)[source]#

Generate a plaquette at site coo0 by stepping first in steps and then the reverse steps.

Parameters:
  • coo0 (tuple) – The coordinate of the first site in the plaquette.

  • steps (tuple) – The steps to take to generate the plaquette. Each element should be one of ('x+', 'x-', 'y+', 'y-', 'z+', 'z-').

Yields:

coo (tuple) – The coordinates of the sites in the plaquette, including the last site which will be the same as the first.

quimb.tensor.tensor_3d.gen_3d_plaquettes(Lx, Ly, Lz, tiling='1')[source]#

Generate a tiling of plaquettes in a cubic 3D lattice.

Parameters:
  • Lx (int) – The length of the lattice in the x direction.

  • Ly (int) – The length of the lattice in the y direction.

  • Lz (int) – The length of the lattice in the z direction.

  • tiling ({'1', '2', '4', 'full'}) –

    The tiling to use:

    • ’1’: plaquettes in a sparse checkerboard pattern, such that each edge

      is covered by a maximum of one plaquette.

    • ’2’: less sparse checkerboard pattern, such that each edge is

      covered by a maximum of two plaquettes.

    • ’4’ or ‘full’: dense tiling of plaquettes. All bulk edges will

      be covered four times.

Yields:

plaquette (tuple[tuple[int]]) – The coordinates of the sites in each plaquette, including the last site which will be the same as the first.

quimb.tensor.tensor_3d.gen_3d_strings(Lx, Ly, Lz)[source]#

Generate all length-wise strings in a cubic 3D lattice.

class quimb.tensor.tensor_3d.Rotator3D(tn, xrange, yrange, zrange, from_which)[source]#

Object for rotating coordinates and various contraction functions so that the core algorithms only have to written once, but nor does the actual TN have to be modified.

property sweep_other#
cyclic_x()#
cyclic_y()#
cyclic_z()#
get_jnext(j)[source]#
get_knext(k)[source]#
quimb.tensor.tensor_3d._canonize_plane_opts#
quimb.tensor.tensor_3d._compress_plane_opts#
quimb.tensor.tensor_3d.BOUNDARY_SEQUENCE_MAP#
quimb.tensor.tensor_3d.parse_boundary_sequence(sequence)[source]#
class quimb.tensor.tensor_3d.TensorNetwork3D(ts=(), *, virtual=False, check_collisions=True)[source]#

Bases: quimb.tensor.tensor_arbgeom.TensorNetworkGen

Mixin class for tensor networks with a cubic lattice three-dimensional structure.

property Lx#

The number of x-slices.

property Ly#

The number of y-slices.

property Lz#

The number of z-slices.

property nsites#

The total number of sites.

property x_tag_id#

The string specifier for tagging each x-slice of this 3D TN.

property x_tags#

A tuple of all of the Lx different x-slice tags.

property y_tag_id#

The string specifier for tagging each y-slice of this 3D TN.

property y_tags#

A tuple of all of the Ly different y-slice tags.

property z_tag_id#

The string specifier for tagging each z-slice of this 3D TN.

property z_tags#

A tuple of all of the Lz different z-slice tags.

_NDIMS = 3#
_EXTRA_PROPS = ('_site_tag_id', '_x_tag_id', '_y_tag_id', '_z_tag_id', '_Lx', '_Ly', '_Lz')#
flatten_[source]#
contract_boundary_from_[source]#
contract_boundary_[source]#
contract_ctmrg_[source]#
coarse_grain_hotrg_[source]#
contract_hotrg_[source]#
_compatible_3d(other)[source]#

Check whether self and other are compatible 3D tensor networks such that they can remain a 3D 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 TensorNetwork3D instance.

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

TensorNetwork3D or TensorNetwork

site_tag(i, j=None, k=None)[source]#

The name of the tag specifiying the tensor at site (i, j, k).

x_tag(i)[source]#
y_tag(j)[source]#
z_tag(k)[source]#
maybe_convert_coo(coo)[source]#

Check if coo is a tuple of three ints and convert to the corresponding site tag if so.

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

This is the function that lets coordinates such as (i, j, k) be used for many ‘tag’ based functions.

gen_site_coos()[source]#

Generate coordinates for all the sites in this 3D TN.

gen_bond_coos()[source]#

Generate pairs of coordinates for all the bonds in this 3D TN.

valid_coo(coo, xrange=None, yrange=None, zrange=None)[source]#

Check whether coo is in-bounds.

Parameters:
  • coo ((int, int, int), optional) – The coordinates to check.

  • xrange ((int, int), optional) – The range of allowed values for the x, y, and z coordinates.

  • yrange ((int, int), optional) – The range of allowed values for the x, y, and z coordinates.

  • zrange ((int, int), optional) – The range of allowed values for the x, y, and z coordinates.

Return type:

bool

get_ranges_present()[source]#

Return the range of site coordinates present in this TN.

Returns:

xrange, yrange, zrange – The minimum and maximum site coordinates present in each direction.

Return type:

tuple[tuple[int, int]]

Examples

>>> tn = qtn.TN3D_rand(4, 4, 4, 2)
>>> tn_sub = tn.select_local('I1,2,3', max_distance=1)
>>> tn_sub.get_ranges_present()
((0, 2), (1, 3), (2, 3))
is_cyclic_x(j=None, k=None, imin=None, imax=None)[source]#

Check if the x dimension is cyclic (periodic), specifically whether a bond exists between (imin, j, k) and (imax, j, k), with default values of imin = 0 and imax = Lx - 1, and j and k the center of the lattice. If imin and imax are adjacent then this is considered False, since there is no ‘extra’ connectivity.

is_cyclic_y(k=None, i=None, jmin=None, jmax=None)[source]#

Check if the y dimension is cyclic (periodic), specifically whether a bond exists between (i, jmin, k) and (i, jmax, k), with default values of jmin = 0 and jmax = Ly - 1, and i and k the center of the lattice. If jmin and jmax are adjacent then this is considered False, since there is no ‘extra’ connectivity.

is_cyclic_z(i=None, j=None, kmin=None, kmax=None)[source]#

Check if the z dimension is cyclic (periodic), specifically whether a bond exists between (i, j, kmin) and (i, j, kmax), with default values of kmin = 0 and kmax = Lz - 1, and i and j the center of the lattice. If kmin and kmax are adjacent then this is considered False, since there is no ‘extra’ connectivity.

__getitem__(key)[source]#

Key based tensor selection, checking for integer based shortcut.

_repr_info()[source]#

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

flatten(fuse_multibonds=True, inplace=False)[source]#

Contract all tensors corresponding to each site into one.

gen_pairs(xrange=None, yrange=None, zrange=None, xreverse=False, yreverse=False, zreverse=False, coordinate_order='xyz', xstep=None, ystep=None, zstep=None, stepping_order='xyz', step_only=None)[source]#

Helper function for generating pairs of cooordinates for all bonds within a certain range, optionally specifying an order.

Parameters:
  • xrange ((int, int), optional) – The range of allowed values for the x, y, and z coordinates.

  • yrange ((int, int), optional) – The range of allowed values for the x, y, and z coordinates.

  • zrange ((int, int), optional) – The range of allowed values for the x, y, and z coordinates.

  • xreverse (bool, optional) – Whether to reverse the order of the x, y, and z sweeps.

  • yreverse (bool, optional) – Whether to reverse the order of the x, y, and z sweeps.

  • zreverse (bool, optional) – Whether to reverse the order of the x, y, and z sweeps.

  • coordinate_order (str, optional) – The order in which to sweep the x, y, and z coordinates. Earlier dimensions will change slower. If the corresponding range has size 1 then that dimension doesn’t need to be specified.

  • xstep (int, optional) – When generating a bond, step in this direction to yield the neighboring coordinate. By default, these follow xreverse, yreverse, and zreverse respectively.

  • ystep (int, optional) – When generating a bond, step in this direction to yield the neighboring coordinate. By default, these follow xreverse, yreverse, and zreverse respectively.

  • zstep (int, optional) – When generating a bond, step in this direction to yield the neighboring coordinate. By default, these follow xreverse, yreverse, and zreverse respectively.

  • stepping_order (str, optional) – The order in which to step the x, y, and z coordinates to generate bonds. Does not need to include all dimensions.

  • step_only (int, optional) – Only perform the ith steps in stepping_order, used to interleave canonizing and compressing for example.

Yields:

coo_a, coo_b (((int, int, int), (int, int, int)))

canonize_plane(xrange, yrange, zrange, equalize_norms=False, canonize_opts=None, **gen_pair_opts)[source]#

Canonize every pair of tensors within a subrange, optionally specifying a order to visit those pairs in.

compress_plane(xrange, yrange, zrange, max_bond=None, cutoff=1e-10, equalize_norms=False, compress_opts=None, **gen_pair_opts)[source]#

Compress every pair of tensors within a subrange, optionally specifying a order to visit those pairs in.

_contract_boundary_core(xrange, yrange, zrange, from_which, max_bond, cutoff=1e-10, canonize=True, canonize_interleave=True, layer_tags=None, compress_late=True, equalize_norms=False, compress_opts=None, canonize_opts=None)[source]#
_contract_boundary_projector(xrange, yrange, zrange, from_which, max_bond=None, cutoff=1e-10, canonize=False, layer_tags=None, lazy=False, equalize_norms=False, optimize='auto-hq', compress_opts=None, canonize_opts=None)[source]#
contract_boundary_from(xrange, yrange, zrange, from_which, max_bond=None, *, cutoff=1e-10, mode='peps', equalize_norms=False, compress_opts=None, canonize_opts=None, inplace=False, **contract_boundary_opts)[source]#

Unified entrypoint for contracting any rectangular patch of tensors from any direction, with any boundary method.

_contract_interleaved_boundary_sequence(*, contract_boundary_opts=None, sequence=None, xmin=None, xmax=None, ymin=None, ymax=None, zmin=None, zmax=None, max_separation=1, max_unfinished=1, around=None, equalize_norms=False, canonize=False, canonize_opts=None, final_contract=True, final_contract_opts=None, optimize='auto-hq', progbar=False, inplace=False)[source]#

Unified handler for performing iterleaved contractions in a sequence of inwards boundary directions.

contract_boundary(max_bond=None, *, cutoff=1e-10, mode='peps', canonize=True, compress_opts=None, sequence=None, xmin=None, xmax=None, ymin=None, ymax=None, zmin=None, zmax=None, max_separation=1, max_unfinished=1, around=None, equalize_norms=False, final_contract=True, final_contract_opts=None, optimize='auto-hq', progbar=False, inplace=False, **contract_boundary_opts)[source]#
contract_ctmrg(max_bond=None, *, cutoff=1e-10, canonize=False, canonize_opts=None, lazy=False, mode='projector', compress_opts=None, sequence=('xmin', 'xmax', 'ymin', 'ymax', 'zmin', 'zmax'), xmin=None, xmax=None, ymin=None, ymax=None, zmin=None, zmax=None, max_separation=1, max_unfinished=1, around=None, equalize_norms=False, final_contract=True, final_contract_opts=None, progbar=False, inplace=False, **contract_boundary_opts)[source]#
_compute_plane_envs(xrange, yrange, zrange, from_which, envs=None, storage_factory=None, **contract_boundary_opts)[source]#

Compute all ‘plane’ environments for the cube given by xrange, yrange, zrange, with direction given by from_which.

_maybe_compute_cell_env(key, envs=None, storage_factory=None, boundary_order=None, **contract_boundary_opts)[source]#

Recursively compute the necessary 2D, 1D, and 0D environments.

coarse_grain_hotrg(direction, max_bond=None, cutoff=1e-10, lazy=False, equalize_norms=False, optimize='auto-hq', compress_opts=None, inplace=False)[source]#

Coarse grain this tensor network in direction using HOTRG. This inserts oblique projectors between tensor pairs and then optionally contracts them into new sites for form a lattice half the size.

Parameters:
  • direction ({'x', 'y', 'z'}) – The direction to coarse grain in.

  • max_bond (int, optional) – The maximum bond dimension of the projector pairs inserted.

  • cutoff (float, optional) – The cutoff for the singular values of the projector pairs.

  • lazy (bool, optional) – Whether to contract the coarse graining projectors or leave them in the tensor network lazily. Default is to contract them.

  • equalize_norms (bool, optional) – Whether to equalize the norms of the tensors in the coarse grained lattice.

  • optimize (str, optional) – The optimization method to use when contracting the coarse grained lattice, if lazy=False.

  • compress_opts (None or dict, optional) – Supplied to insert_compressor_between_regions().

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

Returns:

The coarse grained tensor network, with size halved in direction.

Return type:

TensorNetwork2D

See also

contract_hotrg, insert_compressor_between_regions

contract_hotrg(max_bond=None, cutoff=1e-10, canonize=False, canonize_opts=None, sequence=('x', 'y', 'z'), max_separation=1, max_unfinished=1, lazy=False, equalize_norms=False, final_contract=True, final_contract_opts=None, progbar=False, inplace=False, **coarse_grain_opts)[source]#

Contract this tensor network using the finite version of HOTRG. See https://arxiv.org/abs/1201.1144v4 and https://arxiv.org/abs/1905.02351 for the more optimal computaton of the projectors used here. The TN is contracted sequentially in sequence directions by inserting oblique projectors between tensor pairs, and then optionally contracting these new effective sites. The algorithm stops when only one direction has a length larger than 2, and thus exact contraction can be used.

Parameters:
  • max_bond (int, optional) – The maximum bond dimension of the projector pairs inserted.

  • cutoff (float, optional) – The cutoff for the singular values of the projector pairs.

  • sequence (tuple of str, optional) – The directions to contract in. Default is to contract in all directions.

  • max_separation (int, optional) – The maximum distance between sides (i.e. length - 1) of the tensor network before that direction is considered finished.

  • max_unfinished (int, optional) – The maximum number of directions that can be unfinished (i.e. are still longer than max_separation + 1) before the coarse graining terminates.

  • lazy (bool, optional) – Whether to contract the coarse graining projectors or leave them in the tensor network lazily. Default is to contract them.

  • equalize_norms (bool or float, optional) – Whether to equalize the norms of the tensors in the tensor network after each coarse graining step.

  • final_contract (bool, optional) – Whether to exactly contract the remaining tensor network after the coarse graining contractions.

  • final_contract_opts (None or dict, optional) – Options to pass to contract(), optimize defaults to 'auto-hq'.

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

  • coarse_grain_opts – Additional options to pass to coarse_grain_hotrg().

Returns:

The contracted tensor network, which will have no more than one directino of length > 2.

Return type:

TensorNetwork3D

See also

coarse_grain_hotrg, insert_compressor_between_regions

quimb.tensor.tensor_3d.is_lone_coo(where)[source]#

Check if where has been specified as a single coordinate triplet.

quimb.tensor.tensor_3d.calc_cell_sizes(coo_groups, autogroup=True)[source]#
quimb.tensor.tensor_3d.cell_to_sites(p)[source]#

Turn a cell ((i0, j0, k0), (di, dj, dk)) into the sites it contains.

Examples

>>> cell_to_sites([(3, 4), (2, 2)])
((3, 4), (3, 5), (4, 4), (4, 5))
quimb.tensor.tensor_3d.sites_to_cell(sites)[source]#

Get the minimum covering cell for sites.

Examples

>>> sites_to_cell([(1, 3, 3), (2, 2, 4)])
((1, 2, 3) , (2, 2, 2))
quimb.tensor.tensor_3d.calc_cell_map(cells)[source]#
class quimb.tensor.tensor_3d.TensorNetwork3DVector(ts=(), *, virtual=False, check_collisions=True)[source]#

Bases: TensorNetwork3D, quimb.tensor.tensor_arbgeom.TensorNetworkGenVector

Mixin class for a 3D square lattice vector TN, i.e. one with a single physical index per site.

_EXTRA_PROPS = ('_site_tag_id', '_x_tag_id', '_y_tag_id', '_z_tag_id', '_Lx', '_Ly', '_Lz', '_site_ind_id')#
gate_[source]#
site_ind(i, j=None, k=None)[source]#
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.

phys_dim(i=None, j=None, k=None)[source]#

Get the size of the physical indices / a specific physical index.

gate(G, where, contract=False, tags=None, info=None, inplace=False, **compress_opts)[source]#

Apply a gate G to sites where, preserving the outer site inds.

class quimb.tensor.tensor_3d.TensorNetwork3DFlat(ts=(), *, virtual=False, check_collisions=True)[source]#

Bases: TensorNetwork3D

Mixin class for a 3D square lattice tensor network with a single tensor per site, for example, both PEPS and PEPOs.

_EXTRA_PROPS = ('_site_tag_id', '_x_tag_id', '_y_tag_id', '_z_tag_id', '_Lx', '_Ly', '_Lz')#
bond(coo1, coo2)[source]#

Get the name of the index defining the bond between sites at coo1 and coo2.

bond_size(coo1, coo2)[source]#

Return the size of the bond between sites at coo1 and coo2.

class quimb.tensor.tensor_3d.PEPS3D(arrays, *, shape='urfdlbp', tags=None, site_ind_id='k{},{},{}', site_tag_id='I{},{},{}', x_tag_id='X{}', y_tag_id='Y{}', z_tag_id='Z{}', **tn_opts)[source]#

Bases: TensorNetwork3DVector, TensorNetwork3DFlat

Projected Entangled Pair States object (3D).

Parameters:
  • arrays (sequence of sequence of sequence of array) – The core tensor data arrays.

  • shape (str, optional) – Which order the dimensions of the arrays are stored in, the default 'urfdlbp' stands for (‘up’, ‘right’, ‘front’, down’, ‘left’, ‘behind’, ‘physical’) meaning (x+, y+, z+, x-, y-, z-, physical) respectively. Arrays on the edge of lattice are assumed to be missing the corresponding dimension.

  • tags (set[str], optional) – Extra global tags to add to the tensor network.

  • site_ind_id (str, optional) – String specifier for naming convention of site indices.

  • site_tag_id (str, optional) – String specifier for naming convention of site tags.

  • x_tag_id (str, optional) – String specifier for naming convention of x-slice tags.

  • y_tag_id (str, optional) – String specifier for naming convention of y-slice tags.

  • z_tag_id (str, optional) – String specifier for naming convention of z-slice tags.

_EXTRA_PROPS = ('_site_tag_id', '_x_tag_id', '_y_tag_id', '_z_tag_id', '_Lx', '_Ly', '_Lz', '_site_ind_id')#
permute_arrays(shape='urfdlbp')[source]#

Permute the indices of each tensor in this PEPS3D 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.

classmethod from_fill_fn(fill_fn, Lx, Ly, Lz, bond_dim, phys_dim=2, **peps3d_opts)[source]#

Create a 3D PEPS from a filling function with signature fill_fn(shape).

Parameters:
  • Lx (int) – The number of x-slices.

  • Ly (int) – The number of y-slices.

  • Lz (int) – The number of z-slices.

  • bond_dim (int) – The bond dimension.

  • physical (int, optional) – The physical index dimension.

  • peps_opts – Supplied to PEPS3D.

Returns:

psi

Return type:

PEPS3D

classmethod empty(Lx, Ly, Lz, bond_dim, phys_dim=2, like='numpy', **peps3d_opts)[source]#

Create an empty 3D PEPS.

Parameters:
  • Lx (int) – The number of x-slices.

  • Ly (int) – The number of y-slices.

  • Lz (int) – The number of z-slices.

  • bond_dim (int) – The bond dimension.

  • physical (int, optional) – The physical index dimension.

  • peps3d_opts – Supplied to PEPS3D.

Returns:

psi

Return type:

PEPS3D

classmethod ones(Lx, Ly, Lz, bond_dim, phys_dim=2, like='numpy', **peps3d_opts)[source]#

Create a 3D PEPS whose tensors are filled with ones.

Parameters:
  • Lx (int) – The number of x-slices.

  • Ly (int) – The number of y-slices.

  • Lz (int) – The number of z-slices.

  • bond_dim (int) – The bond dimension.

  • physical (int, optional) – The physical index dimension.

  • peps3d_opts – Supplied to PEPS3D.

Returns:

psi

Return type:

PEPS3D

classmethod rand(Lx, Ly, Lz, bond_dim, phys_dim=2, dtype='float64', seed=None, **peps3d_opts)[source]#

Create a random (un-normalized) 3D PEPS.

Parameters:
  • Lx (int) – The number of x-slices.

  • Ly (int) – The number of y-slices.

  • Lz (int) – The number of z-slices.

  • bond_dim (int) – The bond dimension.

  • physical (int, optional) – The physical index dimension.

  • dtype (dtype, optional) – The dtype to create the arrays with, default is real double.

  • seed (int, optional) – A random seed.

  • peps_opts – Supplied to PEPS3D.

Returns:

psi

Return type:

PEPS3D

partial_trace_cluster(keep, max_bond=None, *, cutoff=1e-10, max_distance=0, fillin=0, gauges=False, flatten=False, normalized=True, symmetrized='auto', get=None, **contract_boundary_opts)[source]#
partial_trace(keep, max_bond=None, *, cutoff=1e-10, canonize=True, flatten=False, normalized=True, symmetrized='auto', envs=None, storage_factory=None, boundary_order=None, contract_cell_optimize='auto-hq', contract_cell_method='boundary', contract_cell_opts=None, get=None, **contract_boundary_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

compute_local_expectation(terms, max_bond=None, *, cutoff=1e-10, canonize=True, flatten=False, normalized=True, symmetrized='auto', return_all=False, envs=None, storage_factory=None, progbar=False, **contract_boundary_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]