quimb.tensor.tensor_2d¶
Classes and algorithms related to 2D tensor networks.
Attributes¶
Classes¶
Object for rotating coordinates and various contraction functions so |
|
Mixin class for tensor networks with a square lattice two-dimensional |
|
Mixin class for a 2D square lattice vector TN, i.e. one with a single |
|
Mixin class for a 2D square lattice TN operator, i.e. one with both |
|
Mixin class for a 2D square lattice tensor network with a single tensor |
|
Projected Entangled Pair States object (2D): |
|
Projected Entangled Pair Operator object: |
Functions¶
|
|
|
|
|
Convenience function for tiling pairs of bond coordinates on a 2D |
|
Generate a plaquette at site |
|
Generate a tiling of plaquettes in a square 2D lattice. |
|
Generate all length-wise strings in a square 2D lattice. |
|
Ensure |
|
Check if |
|
|
|
|
|
|
|
Base function for printing a unicode schematic of flat 2D TNs. |
|
Find a sequence of plaquette blocksizes that will cover all the terms |
Turn a plaquette |
|
|
Generate a dictionary of all the coordinate pairs in |
|
Generate a string of coordinates, in order, from |
|
Generate the coordinates or a series of swaps that would bring |
|
Generates the ordered long-range path - a sequence of coordinates - from |
|
Module Contents¶
- quimb.tensor.tensor_2d.gen_2d_bonds(Lx, Ly, steppers=None, coo_filter=None, cyclic=False)[source]¶
Convenience function for tiling pairs of bond coordinates on a 2D lattice given a function like
lambda i, j: (i + 1, j + 1)
.- Parameters:
Lx (int) – The number of rows.
Ly (int) – The number of columns.
steppers (callable or sequence of callable, optional) – Function(s) that take args
(i, j)
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)
and only returnsTrue
if this is to be a valid starting coordinate.
- Yields:
bond (tuple[tuple[int, int], tuple[int, int]]) – A pair of coordinates.
Examples
Generate nearest neighbor bonds:
>>> for bond in gen_2d_bonds(2, 2, [lambda i, j: (i, j + 1), >>> lambda i, j: (i + 1, j)]): >>> print(bond) ((0, 0), (0, 1)) ((0, 0), (1, 0)) ((0, 1), (1, 1)) ((1, 0), (1, 1))
Generate next nearest neighbor digonal bonds:
>>> for bond in gen_2d_bonds(2, 2, [lambda i, j: (i + 1, j + 1), >>> lambda i, j: (i + 1, j - 1)]): >>> print(bond) ((0, 0), (1, 1)) ((0, 1), (1, 0))
- quimb.tensor.tensor_2d.gen_2d_plaquette(coo0, steps)[source]¶
Generate a plaquette at site
coo0
by stepping first insteps
and then the reverse steps.- Parameters:
- 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_2d.gen_2d_plaquettes(Lx, Ly, tiling)[source]¶
Generate a tiling of plaquettes in a square 2D lattice.
- Parameters:
Lx (int) – The length of the lattice in the x direction.
Ly (int) – The length of the lattice in the y direction.
tiling ({'1', '2', 'full'}) –
The tiling to use:
- ’1’: plaquettes in a checkerboard pattern, such that each edge
is covered by a maximum of one plaquette.
- ’2’ or ‘full’: dense tiling of plaquettes. All bulk edges will
be covered twice.
- 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_2d.gen_2d_strings(Lx, Ly)[source]¶
Generate all length-wise strings in a square 2D lattice.
- class quimb.tensor.tensor_2d.Rotator2D(tn, xrange, yrange, from_which, stepsize=1)[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.
- Parameters:
tn (TensorNetwork2D) – The tensor network to rotate coordinates for.
xrange (tuple[int, int]) – The range of x-coordinates to range over.
yrange (tuple[int, int]) – The range of y-coordinates to range over.
from_which ({'xmin', 'xmax', 'ymin', 'ymax'}) – The direction to sweep from.
stepsize (int, optional) – The step size to use when sweeping.
- tn¶
- xrange¶
- yrange¶
- from_which¶
- plane¶
- sweep_other()¶
- cyclic_x()¶
- cyclic_y()¶
- quimb.tensor.tensor_2d.BOUNDARY_SEQUENCE_VALID¶
- quimb.tensor.tensor_2d.BOUNDARY_SEQUENCE_MAP¶
- quimb.tensor.tensor_2d.parse_boundary_sequence(sequence)[source]¶
Ensure
sequence
is a tuple of boundary sequence strings from{'xmin', 'xmax', 'ymin', 'ymax'}
- class quimb.tensor.tensor_2d.TensorNetwork2D(ts=(), *, virtual=False, check_collisions=True)[source]¶
Bases:
quimb.tensor.tensor_arbgeom.TensorNetworkGen
Mixin class for tensor networks with a square lattice two-dimensional structure, indexed by
[{row},{column}]
so that:'Y{j}' v i=Lx-1 ●──●──●──●──●──●── ──● | | | | | | | ... | | | | | | 'I{i},{j}' = 'I3,5' e.g. i=3 ●──●──●──●──●──●── | | | | | | | i=2 ●──●──●──●──●──●── ──● <== 'X{i}' | | | | | | ... | i=1 ●──●──●──●──●──●── ──● | | | | | | | i=0 ●──●──●──●──●──●── ──● j=0, 1, 2, 3, 4, 5 j=Ly-1
This implies the following conventions:
the ‘up’ bond is coordinates
(i, j), (i + 1, j)
the ‘down’ bond is coordinates
(i, j), (i - 1, j)
the ‘right’ bond is coordinates
(i, j), (i, j + 1)
the ‘left’ bond is coordinates
(i, j), (i, j - 1)
- _NDIMS = 2¶
- _EXTRA_PROPS = ('_site_tag_id', '_x_tag_id', '_y_tag_id', '_Lx', '_Ly')¶
- _compatible_2d(other)[source]¶
Check whether
self
andother
are compatible 2D tensor networks such that they can remain a 2D tensor network when combined.
- combine(other, *, virtual=False, check_collisions=True)[source]¶
Combine this tensor network with another, returning a new tensor network. If the two are compatible, cast the resulting tensor network to a
TensorNetwork2D
instance.- Parameters:
other (TensorNetwork2D 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:
- property Lx¶
- The number of rows.
- property Ly¶
- The number of columns.
- property nsites¶
- The total number of sites.
- property x_tag_id¶
- The string specifier for tagging each row of this 2D TN.
- property x_tags¶
- A tuple of all of the ``Lx`` different row tags.
- row_tags¶
- property y_tag_id¶
- The string specifier for tagging each column of this 2D TN.
- property y_tags¶
- A tuple of all of the ``Ly`` different column tags.
- col_tags¶
- maybe_convert_coo(x)[source]¶
Check if
x
is a tuple of two 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)
be used for many ‘tag’ based functions.
- gen_horizontal_even_bond_coos()[source]¶
Generate all coordinate pairs like
(i, j), (i, j + 1)
wherej
is even, which thus don’t overlap at all.
- gen_horizontal_odd_bond_coos()[source]¶
Generate all coordinate pairs like
(i, j), (i, j + 1)
wherej
is odd, which thus don’t overlap at all.
- gen_vertical_even_bond_coos()[source]¶
Generate all coordinate pairs like
(i, j), (i + 1, j)
wherei
is even, which thus don’t overlap at all.
- gen_vertical_odd_bond_coos()[source]¶
Generate all coordinate pairs like
(i, j), (i + 1, j)
wherei
is odd, which thus don’t overlap at all.
- gen_diagonal_left_even_bond_coos()[source]¶
Generate all coordinate pairs like
(i, j), (i + 1, j - 1)
wherej
is even, which thus don’t overlap at all.
- gen_diagonal_left_odd_bond_coos()[source]¶
Generate all coordinate pairs like
(i, j), (i + 1, j - 1)
wherej
is odd, which thus don’t overlap at all.
- gen_diagonal_right_even_bond_coos()[source]¶
Generate all coordinate pairs like
(i, j), (i + 1, j + 1)
wherei
is even, which thus don’t overlap at all.
- gen_diagonal_right_odd_bond_coos()[source]¶
Generate all coordinate pairs like
(i, j), (i + 1, j + 1)
wherei
is odd, which thus don’t overlap at all.
- is_cyclic_x(j=None, imin=None, imax=None)[source]¶
Check if the x dimension is cyclic (periodic), specifically whether a bond exists between
(imin, j)
and(imax, j)
, with default values ofimin = 0
andimax = Lx - 1
, andj
at the center of the lattice. Ifimin
andimax
are adjacent then this is considered False, since there is no ‘extra’ connectivity.
- is_cyclic_y(i=None, jmin=None, jmax=None)[source]¶
Check if the y dimension is cyclic (periodic), specifically whether a bond exists between
(i, jmin)
and(i, jmax)
, with default values ofjmin = 0
andjmax = Ly - 1
, andi
at the center of the lattice. Ifjmin
andjmax
are adjacent then this is considered False, since there is no ‘extra’ connectivity.
- _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, xreverse=False, yreverse=False, coordinate_order='xy', xstep=None, ystep=None, stepping_order='xy', 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 and y coordinates.
yrange ((int, int), optional) – The range of allowed values for the x and y coordinates.
xreverse (bool, optional) – Whether to reverse the order of the x and y sweeps.
yreverse (bool, optional) – Whether to reverse the order of the x and y sweeps.
coordinate_order (str, optional) – The order in which to sweep the x and y 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
andyreverse
respectively.ystep (int, optional) – When generating a bond, step in this direction to yield the neighboring coordinate. By default, these follow
xreverse
andyreverse
respectively.stepping_order (str, optional) – The order in which to step the x and y 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)))
- canonize_plane(xrange, yrange, 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.
- canonize_row(i, sweep, yrange=None, **canonize_opts)[source]¶
Canonize all or part of a row.
If
sweep == 'right'
then:| | | | | | | | | | | | | | ─●──●──●──●──●──●──●─ ─●──●──●──●──●──●──●─ | | | | | | | | | | | | | | ─●──●──●──●──●──●──●─ ==> ─●──>──>──>──>──o──●─ row=i | | | | | | | | | | | | | | ─●──●──●──●──●──●──●─ ─●──●──●──●──●──●──●─ | | | | | | | | | | | | | | . . . . jstart jstop jstart jstop
If
sweep == 'left'
then:| | | | | | | | | | | | | | ─●──●──●──●──●──●──●─ ─●──●──●──●──●──●──●─ | | | | | | | | | | | | | | ─●──●──●──●──●──●──●─ ==> ─●──o──<──<──<──<──●─ row=i | | | | | | | | | | | | | | ─●──●──●──●──●──●──●─ ─●──●──●──●──●──●──●─ | | | | | | | | | | | | | | . . . . jstop jstart jstop jstart
Does not yield an orthogonal form in the same way as in 1D.
- canonize_column(j, sweep, xrange=None, **canonize_opts)[source]¶
Canonize all or part of a column.
If
sweep='up'
then:| | | | | | ─●──●──●─ ─●──●──●─ | | | | | | ─●──●──●─ ─●──o──●─ istop | | | ==> | | | ─●──●──●─ ─●──^──●─ | | | | | | ─●──●──●─ ─●──^──●─ istart | | | | | | ─●──●──●─ ─●──●──●─ | | | | | | . . j j
If
sweep='down'
then:| | | | | | ─●──●──●─ ─●──●──●─ | | | | | | ─●──●──●─ ─●──v──●─ istart | | | ==> | | | ─●──●──●─ ─●──v──●─ | | | | | | ─●──●──●─ ─●──o──●─ istop | | | | | | ─●──●──●─ ─●──●──●─ | | | | | | . . j j
Does not yield an orthogonal form in the same way as in 1D.
- compress_plane(xrange, yrange, 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.
- compress_row(i, sweep, yrange=None, max_bond=None, cutoff=1e-10, equalize_norms=False, compress_opts=None)[source]¶
Compress all or part of a row.
If
sweep == 'right'
then:| | | | | | | | | | | | | | ━●━━●━━●━━●━━●━━●━━●━ ━●━━●━━●━━●━━●━━●━━●━ | | | | | | | | | | | | | | ━●━━●━━●━━●━━●━━●━━●━ ━━> ━●━━>──>──>──>──o━━●━ row=i | | | | | | | | | | | | | | ━●━━●━━●━━●━━●━━●━━●━ ━●━━●━━●━━●━━●━━●━━●━ | | | | | | | | | | | | | | . . . . jstart jstop jstart jstop
If
sweep == 'left'
then:| | | | | | | | | | | | | | ━●━━●━━●━━●━━●━━●━━●━ ━●━━●━━●━━●━━●━━●━━●━ | | | | | | | | | | | | | | ━●━━●━━●━━●━━●━━●━━●━ ━━> ━●━━o──<──<──<──<━━●━ row=i | | | | | | | | | | | | | | ━●━━●━━●━━●━━●━━●━━●━ ━●━━●━━●━━●━━●━━●━━●━ | | | | | | | | | | | | | | . . . . jstop jstart jstop jstart
Does not yield an orthogonal form in the same way as in 1D.
- Parameters:
i (int) – Which row to compress.
sweep ({'right', 'left'}) – Which direction to sweep in.
yrange (tuple[int, int] or None) – The range of columns to compress.
max_bond (int, optional) – The maximum boundary dimension, AKA ‘chi’. The default of
None
means truncation is left purely tocutoff
and is not recommended in 2D.cutoff (float, optional) – Cut-off value to used to truncate singular values in the boundary contraction.
compress_opts (None or dict, optional) – Supplied to
compress_between()
.
- compress_column(j, sweep, xrange=None, max_bond=None, cutoff=1e-10, equalize_norms=False, compress_opts=None)[source]¶
Compress all or part of a column.
If
sweep='up'
then:┃ ┃ ┃ ┃ ┃ ┃ ─●──●──●─ ─●──●──●─ ┃ ┃ ┃ ┃ ┃ ┃ ─●──●──●─ ─●──o──●─ . ┃ ┃ ┃ ==> ┃ | ┃ . ─●──●──●─ ─●──^──●─ . xrange ┃ ┃ ┃ ┃ | ┃ . ─●──●──●─ ─●──^──●─ . ┃ ┃ ┃ ┃ ┃ ┃ ─●──●──●─ ─●──●──●─ ┃ ┃ ┃ ┃ ┃ ┃ . . j j
If
sweep='down'
then:┃ ┃ ┃ ┃ ┃ ┃ ─●──●──●─ ─●──●──●─ ┃ ┃ ┃ ┃ ┃ ┃ ─●──●──●─ ─●──v──●─ . ┃ ┃ ┃ ==> ┃ | ┃ . ─●──●──●─ ─●──v──●─ . xrange ┃ ┃ ┃ ┃ | ┃ . ─●──●──●─ ─●──o──●─ . ┃ ┃ ┃ ┃ ┃ ┃ ─●──●──●─ ─●──●──●─ ┃ ┃ ┃ ┃ ┃ ┃ . . j j
Does not yield an orthogonal form in the same way as in 1D.
- Parameters:
j (int) – Which column to compress.
sweep ({'up', 'down'}) – Which direction to sweep in.
xrange (None or (int, int), optional) – The range of rows to compress.
max_bond (int, optional) – The maximum boundary dimension, AKA ‘chi’. The default of
None
means truncation is left purely tocutoff
and is not recommended in 2D.cutoff (float, optional) – Cut-off value to used to truncate singular values in the boundary contraction.
compress_opts (None or dict, optional) – Supplied to
compress_between()
.
- _contract_boundary_core_via_1d(xrange, yrange, from_which, max_bond, cutoff=1e-10, method='dm', layer_tags=None, **compress_opts)[source]¶
- _contract_boundary_core(xrange, yrange, from_which, max_bond, cutoff=1e-10, canonize=True, layer_tags=None, compress_late=True, sweep_reverse=False, equalize_norms=False, compress_opts=None, canonize_opts=None)[source]¶
- _contract_boundary_full_bond(xrange, yrange, from_which, max_bond, cutoff=0.0, method='eigh', renorm=False, optimize='auto-hq', opposite_envs=None, equalize_norms=False, contract_boundary_opts=None)[source]¶
Contract the boundary of this 2D TN using the ‘full bond’ environment information obtained from a boundary contraction in the opposite direction.
- Parameters:
xrange ((int, int) or None, optional) – The range of rows to contract and compress.
yrange ((int, int)) – The range of columns to contract and compress.
from_which ({'xmin', 'ymin', 'xmax', 'ymax'}) – Which direction to contract the rectangular patch from.
max_bond (int) – The maximum boundary dimension, AKA ‘chi’. By default used for the opposite direction environment contraction as well.
cutoff (float, optional) – Cut-off value to used to truncate singular values in the boundary contraction - only for the opposite direction environment contraction.
method ({'eigh', 'eig', 'svd', 'biorthog'}, optional) – Which similarity decomposition method to use to compress the full bond environment.
renorm (bool, optional) – Whether to renormalize the isometric projection or not.
optimize (str or PathOptimize, optimize) – Contraction optimizer to use for the exact contractions.
opposite_envs (dict, optional) – If supplied, the opposite environments will be fetched or lazily computed into this dict depending on whether they are missing.
contract_boundary_opts – Other options given to the opposite direction environment contraction.
- _contract_boundary_projector(xrange, yrange, from_which, max_bond=None, cutoff=1e-10, lazy=False, equalize_norms=False, optimize='auto-hq', compress_opts=None)[source]¶
Contract the boundary of this 2D tensor network by explicitly computing and inserting explicit local projector tensors, which can optionally be left uncontracted. Multilayer networks are naturally supported.
- Parameters:
xrange (tuple) – The range of x indices to contract.
yrange (tuple) – The range of y indices to contract.
from_which ({'xmin', 'xmax', 'ymin', 'ymax'}) – From which boundary to contract.
max_bond (int, optional) – The maximum bond dimension to contract to. If
None
(default), compression is left tocutoff
.cutoff (float, optional) – The cutoff to use for boundary compression.
lazy (bool, optional) – Whether to leave the boundary tensors uncontracted. If
False
(the default), the boundary tensors are contracted and the resulting boundary has a single tensor per site.equalize_norms (bool, optional) – Whether to actively absorb the norm of modified tensors into
self.exponent
.optimize (str or PathOptimizer, optional) – The contract path optimization to use when forming the projector tensors.
compress_opts (dict, optional) – Other options to pass to
svd_truncated()
.
See also
TensorNetwork.insert_compressor_between_regions
- contract_boundary_from(xrange, yrange, from_which, max_bond=None, *, cutoff=1e-10, canonize=True, mode='mps', layer_tags=None, sweep_reverse=False, compress_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_boundary_from_xmin(xrange, yrange=None, max_bond=None, *, cutoff=1e-10, canonize=True, mode='mps', layer_tags=None, sweep_reverse=False, compress_opts=None, inplace=False, **contract_boundary_opts)[source]¶
Contract a 2D tensor network inwards from the bottom, canonizing and compressing (left to right) along the way. If
layer_tags is None
this looks like:a) contract │ │ │ │ │ ●──●──●──●──● │ │ │ │ │ │ │ │ │ │ --> ●══●══●══●══● ●──●──●──●──● b) optionally canonicalize │ │ │ │ │ ●══●══<══<══< c) compress in opposite direction │ │ │ │ │ --> │ │ │ │ │ --> │ │ │ │ │ >──●══●══●══● --> >──>──●══●══● --> >──>──>──●══● . . --> . . --> . .
If
layer_tags
is specified, each then each layer is contracted in and compressed separately, resulting generally in a lower memory scaling. For two layer tags this looks like:a) first flatten the outer boundary only │ ││ ││ ││ ││ │ │ ││ ││ ││ ││ │ ●─○●─○●─○●─○●─○ ●─○●─○●─○●─○●─○ │ ││ ││ ││ ││ │ ==> ╲│ ╲│ ╲│ ╲│ ╲│ ●─○●─○●─○●─○●─○ ●══●══●══●══● b) contract and compress a single layer only │ ││ ││ ││ ││ │ │ ○──○──○──○──○ │╱ │╱ │╱ │╱ │╱ ●══<══<══<══< c) contract and compress the next layer ╲│ ╲│ ╲│ ╲│ ╲│ >══>══>══>══●
- Parameters:
xrange ((int, int)) – The range of rows to compress (inclusive).
yrange ((int, int) or None, optional) – The range of columns to compress (inclusive), sweeping along with canonization and compression. Defaults to all columns.
max_bond (int, optional) – The maximum boundary dimension, AKA ‘chi’. The default of
None
means truncation is left purely tocutoff
and is not recommended in 2D.cutoff (float, optional) – Cut-off value to used to truncate singular values in the boundary contraction.
canonize (bool, optional) – Whether to sweep one way with canonization before compressing.
mode ({'mps', 'full-bond'}, optional) – How to perform the compression on the boundary.
layer_tags (None or sequence[str], optional) – If
None
, all tensors at each coordinate pair[(i, j), (i + 1, j)]
will be first contracted. If specified, then the outer tensor at(i, j)
will be contracted with the tensor specified by[(i + 1, j), layer_tag]
, for eachlayer_tag
inlayer_tags
.sweep_reverse (bool, optional) – Which way to perform the compression sweep, which has an effect on which tensors end up being canonized. Setting this to true sweeps the compression from largest to smallest coordinates.
compress_opts (None or dict, optional) – Supplied to
compress_between()
.inplace (bool, optional) – Whether to perform the contraction inplace or not.
- contract_boundary_from_xmax(xrange, yrange=None, max_bond=None, *, cutoff=1e-10, canonize=True, mode='mps', layer_tags=None, inplace=False, sweep_reverse=False, compress_opts=None, **contract_boundary_opts)[source]¶
Contract a 2D tensor network inwards from the top, canonizing and compressing (right to left) along the way. If
layer_tags is None
this looks like:a) contract ●──●──●──●──● | | | | | --> ●══●══●══●══● ●──●──●──●──● | | | | | | | | | | b) optionally canonicalize ●══●══<══<══< | | | | | c) compress in opposite direction >──●══●══●══● --> >──>──●══●══● --> >──>──>──●══● | | | | | --> | | | | | --> | | | | | . . --> . . --> . .
If
layer_tags
is specified, each then each layer is contracted in and compressed separately, resulting generally in a lower memory scaling. For two layer tags this looks like:a) first flatten the outer boundary only ●─○●─○●─○●─○●─○ ●══●══●══●══● │ ││ ││ ││ ││ │ ==> ╱│ ╱│ ╱│ ╱│ ╱│ ●─○●─○●─○●─○●─○ ●─○●─○●─○●─○●─○ │ ││ ││ ││ ││ │ │ ││ ││ ││ ││ │ b) contract and compress a single layer only ●══<══<══<══< │╲ │╲ │╲ │╲ │╲ │ ○──○──○──○──○ │ ││ ││ ││ ││ │ c) contract and compress the next layer ●══●══●══●══● ╱│ ╱│ ╱│ ╱│ ╱│
- Parameters:
xrange ((int, int)) – The range of rows to compress (inclusive).
yrange ((int, int) or None, optional) – The range of columns to compress (inclusive), sweeping along with canonization and compression. Defaults to all columns.
max_bond (int, optional) – The maximum boundary dimension, AKA ‘chi’. The default of
None
means truncation is left purely tocutoff
and is not recommended in 2D.cutoff (float, optional) – Cut-off value to used to truncate singular values in the boundary contraction.
canonize (bool, optional) – Whether to sweep one way with canonization before compressing.
mode ({'mps', 'full-bond'}, optional) – How to perform the compression on the boundary.
layer_tags (None or str, optional) – If
None
, all tensors at each coordinate pair[(i, j), (i - 1, j)]
will be first contracted. If specified, then the outer tensor at(i, j)
will be contracted with the tensor specified by[(i - 1, j), layer_tag]
, for eachlayer_tag
inlayer_tags
.sweep_reverse (bool, optional) – Which way to perform the compression sweep, which has an effect on which tensors end up being canonized. Setting this to true sweeps the compression from largest to smallest coordinates.
compress_opts (None or dict, optional) – Supplied to
compress_between()
.inplace (bool, optional) – Whether to perform the contraction inplace or not.
- contract_boundary_from_ymin(yrange, xrange=None, max_bond=None, *, cutoff=1e-10, canonize=True, mode='mps', layer_tags=None, sweep_reverse=False, compress_opts=None, inplace=False, **contract_boundary_opts)[source]¶
Contract a 2D tensor network inwards from the left, canonizing and compressing (bottom to top) along the way. If
layer_tags is None
this looks like:a) contract ●──●── ●── │ │ ║ ●──●── ==> ●── │ │ ║ ●──●── ●── b) optionally canonicalize ●── v── ║ ║ ●── ==> v── ║ ║ ●── ●── c) compress in opposite direction v── ●── ║ │ v── ==> ^── ║ │ ●── ^──
If
layer_tags
is specified, each then each layer is contracted in and compressed separately, resulting generally in a lower memory scaling. For two layer tags this looks like:a) first flatten the outer boundary only ○──○── ●──○── │╲ │╲ │╲ │╲ ●─○──○── ╰─●──○── ╲│╲╲│╲ ==> │╲╲│╲ ●─○──○── ╰─●──○── ╲│ ╲│ │ ╲│ ●──●── ╰──●── b) contract and compress a single layer only ○── ╱╱ ╲ ●─── ○── ╲ ╱╱ ╲ ^─── ○── ╲ ╱╱ ^───── c) contract and compress the next layer ●── │╲ ╰─●── │╲ ╰─●── │ ╰──
- Parameters:
yrange ((int, int)) – The range of columns to compress (inclusive).
xrange ((int, int) or None, optional) – The range of rows to compress (inclusive), sweeping along with canonization and compression. Defaults to all rows.
max_bond (int, optional) – The maximum boundary dimension, AKA ‘chi’. The default of
None
means truncation is left purely tocutoff
and is not recommended in 2D.cutoff (float, optional) – Cut-off value to used to truncate singular values in the boundary contraction.
canonize (bool, optional) – Whether to sweep one way with canonization before compressing.
mode ({'mps', 'full-bond'}, optional) – How to perform the compression on the boundary.
layer_tags (None or str, optional) – If
None
, all tensors at each coordinate pair[(i, j), (i, j + 1)]
will be first contracted. If specified, then the outer tensor at(i, j)
will be contracted with the tensor specified by[(i + 1, j), layer_tag]
, for eachlayer_tag
inlayer_tags
.sweep_reverse (bool, optional) – Which way to perform the compression sweep, which has an effect on which tensors end up being canonized. Setting this to true sweeps the compression from largest to smallest coordinates.
compress_opts (None or dict, optional) – Supplied to
compress_between()
.inplace (bool, optional) – Whether to perform the contraction inplace or not.
- contract_boundary_from_ymax(yrange, xrange=None, max_bond=None, *, cutoff=1e-10, canonize=True, mode='mps', layer_tags=None, sweep_reverse=False, compress_opts=None, inplace=False, **contract_boundary_opts)[source]¶
Contract a 2D tensor network inwards from the left, canonizing and compressing (top to bottom) along the way. If
layer_tags is None
this looks like:a) contract ──●──● ──● │ │ ║ ──●──● ==> ──● │ │ ║ ──●──● ──● b) optionally canonicalize ──● ──v ║ ║ ──● ==> ──v ║ ║ ──● ──● c) compress in opposite direction ──v ──● ║ │ ──v ==> ──^ ║ │ ──● ──^
If
layer_tags
is specified, each then each layer is contracted in and compressed separately, resulting generally in a lower memory scaling. For two layer tags this looks like:a) first flatten the outer boundary only ──○──○ ──○──● ╱│ ╱│ ╱│ ╱│ ──○──○─● ──○──●─╯ ╱│╱╱│╱ ==> ╱│╱╱│ ──○──○─● ──○──●─╯ │╱ │╱ │╱ │ ──●──● ──●──╯ b) contract and compress a single layer only ──○ ╱ ╲╲ ──○────v ╱ ╲╲ ╱ ──○────v ╲╲ ╱ ─────● c) contract and compress the next layer ╲ ────v ╲ ╱ ────v ╲ ╱ ────●
- Parameters:
yrange ((int, int)) – The range of columns to compress (inclusive).
xrange ((int, int) or None, optional) – The range of rows to compress (inclusive), sweeping along with canonization and compression. Defaults to all rows.
max_bond (int, optional) – The maximum boundary dimension, AKA ‘chi’. The default of
None
means truncation is left purely tocutoff
and is not recommended in 2D.cutoff (float, optional) – Cut-off value to used to truncate singular values in the boundary contraction.
canonize (bool, optional) – Whether to sweep one way with canonization before compressing.
mode ({'mps', 'full-bond'}, optional) – How to perform the compression on the boundary.
layer_tags (None or str, optional) – If
None
, all tensors at each coordinate pair[(i, j), (i, j - 1)]
will be first contracted. If specified, then the outer tensor at(i, j)
will be contracted with the tensor specified by[(i + 1, j), layer_tag]
, for eachlayer_tag
inlayer_tags
.sweep_reverse (bool, optional) – Which way to perform the compression sweep, which has an effect on which tensors end up being canonized. Setting this to true sweeps the compression from largest to smallest coordinates.
compress_opts (None or dict, optional) – Supplied to
compress_between()
.inplace (bool, optional) – Whether to perform the contraction inplace or not.
- _contract_interleaved_boundary_sequence(*, contract_boundary_opts=None, sequence=None, xmin=None, xmax=None, ymin=None, ymax=None, max_separation=1, max_unfinished=1, around=None, equalize_norms=False, canonize=False, canonize_opts=None, final_contract=True, final_contract_opts=None, 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, canonize=True, mode='mps', layer_tags=None, compress_opts=None, sequence=None, xmin=None, xmax=None, ymin=None, ymax=None, max_separation=1, max_unfinished=1, around=None, equalize_norms=False, final_contract=True, final_contract_opts=None, progbar=None, inplace=False, **contract_boundary_opts)[source]¶
Contract the boundary of this 2D tensor network inwards:
●──●──●──● ●──●──●──● ●──●──● │ │ │ │ │ │ │ │ ║ │ │ ●──●──●──● ●──●──●──● ^──●──● >══>══● >──v │ │ij│ │ ==> │ │ij│ │ ==> ║ij│ │ ==> │ij│ │ ==> │ij║ ●──●──●──● ●══<══<══< ^──<──< ^──<──< ^──< │ │ │ │ ●──●──●──●
Optionally from any or all of the boundary, in multiple layers, and stopping around a region. The default is to contract the boundary from the two shortest opposing sides.
- Parameters:
around (None or sequence of (int, int), optional) – If given, don’t contract the square of sites bounding these coordinates.
max_bond (int, optional) – The maximum boundary dimension, AKA ‘chi’. The default of
None
means truncation is left purely tocutoff
and is not recommended in 2D.cutoff (float, optional) – Cut-off value to used to truncate singular values in the boundary contraction.
canonize (bool, optional) – Whether to sweep one way with canonization before compressing.
mode ({'mps', 'full-bond'}, optional) – How to perform the compression on the boundary.
layer_tags (None or sequence of str, optional) – If given, perform a multilayer contraction, contracting the inner sites in each layer into the boundary individually.
compress_opts (None or dict, optional) – Other low level options to pass to
compress_between()
.sequence (sequence of {'xmin', 'xmax', 'ymin', 'ymax'}, optional) – Which directions to cycle throught when performing the inwards contractions, i.e. from that direction. If
around
is specified you will likely need all of these! Default is to contract from the two shortest opposing sides.xmin (int, optional) – The initial bottom boundary row, defaults to 0.
xmax (int, optional) – The initial top boundary row, defaults to
Lx - 1
.ymin (int, optional) – The initial left boundary column, defaults to 0.
ymax (int, optional) – The initial right boundary column, defaults to
Ly - 1
..max_separation (int, optional) – If
around is None
, when any two sides become this far apart simply contract the remaining tensor network.max_unfinished (int, optional) – If
around is None
, when this many sides are still not withinmax_separation
simply contract the remaining tensor network.around – If given, don’t contract the square of sites bounding these coordinates.
equalize_norms (bool or float, optional) – Whether to equalize the norms of the boundary tensors after each contraction, gathering the overall scaling coefficient, log10, in
tn.exponent
.final_contract (bool, optional) – Whether to exactly contract the remaining tensor network after the boundary contraction.
final_contract_opts (None or dict, optional) – Options to pass to
contract()
,optimize
defaults to'auto-hq'
.progbar (bool, optional) – Whether to show a progress bar.
inplace (bool, optional) – Whether to perform the contraction in place or not.
contract_boundary_opts – Supplied to
contract_boundary_from()
, including compression and canonization options.
- contract_mps_sweep(max_bond=None, *, cutoff=1e-10, canonize=True, direction=None, **contract_boundary_opts)[source]¶
Contract this 2D tensor network by sweeping an MPS across from one side to the other.
- Parameters:
max_bond (int, optional) – The maximum boundary dimension, AKA ‘chi’. The default of
None
means truncation is left purely tocutoff
and is not recommended in 2D.cutoff (float, optional) – Cut-off value to used to truncate singular values in the boundary contraction.
canonize (bool, optional) – Whether to sweep one way with canonization before compressing.
direction ({'xmin', 'xmax', 'ymin', 'ymax'}, optional) – Which direction to sweep from. If
None
(default) then the shortest boundary is chosen.contract_boundary_opts – Supplied to
contract_boundary_from()
, including compression and canonization options.
- compute_environments(from_which, xrange=None, yrange=None, max_bond=None, *, cutoff=1e-10, canonize=True, mode='mps', layer_tags=None, dense=False, compress_opts=None, envs=None, **contract_boundary_opts)[source]¶
Compute the 1D boundary tensor networks describing the environments of rows and columns.
- Parameters:
from_which ({'xmin', 'xmax', 'ymin', 'ymax'}) – Which boundary to compute the environments from.
xrange (tuple[int], optional) – The range of rows to compute the environments for.
yrange (tuple[int], optional) – The range of columns to compute the environments for.
max_bond (int, optional) – The maximum bond dimension of the environments.
cutoff (float, optional) – The cutoff for the singular values of the environments.
canonize (bool, optional) – Whether to canonicalize along each MPS environment before compressing.
mode ({'mps', 'projector', 'full-bond'}, optional) – Which contraction method to use for the environments.
layer_tags (str or iterable[str], optional) – If this 2D TN is multi-layered (e.g. a bra and a ket), and
mode == 'mps'
, contract and compress each specified layer separately, for a cheaper contraction.dense (bool, optional) – Whether to use dense tensors for the environments.
compress_opts (dict, optional) – Other options to pass to
tensor_compress_bond()
.envs (dict, optional) – An existing dictionary to store the environments in.
contract_boundary_opts – Other options to pass to
contract_boundary_from()
.
- Returns:
envs – A dictionary of the environments, with keys of the form
(from_which, row_or_col_index)
.- Return type:
- compute_xmin_environments[source]¶
Compute the
self.Lx
1D boundary tensor networks describing the lower environments of each row in this 2D tensor network. Seecompute_x_environments()
for full details.
- compute_xmax_environments[source]¶
Compute the
self.Lx
1D boundary tensor networks describing the upper environments of each row in this 2D tensor network. Seecompute_x_environments()
for full details.
- compute_ymin_environments[source]¶
Compute the
self.Ly
1D boundary tensor networks describing the left environments of each column in this 2D tensor network. Seecompute_y_environments()
for full details.
- compute_ymax_environments[source]¶
Compute the
self.Ly
1D boundary tensor networks describing the right environments of each column in this 2D tensor network. Seecompute_y_environments()
for full details.
- compute_x_environments(max_bond=None, *, cutoff=1e-10, canonize=True, dense=False, mode='mps', layer_tags=None, compress_opts=None, envs=None, **contract_boundary_opts)[source]¶
Compute the
2 * self.Lx
1D boundary tensor networks describing the lower and upper environments of each row in this 2D tensor network, assumed to represent the norm.The top or ‘xmax’ environment for row
i
will be a contraction of all rowsi + 1, i + 2, ...
etc:●━━━●━━━●━━━●━━━●━━━●━━━●━━━●━━━●━━━● ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲
The bottom or ‘xmin’ environment for row
i
will be a contraction of all rowsi - 1, i - 2, ...
etc:╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ●━━━●━━━●━━━●━━━●━━━●━━━●━━━●━━━●━━━●
Such that
envs['xmax', i] & self.select(self.x_tag(i)) & envs['xmin', i]
would look like:●━━━●━━━●━━━●━━━●━━━●━━━●━━━●━━━●━━━● ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ o─o─o─o─o─o─o─o─o─o─o─o─o─o─o─o─o─o─o─o ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ╲ ╱ ●━━━●━━━●━━━●━━━●━━━●━━━●━━━●━━━●━━━●
And be (an approximation of) the norm centered around row
i
- Parameters:
max_bond (int, optional) – The maximum boundary dimension, AKA ‘chi’. The default of
None
means truncation is left purely tocutoff
and is not recommended in 2D.cutoff (float, optional) – Cut-off value to used to truncate singular values in the boundary contraction.
canonize (bool, optional) – Whether to sweep one way with canonization before compressing.
dense (bool, optional) – If true, contract the boundary in as a single dense tensor.
mode ({'mps', 'full-bond'}, optional) – How to perform the boundary compression.
layer_tags (None or sequence[str], optional) – If
None
, all tensors at each coordinate pair[(i, j), (i + 1, j)]
will be first contracted. If specified, then the outer tensor at(i, j)
will be contracted with the tensor specified by[(i + 1, j), layer_tag]
, for eachlayer_tag
inlayer_tags
.compress_opts (None or dict, optional) – Supplied to
compress_between()
.envs (dict, optional) – Supply an existing dictionary to store the environments in.
contract_boundary_opts – Supplied to
contract_boundary_from_xmin()
andcontract_boundary_from_xmax()
.
- Returns:
x_envs – The two environment tensor networks of row
i
will be stored inx_envs['xmin', i]
andx_envs['xmax', i]
.- Return type:
dict[(str, int), TensorNetwork]
- compute_y_environments(max_bond=None, *, cutoff=1e-10, canonize=True, dense=False, mode='mps', layer_tags=None, compress_opts=None, envs=None, **contract_boundary_opts)[source]¶
Compute the
2 * self.Ly
1D boundary tensor networks describing the left (‘ymin’) and right (‘ymax’) environments of each column in this 2D tensor network, assumed to represent the norm.The left or ‘ymin’ environment for column
j
will be a contraction of all columnsj - 1, j - 2, ...
etc:●< ┃ ●< ┃ ●< ┃ ●<
The right or ‘ymax’ environment for row
j
will be a contraction of all rowsj + 1, j + 2, ...
etc:>● ┃ >● ┃ >● ┃ >●
Such that
envs['ymin', j] & self.select(self.y_tag(j)) & envs['ymax', j]
would look like:╱o ●< o| >● ┃ |o ┃ ●< o| >● ┃ |o ┃ ●< o| >● ┃ |o ┃ ●< o╱ >●
And be (an approximation of) the norm centered around column
j
- Parameters:
max_bond (int, optional) – The maximum boundary dimension, AKA ‘chi’. The default of
None
means truncation is left purely tocutoff
and is not recommended in 2D.cutoff (float, optional) – Cut-off value to used to truncate singular values in the boundary contraction.
canonize (bool, optional) – Whether to sweep one way with canonization before compressing.
dense (bool, optional) – If true, contract the boundary in as a single dense tensor.
mode ({'mps', 'full-bond'}, optional) – How to perform the boundary compression.
layer_tags (None or sequence[str], optional) – If
None
, all tensors at each coordinate pair[(i, j), (i + 1, j)]
will be first contracted. If specified, then the outer tensor at(i, j)
will be contracted with the tensor specified by[(i + 1, j), layer_tag]
, for eachlayer_tag
inlayer_tags
.compress_opts (None or dict, optional) – Supplied to
compress_between()
.contract_boundary_opts – Supplied to
contract_boundary_from_ymin()
andcontract_boundary_from_ymax()
.
- Returns:
y_envs – The two environment tensor networks of column
j
will be stored iny_envs['ymin', j]
andy_envs['ymax', j]
.- Return type:
dict[(str, int), TensorNetwork]
- _compute_plaquette_environments_x_first(x_bsz, y_bsz, max_bond=None, cutoff=1e-10, canonize=True, layer_tags=None, second_dense=None, x_envs=None, **compute_environment_opts)[source]¶
- _compute_plaquette_environments_y_first(x_bsz, y_bsz, max_bond=None, cutoff=1e-10, canonize=True, layer_tags=None, second_dense=None, y_envs=None, **compute_environment_opts)[source]¶
- compute_plaquette_environments(x_bsz=2, y_bsz=2, max_bond=None, *, cutoff=1e-10, canonize=True, mode='mps', layer_tags=None, first_contract=None, second_dense=None, compress_opts=None, **compute_environment_opts)[source]¶
Compute all environments like:
second_dense=False second_dense=True (& first_contract='columns') ●──● ╭───●───╮ ╱│ │╲ │ ╱ ╲ │ ●─. .─● ┬ ●─ . . ─● ┬ │ │ ┊ x_bsz │ │ ┊ x_bsz ●─. .─● ┴ ●─ . . ─● ┴ ╲│ │╱ │ ╲ ╱ │ ●──● ╰───●───╯ <--> <-> y_bsz y_bsz
Use two boundary contractions sweeps.
- Parameters:
x_bsz (int, optional) – The size of the plaquettes in the x-direction (number of rows).
y_bsz (int, optional) – The size of the plaquettes in the y-direction (number of columns).
max_bond (int, optional) – The maximum boundary dimension, AKA ‘chi’. The default of
None
means truncation is left purely tocutoff
and is not recommended in 2D.cutoff (float, optional) – Cut-off value to used to truncate singular values in the boundary contraction.
canonize (bool, optional) – Whether to sweep one way with canonization before compressing.
mode ({'mps', 'full-bond'}, optional) – How to perform the boundary compression.
layer_tags (None or sequence[str], optional) – If
None
, all tensors at each coordinate pair[(i, j), (i + 1, j)]
will be first contracted. If specified, then the outer tensor at(i, j)
will be contracted with the tensor specified by[(i + 1, j), layer_tag]
, for eachlayer_tag
inlayer_tags
.first_contract ({None, 'x', 'y'}, optional) – The environments can either be generated with initial sweeps in the row (‘x’) or column (‘y’) direction. Generally it makes sense to perform this approximate step in whichever is smaller (the default).
second_dense (None or bool, optional) – Whether to perform the second set of contraction sweeps (in the rotated direction from whichever
first_contract
is) using a dense tensor or boundary method. By default this is only turned on if thebsz
in the corresponding direction is 1.compress_opts (None or dict, optional) – Supplied to
compress_between()
.compute_environment_opts – Supplied to
compute_y_environments()
orcompute_x_environments()
.
- Returns:
The plaquette environments. The key is two tuples of ints, the startings coordinate of the plaquette being the first and the size of the plaquette being the second pair.
- Return type:
- 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'}) – 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:
See also
contract_hotrg
,TensorNetwork.insert_compressor_between_regions
- contract_hotrg(max_bond=None, *, cutoff=1e-10, canonize=False, canonize_opts=None, sequence=('x', 'y'), 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 plaquettes, 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.
canonize (bool, optional) – Whether to canonize all tensors before each contraction, via
gauge_all()
.canonize_opts (None or dict, optional) – Additional options to pass to
gauge_all()
.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 all tensors after each contraction, gathering the overall scaling coefficient, log10, in
tn.exponent
.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'
.progbar (bool, optional) – Whether to show a progress bar.
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 direction of length > 2.
- Return type:
See also
coarse_grain_hotrg
,contract_ctmrg
,TensorNetwork.insert_compressor_between_regions
- 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, max_separation=1, around=None, equalize_norms=False, final_contract=True, final_contract_opts=None, progbar=False, inplace=False, **contract_boundary_opts)[source]¶
Contract this 2D tensor network using the finite analog of the CTMRG algorithm - https://arxiv.org/abs/cond-mat/9507087. The TN is contracted sequentially in
sequence
directions by inserting oblique projectors between boundary 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.
canonize (bool, optional) – Whether to canonize the boundary tensors before each contraction, via
gauge_all()
.canonize_opts (None or dict, optional) – Additional options to pass to
gauge_all()
.lazy (bool, optional) – Whether to contract the coarse graining projectors or leave them in the tensor network lazily. Default is to contract them.
mode (str, optional) – The method to perform the boundary contraction. Defaults to
'projector'
.compress_opts (None or dict, optional) – Other low level options to pass to
insert_compressor_between_regions()
.sequence (sequence of {'xmin', 'xmax', 'ymin', 'ymax'}, optional) – Which directions to cycle throught when performing the inwards contractions, i.e. from that direction. If
around
is specified you will likely need all of these! Default is to contract in all directions.xmin (int, optional) – The initial bottom boundary row, defaults to 0.
xmax (int, optional) – The initial top boundary row, defaults to
Lx - 1
.ymin (int, optional) – The initial left boundary column, defaults to 0.
ymax (int, optional) – The initial right boundary column, defaults to
Ly - 1
..max_separation (int, optional) – If
around is None
, when any two sides become this far apart simply contract the remaining tensor network.around (None or sequence of (int, int), optional) – If given, don’t contract the square of sites bounding these coordinates.
equalize_norms (bool or float, optional) – Whether to equalize the norms of the boundary tensors after each contraction, gathering the overall scaling coefficient, log10, in
tn.exponent
.final_contract (bool, optional) – Whether to exactly contract the remaining tensor network after the boundary contraction.
final_contract_opts (None or dict, optional) – Options to pass to
contract()
,optimize
defaults to'auto-hq'
.progbar (bool, optional) – Whether to show a progress bar.
inplace (bool, optional) – Whether to perform the boundary contraction in place.
contract_boundary_opts – Additional options to pass to
contract_boundary_from()
.
- Returns:
Either the fully contracted scalar (if
final_contract=True
andaround=None
) or the partially contracted tensor network.- Return type:
scalar or TensorNetwork2D
See also
contract_boundary_from
,contract_hotrg
,TensorNetwork.insert_compressor_between_regions
- quimb.tensor.tensor_2d.is_lone_coo(where)[source]¶
Check if
where
has been specified as a single coordinate pair.
- quimb.tensor.tensor_2d.gate_string_split_(TG, where, string, original_ts, bonds_along, reindex_map, site_ix, info, **compress_opts)[source]¶
- quimb.tensor.tensor_2d.gate_string_reduce_split_(TG, where, string, original_ts, bonds_along, reindex_map, site_ix, info, **compress_opts)[source]¶
- quimb.tensor.tensor_2d.gate_2d_long_range(self, G, where, contract=False, tags=None, inplace=False, info=None, long_range_use_swaps=False, long_range_path_sequence=None, **compress_opts)[source]¶
- class quimb.tensor.tensor_2d.TensorNetwork2DVector(ts=(), *, virtual=False, check_collisions=True)[source]¶
Bases:
TensorNetwork2D
,quimb.tensor.tensor_arbgeom.TensorNetworkGenVector
Mixin class for a 2D 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', '_Lx', '_Ly', '_site_ind_id')¶
- 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
).
- gate(G, where, contract=False, tags=None, propagate_tags='sites', inplace=False, info=None, long_range_use_swaps=False, long_range_path_sequence=None, **compress_opts)[source]¶
Apply the dense gate
G
, maintaining the physical indices of this 2D vector tensor network.- Parameters:
G (array_like) – The gate array to apply, should match or be factorable into the shape
(phys_dim,) * (2 * len(where))
.where (sequence of tuple[int, int] or tuple[int, int]) – Which site coordinates to apply the gate to.
contract ({'reduce-split', 'split', False, True}, optional) –
How to contract the gate into the 2D tensor network:
False: gate is added to network and nothing is contracted, tensor network structure is thus not maintained.
True: gate is contracted with all tensors involved, tensor network structure is thus only maintained if gate acts on a single site only.
’split’: contract all involved tensors then split the result back into two.
’reduce-split’: factor the two physical indices into ‘R-factors’ using QR decompositions on the original site tensors, then contract the gate, split it and reabsorb each side. Much cheaper than
'split'
.
The final two methods are relevant for two site gates only, for single site gates they use the
contract=True
option which also maintains the structure of the TN. See below for a pictorial description of each method.tags (str or sequence of str, optional) – Tags to add to the new gate tensor.
propagate_tags ({'sites', 'register', True, False}, optional) –
If
contract==False
, which tags to propagate to the new gate tensor from the tensors it was applied to:If
'sites'
, then only propagate tags matching e.g. ‘I{},{}’ and ignore all others. I.e. assuming unitary gates just propagate the causal lightcone.If
'register'
, then only propagate tags matching the sites of where this gate was actually applied. I.e. ignore the lightcone, just keep track of which ‘registers’ the gate was applied to.If
False
, propagate nothing.If
True
, propagate all tags.
inplace (bool, optional) – Whether to perform the gate operation inplace on the tensor network or not.
info (None or dict, optional) – Used to store extra optional information such as the singular values if not absorbed.
compress_opts – Supplied to
tensor_split()
for anycontract
methods that involve splitting. Ignored otherwise.
- Returns:
G_psi – The new 2D vector TN like
IIIGII @ psi
etc.- Return type:
Notes
The
contract
options look like the following (for two site gates).contract=False
:│ │ GGGGG │╱ │╱ ──●───●── ╱ ╱
contract=True
:│╱ │╱ ──GGGGG── ╱ ╱
contract='split'
:│╱ │╱ │╱ │╱ ──GGGGG── ==> ──G┄┄┄G── ╱ ╱ ╱ ╱ <SVD>
contract='reduce-split'
:│ │ │ │ GGGGG GGG │ │ │╱ │╱ ==> ╱│ │ ╱ ==> ╱│ │ ╱ │╱ │╱ ──●───●── ──>─●─●─<── ──>─GGG─<── ==> ──G┄┄┄G── ╱ ╱ ╱ ╱ ╱ ╱ ╱ ╱ <QR> <LQ> <SVD>
For one site gates when one of the ‘split’ methods is supplied
contract=True
is assumed.
- compute_norm(layer_tags=('KET', 'BRA'), **contract_opts)[source]¶
Compute the norm of this vector via boundary contraction.
- compute_local_expectation(terms, max_bond=None, *, cutoff=1e-10, canonize=True, mode='mps', layer_tags=('KET', 'BRA'), normalized=False, autogroup=True, contract_optimize='auto-hq', return_all=False, plaquette_envs=None, plaquette_map=None, **plaquette_env_options)[source]¶
Compute the sum of many local expecations by essentially forming the reduced density matrix of all required plaquettes. If you supply
normalized=True
each expecation is locally normalized, which a) is usually more accurate and b) doesn’t require a separate normalization boundary contraction.- Parameters:
terms (dict[tuple[tuple[int], array]) – A dictionary mapping site coordinates to raw operators, which will be supplied to
gate()
. The keys should either be a single coordinate -(i, j)
- describing a single site operator, or a pair of coordinates -((i_a, j_a), (i_b, j_b))
describing a two site operator.max_bond (int, optional) – The maximum boundary dimension, AKA ‘chi’. The default of
None
means truncation is left purely tocutoff
and is not recommended in 2D.cutoff (float, optional) – Cut-off value to used to truncate singular values in the boundary contraction.
canonize (bool, optional) – Whether to sweep one way with canonization before compressing.
mode ({'mps', 'full-bond'}, optional) – How to perform the compression on the boundary.
layer_tags (None or sequence of str, optional) – If given, perform a multilayer contraction, contracting the inner sites in each layer into the boundary individually.
normalized (bool, optional) – If True, normalize the value of each local expectation by the local norm: $langle O_i rangle = Tr[rho_p O_i] / Tr[rho_p]$.
autogroup (bool, optional) – If
True
(the default), group terms into horizontal and vertical sets to be computed separately (usually more efficient) if possible.contract_optimize (str, optional) – Contraction path finder to use for contracting the local plaquette expectation (and optionally normalization).
return_all (bool, optional) – Whether to the return all the values individually as a dictionary of coordinates to tuple[local_expectation, local_norm].
plaquette_envs (None or dict, optional) – Supply precomputed plaquette environments.
plaquette_map (None, dict, optional) – Supply the mapping of which plaquettes (denoted by
((x0, y0), (dx, dy))
) to use for which coordinates, it will be calculated automatically otherwise.plaquette_env_options – Supplied to
compute_plaquette_environments()
to generate the plaquette environments, equivalent to approximately performing the partial trace.
- Return type:
scalar or dict
- normalize(max_bond=None, *, cutoff=1e-10, canonize=True, mode='mps', layer_tags=('KET', 'BRA'), balance_bonds=False, equalize_norms=False, inplace=False, **contract_boundary_opts)[source]¶
Normalize this PEPS.
- Parameters:
max_bond (int, optional) – The maximum boundary dimension, AKA ‘chi’. The default of
None
means truncation is left purely tocutoff
and is not recommended in 2D.cutoff (float, optional) – Cut-off value to used to truncate singular values in the boundary contraction.
canonize (bool, optional) – Whether to sweep one way with canonization before compressing.
mode ({'mps', 'full-bond'}, optional) – How to perform the compression on the boundary.
layer_tags (None or sequence of str, optional) – If given, perform a multilayer contraction, contracting the inner sites in each layer into the boundary individually.
balance_bonds (bool, optional) – Whether to balance the bonds after normalization, a form of conditioning.
equalize_norms (bool, optional) – Whether to set all the tensor norms to the same value after normalization, another form of conditioning.
inplace (bool, optional) – Whether to perform the normalization inplace or not.
contract_boundary_opts – Supplied to
contract_boundary()
, by default, two layer contraction will be used.
- class quimb.tensor.tensor_2d.TensorNetwork2DOperator(ts=(), *, virtual=False, check_collisions=True)[source]¶
Bases:
TensorNetwork2D
,quimb.tensor.tensor_arbgeom.TensorNetworkGenOperator
Mixin class for a 2D square lattice TN operator, i.e. one with both ‘upper’ and ‘lower’ site (physical) indices.
- _EXTRA_PROPS = ('_site_tag_id', '_x_tag_id', '_y_tag_id', '_Lx', '_Ly', '_upper_ind_id', '_lower_ind_id')¶
- reindex_lower_sites(new_id, where=None, inplace=False)[source]¶
Update the lower site index labels to a new string specifier.
- class quimb.tensor.tensor_2d.TensorNetwork2DFlat(ts=(), *, virtual=False, check_collisions=True)[source]¶
Bases:
TensorNetwork2D
Mixin class for a 2D 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', '_Lx', '_Ly')¶
- bond(coo1, coo2)[source]¶
Get the name of the index defining the bond between sites at
coo1
andcoo2
.
- bond_size(coo1, coo2)[source]¶
Return the (combined) size of the bond(s) between sites at
coo1
andcoo2
.
- expand_bond_dimension(new_bond_dim, inplace=True, bra=None, rand_strength=0.0)[source]¶
Increase the bond dimension of this flat, 2D, tensor network, padding the tensor data with either zeros or random entries.
- Parameters:
new_bond_dim (int) – The new dimension. If smaller or equal to the current bond dimension nothing will happend.
inplace (bool, optional) – Whether to expand in place (the default), or return a new TN.
bra (TensorNetwork2DFlat, optional) – Expand this TN with the same data also, assuming it to be the conjugate, bra, TN.
rand_strength (float, optional) – If greater than zero, pad the data arrays with gaussian noise of this strength.
- Returns:
tn
- Return type:
- compress(max_bond=None, cutoff=1e-10, equalize_norms=False, row_sweep='right', col_sweep='up', **compress_opts)[source]¶
Compress all bonds in this flat 2D tensor network.
- Parameters:
max_bond (int, optional) – The maximum boundary dimension, AKA ‘chi’. The default of
None
means truncation is left purely tocutoff
and is not recommended in 2D.cutoff (float, optional) – Cut-off value to used to truncate singular values in the boundary contraction.
compress_opts (None or dict, optional) – Supplied to
compress_between()
.
- class quimb.tensor.tensor_2d.PEPS(arrays, *, shape='urdlp', tags=None, site_ind_id='k{},{}', site_tag_id='I{},{}', x_tag_id='X{}', y_tag_id='Y{}', **tn_opts)[source]¶
Bases:
TensorNetwork2DVector
,TensorNetwork2DFlat
Projected Entangled Pair States object (2D):
... │ │ │ │ │ │ ●────●────●────●────●────●── ╱│ ╱│ ╱│ ╱│ ╱│ ╱│ │ │ │ │ │ │ ●────●────●────●────●────●── ╱│ ╱│ ╱│ ╱│ ╱│ ╱│ │ │ │ │ │ │ ... ●────●────●────●────●────●── ╱│ ╱│ ╱│ ╱│ ╱│ ╱│ │ │ │ │ │ │ ●────●────●────●────●────●── ╱ ╱ ╱ ╱ ╱ ╱
- Parameters:
arrays (sequence of sequence of array_like) – The core tensor data arrays.
shape (str, optional) – Which order the dimensions of the arrays are stored in, the default
'urdlp'
stands for (‘up’, ‘right’, ‘down’, ‘left’, ‘physical’). 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 row (‘x’) tags.
y_tag_id (str, optional) – String specifier for naming convention of column (‘y’) tags.
- _EXTRA_PROPS = ('_site_tag_id', '_x_tag_id', '_y_tag_id', '_Lx', '_Ly', '_site_ind_id')¶
- tags¶
- _site_ind_id¶
- _site_tag_id¶
- _x_tag_id¶
- _y_tag_id¶
- arrays¶
Get the tuple of raw arrays containing all the tensor network data.
- _Lx¶
- _Ly¶
- cyclicx¶
- cyclicy¶
- ix¶
- tensors = []¶
Get the tuple of tensors in this tensor network.
- classmethod from_fill_fn(fill_fn, Lx, Ly, bond_dim, phys_dim=2, cyclic=False, shape='urdlp', **peps_opts)[source]¶
Create a 2D PEPS from a filling function with signature
fill_fn(shape)
.- Parameters:
Lx (int) – The number of rows.
Ly (int) – The number of columns.
bond_dim (int) – The bond dimension.
phys_dim (int, optional) – The physical index dimension.
cyclic (bool or tuple[bool, bool], optional) – Whether the lattice is cyclic in the x and y directions.
shape (str, optional) – How to layout the indices of the tensors, the default is
(up, right, down, left, phys) == 'urdlbk'
. This is the order of the shape supplied to the filling function.peps_opts – Supplied to
PEPS
.
- Returns:
psi
- Return type:
- classmethod empty(Lx, Ly, bond_dim, phys_dim=2, like='numpy', **peps_opts)[source]¶
Create an empty 2D PEPS.
- Parameters:
- Returns:
psi
- Return type:
See also
- classmethod ones(Lx, Ly, bond_dim, phys_dim=2, like='numpy', **peps_opts)[source]¶
Create a 2D PEPS whose tensors are filled with ones.
- Parameters:
- Returns:
psi
- Return type:
See also
- classmethod rand(Lx, Ly, bond_dim, phys_dim=2, dist='normal', loc=0.0, dtype='float64', seed=None, **peps_opts)[source]¶
Create a random (un-normalized) PEPS.
- Parameters:
Lx (int) – The number of rows.
Ly (int) – The number of columns.
bond_dim (int) – The bond dimension.
physical (int, optional) – The physical index dimension.
dist ({'normal', 'uniform', 'rademacher', 'exp'}, optional) – Type of random number to generate, defaults to ‘normal’.
loc (float, optional) – An additive offset to add to the random numbers.
dtype (dtype, optional) – The dtype to create the arrays with, default is real double.
seed (int, optional) – A random seed.
peps_opts – Supplied to
PEPS
.
- Returns:
psi
- Return type:
See also
- class quimb.tensor.tensor_2d.PEPO(arrays, *, shape='urdlbk', tags=None, upper_ind_id='k{},{}', lower_ind_id='b{},{}', site_tag_id='I{},{}', x_tag_id='X{}', y_tag_id='Y{}', **tn_opts)[source]¶
Bases:
TensorNetwork2DOperator
,TensorNetwork2DFlat
Projected Entangled Pair Operator object:
... │╱ │╱ │╱ │╱ │╱ │╱ ●────●────●────●────●────●── ╱│ ╱│ ╱│ ╱│ ╱│ ╱│ │╱ │╱ │╱ │╱ │╱ │╱ ●────●────●────●────●────●── ╱│ ╱│ ╱│ ╱│ ╱│ ╱│ │╱ │╱ │╱ │╱ │╱ │╱ ... ●────●────●────●────●────●── ╱│ ╱│ ╱│ ╱│ ╱│ ╱│ │╱ │╱ │╱ │╱ │╱ │╱ ●────●────●────●────●────●── ╱ ╱ ╱ ╱ ╱ ╱
- Parameters:
arrays (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
'urdlbk'
stands for (‘up’, ‘right’, ‘down’, ‘left’, ‘bra’, ‘ket’). 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.
upper_ind_id (str, optional) – String specifier for naming convention of upper site indices.
lower_ind_id (str, optional) – String specifier for naming convention of lower 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 row (‘x’) tags.
y_tag_id (str, optional) – String specifier for naming convention of column (‘y’) tags.
- _EXTRA_PROPS = ('_site_tag_id', '_x_tag_id', '_y_tag_id', '_Lx', '_Ly', '_upper_ind_id', '_lower_ind_id')¶
- tags¶
- _upper_ind_id¶
- _lower_ind_id¶
- _site_tag_id¶
- _x_tag_id¶
- _y_tag_id¶
- arrays¶
Get the tuple of raw arrays containing all the tensor network data.
- _Lx¶
- _Ly¶
- ix¶
- tensors = []¶
Get the tuple of tensors in this tensor network.
- cyclicx¶
- cyclicy¶
- classmethod from_fill_fn(fill_fn, Lx, Ly, bond_dim, phys_dim=2, cyclic=False, shape='urdlbk', **pepo_opts)[source]¶
Create a PEPO and fill the tensor entries with a supplied function matching signature
fill_fn(shape) -> array
.- Parameters:
fill_fn (callable) – A function that takes a shape tuple and returns a data array.
Lx (int) – The number of rows.
Ly (int) – The number of columns.
bond_dim (int) – The bond dimension.
phys_dim (int, optional) – The physical indices dimension.
cyclic (bool or tuple[bool, bool], optional) – Whether the lattice is cyclic in the x and y directions.
shape (str, optional) – How to layout the indices of the tensors, the default is
(up, right, down, left bra, ket) == 'urdlbk'
.pepo_opts – Supplied to
PEPO
.
- classmethod rand(Lx, Ly, bond_dim, phys_dim=2, herm=False, dist='normal', loc=0.0, dtype='float64', seed=None, **pepo_opts)[source]¶
Create a random PEPO.
- Parameters:
Lx (int) – The number of rows.
Ly (int) – The number of columns.
bond_dim (int) – The bond dimension.
physical (int, optional) – The physical index dimension.
herm (bool, optional) – Whether to symmetrize the tensors across the physical bonds to make the overall operator hermitian.
dtype (dtype, optional) – The dtype to create the arrays with, default is real double.
seed (int, optional) – A random seed.
pepo_opts – Supplied to
PEPO
.
- Returns:
X
- Return type:
- rand_herm¶
- classmethod zeros(Lx, Ly, bond_dim, phys_dim=2, dtype='float64', backend='numpy', **pepo_opts)[source]¶
Create a PEPO with all zero entries.
- Parameters:
Lx (int) – The number of rows.
Ly (int) – The number of columns.
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.
backend (str, optional) – Which backend to use, default is
'numpy'
.pepo_opts – Supplied to
PEPO
.
- apply(other, compress=False, **compress_opts)[source]¶
Act with this PEPO on
other
, returning a new TN likeother
with the same outer indices.- Parameters:
other (PEPS) – The TN to act on.
compress (bool, optional) – Whether to compress the resulting TN.
compress_opts – Supplied to
compress()
.
- Return type:
- quimb.tensor.tensor_2d.show_2d(tn_2d, show_lower=False, show_upper=False)[source]¶
Base function for printing a unicode schematic of flat 2D TNs.
- quimb.tensor.tensor_2d.calc_plaquette_sizes(coo_groups, autogroup=True)[source]¶
Find a sequence of plaquette blocksizes that will cover all the terms (coordinate pairs) in
pairs
.- Parameters:
coo_groups (sequence of tuple[tuple[int]] or tuple[int]) – The sequence of 2D coordinates pairs describing terms. Each should either be a single 2D coordinate or a sequence of 2D coordinates.
autogroup (bool, optional) – Whether to return the minimal sequence of blocksizes that will cover all terms or merge them into a single
((x_bsz, y_bsz),)
.
- Returns:
bszs – Pairs of blocksizes.
- Return type:
Examples
Some nearest neighbour interactions:
>>> H2 = {None: qu.ham_heis(2)} >>> ham = qtn.LocalHam2D(10, 10, H2) >>> calc_plaquette_sizes(ham.terms.keys()) ((1, 2), (2, 1))
>>> calc_plaquette_sizes(ham.terms.keys(), autogroup=False) ((2, 2),)
If we add any next nearest neighbour interaction then we are going to need the (2, 2) blocksize in any case:
>>> H2[(1, 1), (2, 2)] = 0.5 * qu.ham_heis(2) >>> ham = qtn.LocalHam2D(10, 10, H2) >>> calc_plaquette_sizes(ham.terms.keys()) ((2, 2),)
If we add longer range interactions (non-diagonal next nearest) we again can benefit from multiple plaquette blocksizes:
>>> H2[(1, 1), (1, 3)] = 0.25 * qu.ham_heis(2) >>> H2[(1, 1), (3, 1)] = 0.25 * qu.ham_heis(2) >>> ham = qtn.LocalHam2D(10, 10, H2) >>> calc_plaquette_sizes(ham.terms.keys()) ((1, 3), (2, 2), (3, 1))
Or choose the plaquette blocksize that covers all terms:
>>> calc_plaquette_sizes(ham.terms.keys(), autogroup=False) ((3, 3),)
- quimb.tensor.tensor_2d.plaquette_to_sites(p)[source]¶
Turn a plaquette
((i0, j0), (di, dj))
into the sites it contains.Examples
>>> plaquette_to_sites([(3, 4), (2, 2)]) ((3, 4), (3, 5), (4, 4), (4, 5))
- quimb.tensor.tensor_2d.calc_plaquette_map(plaquettes)[source]¶
Generate a dictionary of all the coordinate pairs in
plaquettes
mapped to the ‘best’ (smallest) rectangular plaquette that contains them.Examples
Consider 4 sites, with one 2x2 plaquette and two vertical (2x1) and horizontal (1x2) plaquettes each:
>>> plaquettes = [ ... # 2x2 plaquette covering all sites ... ((0, 0), (2, 2)), ... # horizontal plaquettes ... ((0, 0), (1, 2)), ... ((1, 0), (1, 2)), ... # vertical plaquettes ... ((0, 0), (2, 1)), ... ((0, 1), (2, 1)), ... ]
>>> calc_plaquette_map(plaquettes) {((0, 0), (0, 1)): ((0, 0), (1, 2)), ((0, 0), (1, 0)): ((0, 0), (2, 1)), ((0, 0), (1, 1)): ((0, 0), (2, 2)), ((0, 1), (1, 0)): ((0, 0), (2, 2)), ((0, 1), (1, 1)): ((0, 1), (2, 1)), ((1, 0), (1, 1)): ((1, 0), (1, 2))}
Now every of the size coordinate pairs is mapped to one of the plaquettes, but to the smallest one that contains it. So the 2x2 plaquette (specified by
((0, 0), (2, 2))
) would only used for diagonal terms here.
- quimb.tensor.tensor_2d.gen_long_range_path(ij_a, ij_b, sequence=None)[source]¶
Generate a string of coordinates, in order, from
ij_a
toij_b
.- Parameters:
sequence (None, iterable of {'v', 'h'}, or 'random', optional) – What order to cycle through and try and perform moves in, ‘v’, ‘h’ standing for move vertically and horizontally respectively. The default is
('v', 'h')
.
- Returns:
The path, each element is a single coordinate.
- Return type:
- quimb.tensor.tensor_2d.gen_long_range_swap_path(ij_a, ij_b, sequence=None)[source]¶
Generate the coordinates or a series of swaps that would bring
ij_a
andij_b
together.- Parameters:
sequence (None, it of {'av', 'bv', 'ah', 'bh'}, or 'random', optional) – What order to cycle through and try and perform moves in, ‘av’, ‘bv’, ‘ah’, ‘bh’ standing for move ‘a’ vertically, ‘b’ vertically, ‘a’ horizontally’, and ‘b’ horizontally respectively. The default is
('av', 'bv', 'ah', 'bh')
.
- Returns:
The path, each element is two coordinates to swap.
- Return type: