quimb.tensor.tensor_3d ====================== .. py:module:: quimb.tensor.tensor_3d .. autoapi-nested-parse:: Classes and algorithms related to 3D tensor networks. Attributes ---------- .. autoapisummary:: quimb.tensor.tensor_3d._canonize_plane_opts quimb.tensor.tensor_3d._compress_plane_opts quimb.tensor.tensor_3d.BOUNDARY_SEQUENCE_MAP Classes ------- .. autoapisummary:: quimb.tensor.tensor_3d.Rotator3D quimb.tensor.tensor_3d.TensorNetwork3D quimb.tensor.tensor_3d.TensorNetwork3DVector quimb.tensor.tensor_3d.TensorNetwork3DFlat quimb.tensor.tensor_3d.PEPS3D Functions --------- .. autoapisummary:: quimb.tensor.tensor_3d.gen_3d_bonds quimb.tensor.tensor_3d.gen_3d_plaquette quimb.tensor.tensor_3d.gen_3d_plaquettes quimb.tensor.tensor_3d.gen_3d_strings quimb.tensor.tensor_3d.parse_boundary_sequence quimb.tensor.tensor_3d.is_lone_coo quimb.tensor.tensor_3d.calc_cell_sizes quimb.tensor.tensor_3d.cell_to_sites quimb.tensor.tensor_3d.sites_to_cell quimb.tensor.tensor_3d.calc_cell_map Module Contents --------------- .. py:function:: gen_3d_bonds(Lx, Ly, Lz, steppers=None, coo_filter=None, cyclic=False) 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)``. :param Lx: The number of x-slices. :type Lx: int :param Ly: The number of y-slices. :type Ly: int :param Lz: The number of z-slices. :type Lz: int :param steppers: 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. :type steppers: callable or sequence of callable :param coo_filter: Function that takes args ``(i, j, k)`` and only returns ``True`` if this is to be a valid starting coordinate. :type coo_filter: callable :Yields: **bond** (*tuple[tuple[int, int, int], tuple[int, int, int]]*) -- A pair of coordinates. .. rubric:: 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)) .. py:function:: gen_3d_plaquette(coo0, steps) Generate a plaquette at site ``coo0`` by stepping first in ``steps`` and then the reverse steps. :param coo0: The coordinate of the first site in the plaquette. :type coo0: tuple :param steps: The steps to take to generate the plaquette. Each element should be one of ``('x+', 'x-', 'y+', 'y-', 'z+', 'z-')``. :type steps: tuple :Yields: **coo** (*tuple*) -- The coordinates of the sites in the plaquette, including the last site which will be the same as the first. .. py:function:: gen_3d_plaquettes(Lx, Ly, Lz, tiling='1') Generate a tiling of plaquettes in a cubic 3D lattice. :param Lx: The length of the lattice in the x direction. :type Lx: int :param Ly: The length of the lattice in the y direction. :type Ly: int :param Lz: The length of the lattice in the z direction. :type Lz: int :param tiling: 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. :type tiling: {'1', '2', '4', 'full'} :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. .. py:function:: gen_3d_strings(Lx, Ly, Lz) Generate all length-wise strings in a cubic 3D lattice. .. py:class:: Rotator3D(tn, xrange, yrange, zrange, from_which) 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. .. py:attribute:: xrange .. py:attribute:: yrange .. py:attribute:: zrange .. py:attribute:: from_which .. py:attribute:: plane .. py:property:: sweep_other .. py:method:: cyclic_x() .. py:method:: cyclic_y() .. py:method:: cyclic_z() .. py:method:: get_jnext(j) .. py:method:: get_knext(k) .. py:data:: _canonize_plane_opts .. py:data:: _compress_plane_opts .. py:data:: BOUNDARY_SEQUENCE_MAP .. py:function:: parse_boundary_sequence(sequence) .. py:class:: TensorNetwork3D(ts=(), *, virtual=False, check_collisions=True) Bases: :py:obj:`quimb.tensor.tensor_arbgeom.TensorNetworkGen` Mixin class for tensor networks with a cubic lattice three-dimensional structure. .. py:attribute:: _NDIMS :value: 3 .. py:attribute:: _EXTRA_PROPS :value: ('_site_tag_id', '_x_tag_id', '_y_tag_id', '_z_tag_id', '_Lx', '_Ly', '_Lz') .. py:method:: _compatible_3d(other) Check whether ``self`` and ``other`` are compatible 3D tensor networks such that they can remain a 3D tensor network when combined. .. py:method:: combine(other, *, virtual=False, check_collisions=True) Combine this tensor network with another, returning a new tensor network. If the two are compatible, cast the resulting tensor network to a :class:`TensorNetwork3D` instance. :param other: The other tensor network to combine with. :type other: TensorNetwork3D or TensorNetwork :param virtual: Whether the new tensor network should copy all the incoming tensors (``False``, the default), or view them as virtual (``True``). :type virtual: bool, optional :param check_collisions: 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. :type check_collisions: bool, optional :rtype: TensorNetwork3D or TensorNetwork .. py:property:: Lx The number of x-slices. .. py:property:: Ly The number of y-slices. .. py:property:: Lz The number of z-slices. .. py:property:: nsites The total number of sites. .. py:method:: site_tag(i, j=None, k=None) The name of the tag specifiying the tensor at site ``(i, j, k)``. .. py:property:: x_tag_id The string specifier for tagging each x-slice of this 3D TN. .. py:method:: x_tag(i) .. py:property:: x_tags A tuple of all of the ``Lx`` different x-slice tags. .. py:property:: y_tag_id The string specifier for tagging each y-slice of this 3D TN. .. py:method:: y_tag(j) .. py:property:: y_tags A tuple of all of the ``Ly`` different y-slice tags. .. py:property:: z_tag_id The string specifier for tagging each z-slice of this 3D TN. .. py:method:: z_tag(k) .. py:property:: z_tags A tuple of all of the ``Lz`` different z-slice tags. .. py:method:: maybe_convert_coo(coo) Check if ``coo`` is a tuple of three ints and convert to the corresponding site tag if so. .. py:method:: _get_tids_from_tags(tags, which='all') This is the function that lets coordinates such as ``(i, j, k)`` be used for many 'tag' based functions. .. py:method:: gen_site_coos() Generate coordinates for all the sites in this 3D TN. .. py:method:: gen_bond_coos() Generate pairs of coordinates for all the bonds in this 3D TN. .. py:method:: valid_coo(coo, xrange=None, yrange=None, zrange=None) Check whether ``coo`` is in-bounds. :param coo: The coordinates to check. :type coo: (int, int, int), optional :param xrange: The range of allowed values for the x, y, and z coordinates. :type xrange: (int, int), optional :param yrange: The range of allowed values for the x, y, and z coordinates. :type yrange: (int, int), optional :param zrange: The range of allowed values for the x, y, and z coordinates. :type zrange: (int, int), optional :rtype: bool .. py:method:: get_ranges_present() Return the range of site coordinates present in this TN. :returns: **xrange, yrange, zrange** -- The minimum and maximum site coordinates present in each direction. :rtype: tuple[tuple[int, int]] .. rubric:: 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)) .. py:method:: is_cyclic_x(j=None, k=None, imin=None, imax=None) 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. .. py:method:: is_cyclic_y(k=None, i=None, jmin=None, jmax=None) 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. .. py:method:: is_cyclic_z(i=None, j=None, kmin=None, kmax=None) 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. .. py:method:: __getitem__(key) Key based tensor selection, checking for integer based shortcut. .. py:method:: _repr_info() General info to show in various reprs. Sublasses can add more relevant info to this dict. .. py:method:: flatten(fuse_multibonds=True, inplace=False) Contract all tensors corresponding to each site into one. .. py:attribute:: flatten_ .. py:method:: 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) Helper function for generating pairs of cooordinates for all bonds within a certain range, optionally specifying an order. :param xrange: The range of allowed values for the x, y, and z coordinates. :type xrange: (int, int), optional :param yrange: The range of allowed values for the x, y, and z coordinates. :type yrange: (int, int), optional :param zrange: The range of allowed values for the x, y, and z coordinates. :type zrange: (int, int), optional :param xreverse: Whether to reverse the order of the x, y, and z sweeps. :type xreverse: bool, optional :param yreverse: Whether to reverse the order of the x, y, and z sweeps. :type yreverse: bool, optional :param zreverse: Whether to reverse the order of the x, y, and z sweeps. :type zreverse: bool, optional :param coordinate_order: 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. :type coordinate_order: str, optional :param xstep: When generating a bond, step in this direction to yield the neighboring coordinate. By default, these follow ``xreverse``, ``yreverse``, and ``zreverse`` respectively. :type xstep: int, optional :param ystep: When generating a bond, step in this direction to yield the neighboring coordinate. By default, these follow ``xreverse``, ``yreverse``, and ``zreverse`` respectively. :type ystep: int, optional :param zstep: When generating a bond, step in this direction to yield the neighboring coordinate. By default, these follow ``xreverse``, ``yreverse``, and ``zreverse`` respectively. :type zstep: int, optional :param stepping_order: The order in which to step the x, y, and z coordinates to generate bonds. Does not need to include all dimensions. :type stepping_order: str, optional :param step_only: Only perform the ith steps in ``stepping_order``, used to interleave canonizing and compressing for example. :type step_only: int, optional :Yields: **coo_a, coo_b** (*((int, int, int), (int, int, int))*) .. py:method:: canonize_plane(xrange, yrange, zrange, equalize_norms=False, canonize_opts=None, **gen_pair_opts) Canonize every pair of tensors within a subrange, optionally specifying a order to visit those pairs in. .. py:method:: compress_plane(xrange, yrange, zrange, max_bond=None, cutoff=1e-10, equalize_norms=False, compress_opts=None, **gen_pair_opts) Compress every pair of tensors within a subrange, optionally specifying a order to visit those pairs in. .. py:method:: _contract_boundary_core_via_2d(xrange, yrange, zrange, from_which, max_bond, cutoff=1e-10, method='local-early', layer_tags=None, **compress_opts) .. py:method:: _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) .. py:method:: _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) .. py:method:: 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) Unified entrypoint for contracting any rectangular patch of tensors from any direction, with any boundary method. .. py:attribute:: contract_boundary_from_ .. py: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) Unified handler for performing iterleaved contractions in a sequence of inwards boundary directions. .. py:method:: 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) .. py:attribute:: contract_boundary_ .. py:method:: contract_ctmrg(max_bond=None, *, cutoff=1e-10, canonize=False, canonize_opts=None, lazy=False, mode='projector', 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, progbar=False, inplace=False, **contract_boundary_opts) .. py:attribute:: contract_ctmrg_ .. py:method:: _compute_plane_envs(xrange, yrange, zrange, from_which, envs=None, storage_factory=None, **contract_boundary_opts) Compute all 'plane' environments for the cube given by ``xrange``, ``yrange``, ``zrange``, with direction given by ``from_which``. .. py:method:: _maybe_compute_cell_env(key, envs=None, storage_factory=None, boundary_order=None, **contract_boundary_opts) Recursively compute the necessary 2D, 1D, and 0D environments. .. py:method:: coarse_grain_hotrg(direction, max_bond=None, cutoff=1e-10, lazy=False, equalize_norms=False, optimize='auto-hq', compress_opts=None, inplace=False) 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. :param direction: The direction to coarse grain in. :type direction: {'x', 'y', 'z'} :param max_bond: The maximum bond dimension of the projector pairs inserted. :type max_bond: int, optional :param cutoff: The cutoff for the singular values of the projector pairs. :type cutoff: float, optional :param lazy: Whether to contract the coarse graining projectors or leave them in the tensor network lazily. Default is to contract them. :type lazy: bool, optional :param equalize_norms: Whether to equalize the norms of the tensors in the coarse grained lattice. :type equalize_norms: bool, optional :param optimize: The optimization method to use when contracting the coarse grained lattice, if ``lazy=False``. :type optimize: str, optional :param compress_opts: Supplied to :meth:`~quimb.tensor.tensor_core.TensorNetwork.insert_compressor_between_regions`. :type compress_opts: None or dict, optional :param inplace: Whether to perform the coarse graining in place. :type inplace: bool, optional :returns: The coarse grained tensor network, with size halved in ``direction``. :rtype: TensorNetwork2D .. seealso:: :py:obj:`contract_hotrg`, :py:obj:`insert_compressor_between_regions` .. py:attribute:: coarse_grain_hotrg_ .. py:method:: 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) 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. :param max_bond: The maximum bond dimension of the projector pairs inserted. :type max_bond: int, optional :param cutoff: The cutoff for the singular values of the projector pairs. :type cutoff: float, optional :param sequence: The directions to contract in. Default is to contract in all directions. :type sequence: tuple of str, optional :param max_separation: The maximum distance between sides (i.e. length - 1) of the tensor network before that direction is considered finished. :type max_separation: int, optional :param max_unfinished: The maximum number of directions that can be unfinished (i.e. are still longer than max_separation + 1) before the coarse graining terminates. :type max_unfinished: int, optional :param lazy: Whether to contract the coarse graining projectors or leave them in the tensor network lazily. Default is to contract them. :type lazy: bool, optional :param equalize_norms: Whether to equalize the norms of the tensors in the tensor network after each coarse graining step. :type equalize_norms: bool or float, optional :param final_contract: Whether to exactly contract the remaining tensor network after the coarse graining contractions. :type final_contract: bool, optional :param final_contract_opts: Options to pass to :meth:`contract`, ``optimize`` defaults to ``'auto-hq'``. :type final_contract_opts: None or dict, optional :param inplace: Whether to perform the coarse graining in place. :type inplace: bool, optional :param coarse_grain_opts: Additional options to pass to :meth:`coarse_grain_hotrg`. :returns: The contracted tensor network, which will have no more than one directino of length > 2. :rtype: TensorNetwork3D .. seealso:: :py:obj:`coarse_grain_hotrg`, :py:obj:`insert_compressor_between_regions` .. py:attribute:: contract_hotrg_ .. py:function:: is_lone_coo(where) Check if ``where`` has been specified as a single coordinate triplet. .. py:function:: calc_cell_sizes(coo_groups, autogroup=True) .. py:function:: cell_to_sites(p) Turn a cell ``((i0, j0, k0), (di, dj, dk))`` into the sites it contains. .. rubric:: Examples >>> cell_to_sites([(3, 4), (2, 2)]) ((3, 4), (3, 5), (4, 4), (4, 5)) .. py:function:: sites_to_cell(sites) Get the minimum covering cell for ``sites``. .. rubric:: Examples >>> sites_to_cell([(1, 3, 3), (2, 2, 4)]) ((1, 2, 3) , (2, 2, 2)) .. py:function:: calc_cell_map(cells) .. py:class:: TensorNetwork3DVector(ts=(), *, virtual=False, check_collisions=True) Bases: :py:obj:`TensorNetwork3D`, :py:obj:`quimb.tensor.tensor_arbgeom.TensorNetworkGenVector` Mixin class for a 3D square lattice vector TN, i.e. one with a single physical index per site. .. py:attribute:: _EXTRA_PROPS :value: ('_site_tag_id', '_x_tag_id', '_y_tag_id', '_z_tag_id', '_Lx', '_Ly', '_Lz', '_site_ind_id') .. py:method:: site_ind(i, j=None, k=None) .. py:method:: reindex_sites(new_id, where=None, inplace=False) Modify the site indices for all or some tensors in this vector tensor network (without changing the ``site_ind_id``). :param new_id: A string with a format placeholder to accept a site, e.g. "ket{}". :type new_id: str :param where: Which sites to update the index labels on. If ``None`` (default) all sites. :type where: None or sequence :param inplace: Whether to reindex in place. :type inplace: bool .. py:method:: phys_dim(i=None, j=None, k=None) Get the size of the physical indices / a specific physical index. .. py:method:: gate(G, where, contract=False, tags=None, info=None, inplace=False, **compress_opts) Apply a gate ``G`` to sites ``where``, preserving the outer site inds. .. py:attribute:: gate_ .. py:class:: TensorNetwork3DFlat(ts=(), *, virtual=False, check_collisions=True) Bases: :py:obj:`TensorNetwork3D` Mixin class for a 3D square lattice tensor network with a single tensor per site, for example, both PEPS and PEPOs. .. py:attribute:: _EXTRA_PROPS :value: ('_site_tag_id', '_x_tag_id', '_y_tag_id', '_z_tag_id', '_Lx', '_Ly', '_Lz') .. py:class:: 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) Bases: :py:obj:`TensorNetwork3DVector`, :py:obj:`TensorNetwork3DFlat` Projected Entangled Pair States object (3D). :param arrays: The core tensor data arrays. :type arrays: sequence of sequence of sequence of array :param shape: 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. :type shape: str, optional :param tags: Extra global tags to add to the tensor network. :type tags: set[str], optional :param site_ind_id: String specifier for naming convention of site indices. :type site_ind_id: str, optional :param site_tag_id: String specifier for naming convention of site tags. :type site_tag_id: str, optional :param x_tag_id: String specifier for naming convention of x-slice tags. :type x_tag_id: str, optional :param y_tag_id: String specifier for naming convention of y-slice tags. :type y_tag_id: str, optional :param z_tag_id: String specifier for naming convention of z-slice tags. :type z_tag_id: str, optional .. py:attribute:: _EXTRA_PROPS :value: ('_site_tag_id', '_x_tag_id', '_y_tag_id', '_z_tag_id', '_Lx', '_Ly', '_Lz', '_site_ind_id') .. py:attribute:: _site_ind_id :value: 'k{},{},{}' .. py:attribute:: _site_tag_id :value: 'I{},{},{}' .. py:attribute:: _x_tag_id :value: 'X{}' .. py:attribute:: _y_tag_id :value: 'Y{}' .. py:attribute:: _z_tag_id :value: 'Z{}' .. py:attribute:: _Lx .. py:attribute:: _Ly .. py:attribute:: _Lz .. py:method:: permute_arrays(shape='urfdlbp') 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. :param shape: A permutation of ``'lrp'`` specifying the desired order of the left, right, and physical indices respectively. :type shape: str, optional .. py:method:: from_fill_fn(fill_fn, Lx, Ly, Lz, bond_dim, phys_dim=2, cyclic=False, shape='urfdlbp', **peps3d_opts) :classmethod: Create a 3D PEPS from a filling function with signature ``fill_fn(shape)``. :param Lx: The number of x-slices. :type Lx: int :param Ly: The number of y-slices. :type Ly: int :param Lz: The number of z-slices. :type Lz: int :param bond_dim: The bond dimension. :type bond_dim: int :param phys_dim: The physical index dimension. :type phys_dim: int, optional :param shape: How to layout the indices of the tensors, the default is ``(up, right, front, down, left, back, phys) == 'urfdlbp'``. This is the order of the shape supplied to the filling function. :type shape: str, optional :param peps_opts: Supplied to :class:`~quimb.tensor.tensor_3d.PEPS3D`. :returns: **psi** :rtype: PEPS3D .. py:method:: empty(Lx, Ly, Lz, bond_dim, phys_dim=2, like='numpy', **peps3d_opts) :classmethod: Create an empty 3D PEPS. :param Lx: The number of x-slices. :type Lx: int :param Ly: The number of y-slices. :type Ly: int :param Lz: The number of z-slices. :type Lz: int :param bond_dim: The bond dimension. :type bond_dim: int :param physical: The physical index dimension. :type physical: int, optional :param peps3d_opts: Supplied to :class:`~quimb.tensor.tensor_3d.PEPS3D`. :returns: **psi** :rtype: PEPS3D .. seealso:: :py:obj:`PEPS3D.from_fill_fn` .. py:method:: ones(Lx, Ly, Lz, bond_dim, phys_dim=2, like='numpy', **peps3d_opts) :classmethod: Create a 3D PEPS whose tensors are filled with ones. :param Lx: The number of x-slices. :type Lx: int :param Ly: The number of y-slices. :type Ly: int :param Lz: The number of z-slices. :type Lz: int :param bond_dim: The bond dimension. :type bond_dim: int :param physical: The physical index dimension. :type physical: int, optional :param peps3d_opts: Supplied to :class:`~quimb.tensor.tensor_3d.PEPS3D`. :returns: **psi** :rtype: PEPS3D .. seealso:: :py:obj:`PEPS3D.from_fill_fn` .. py:method:: rand(Lx, Ly, Lz, bond_dim, phys_dim=2, dist='normal', loc=0.0, dtype='float64', seed=None, **peps3d_opts) :classmethod: Create a random (un-normalized) 3D PEPS. :param Lx: The number of x-slices. :type Lx: int :param Ly: The number of y-slices. :type Ly: int :param Lz: The number of z-slices. :type Lz: int :param bond_dim: The bond dimension. :type bond_dim: int :param physical: The physical index dimension. :type physical: int, optional :param dtype: The dtype to create the arrays with, default is real double. :type dtype: dtype, optional :param seed: A random seed. :type seed: int, optional :param peps_opts: Supplied to :class:`~quimb.tensor.tensor_3d.PEPS3D`. :returns: **psi** :rtype: PEPS3D .. seealso:: :py:obj:`PEPS3D.from_fill_fn` .. py:method:: 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) Compute the approximate reduced density matrix at sites ``where`` by contracting a local cluster of tensors, potentially gauging the region with the simple update style bond gauges in ``gauges``. :param where: The sites to keep. :type where: sequence[node] :param gauges: The store of gauge bonds, the keys being indices and the values being the vectors. Only bonds present in this dictionary will be used. :type gauges: dict[str, array_like], optional :param optimize: The contraction path optimizer to use, when exactly contracting the local tensors. :type optimize: str or PathOptimizer, optional :param normalized: Whether to normalize the result. If "return", return the norm separately. :type normalized: bool or "return", optional :param max_distance: 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, 1 includes there direct neighbors, etc. :type max_distance: int, optional :param mode: How to select the local tensors, either by graph distance or by selecting the union of all loopy regions containing ``where``, of size up to ``max_distance``, ensuring no dangling tensors. :type mode: {'graphdistance', 'loopunion'}, optional :param fillin: 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 :meth:`~quimb.tensor.tensor_core.TensorNetwork.select_local`. :type fillin: bool or int, optional :param grow_from: If mode is 'loopunion', whether each loop should contain *all* of the initial tagged tensors, or just *any* of them (generating a larger region). :type grow_from: {"all", "any"}, optional :param smudge: A small value to add to the gauges before multiplying them in and inverting them to avoid numerical issues. :type smudge: float, optional :param power: The power to raise the singular values to before multiplying them in and inverting them. :type power: float, optional :param get: Whether to return the result as a fused matrix (i.e. always 2D), unfused array, or still labeled Tensor. :type get: {'matrix', 'array', 'tensor'}, optional :param rehearse: Whether to perform the computation or not, if ``True`` return a rehearsal info dict. :type rehearse: bool, optional :param contract_opts: Supplied to :meth:`~quimb.tensor.tensor_core.TensorNetwork.contract`. .. py:method:: 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) Partially trace this tensor network state, keeping only the sites in ``keep``, using compressed contraction. :param keep: The sites to keep. :type keep: iterable of hashable :param max_bond: The maximum bond dimensions to use while compressed contracting. :type max_bond: int :param optimize: The contraction path optimizer to use, should specifically generate contractions paths designed for compressed contraction. :type optimize: str or PathOptimizer, optional :param flatten: 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``. :type flatten: {False, True, 'all'}, optional :param reduce: Whether to first 'pull' the physical indices off their respective tensors using QR reduction. Experimental. :type reduce: bool, optional :param normalized: Whether to normalize the reduced density matrix at the end. :type normalized: bool, optional :param symmetrized: Whether to symmetrize the reduced density matrix at the end. This should be unecessary if ``flatten`` is set to ``True``. :type symmetrized: {'auto', True, False}, optional :param rehearse: 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``. :type rehearse: {False, 'tn', 'tree', True}, optional :param contract_compressed_opts: Additional keyword arguments to pass to :meth:`~quimb.tensor.tensor_core.TensorNetwork.contract_compressed`. :type contract_compressed_opts: dict, optional :returns: **rho** -- The reduce density matrix of sites in ``keep``. :rtype: array_like .. py:method:: 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) Compute the local expectations of many local operators, by approximately contracting the full overlap tensor network. :param terms: The terms to compute the expectation of, with keys being the sites and values being the local operators. :type terms: dict[node or (node, node), array_like] :param max_bond: The maximum bond dimension to use during contraction. :type max_bond: int :param optimize: The compressed contraction path optimizer to use. :type optimize: str or PathOptimizer :param method: 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. :type method: {'rho', 'rho-reduced'}, optional :param flatten: 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. :type flatten: bool, optional :param normalized: Whether to locally normalize the result. :type normalized: bool, optional :param symmetrized: Whether to symmetrize the reduced density matrix at the end. This should be unecessary if ``flatten`` is set to ``True``. :type symmetrized: {'auto', True, False}, optional :param return_all: Whether to return all results, or just the summed expectation. If ``rehease is not False``, this is ignored and a dict is always returned. :type return_all: bool, optional :param rehearse: 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. :type rehearse: {False, 'tn', 'tree', True}, optional :param executor: If supplied compute the terms in parallel using this executor. :type executor: Executor, optional :param progbar: Whether to show a progress bar. :type progbar: bool, optional :param contract_compressed_opts: Supplied to :meth:`~quimb.tensor.tensor_core.TensorNetwork.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. :rtype: float or dict[node or (node, node), float]