quimb.experimental.merabuilder¶
Tools for constructing MERA for arbitrary geometry.
TODO:
- [ ] 2D, 3D MERA classes
- [ ] general strategies for arbitrary geometries
- [ ] layer_tag? and hanling of other attributes
- [ ] handle dangling case
- [ ] invariant generators?
DONE::
- [x] layer_gate methods for arbitrary geometry
- [x] 1D: generic way to handle finite and open boundary conditions
- [x] hook into other arbgeom infrastructure for computing rdms etc
Classes¶
A class for building generic 'isometric' or MERA like tensor network |
|
Replacement class for |
Functions¶
|
Given |
|
Return a randomly constructed tree tensor network. |
Module Contents¶
- class quimb.experimental.merabuilder.TensorNetworkGenIso(ts=(), *, virtual=False, check_collisions=True)¶
Bases:
quimb.tensor.tensor_arbgeom.TensorNetworkGenVector
A class for building generic ‘isometric’ or MERA like tensor network states with arbitrary geometry. After supplying the underyling sites of the problem - which can be an arbitrary sequence of hashable objects - one places either unitaries, isometries or tree tensors layered above groups of sites. The isometric and tree tensors effectively coarse grain blocks into a single new site, and the unitaries generally ‘disentangle’ between blocks.
- _EXTRA_PROPS = ('_site_tag_id', '_sites', '_site_ind_id', '_layer_ind_id')¶
- classmethod empty(sites, phys_dim=2, site_tag_id='I{}', site_ind_id='k{}', layer_ind_id='l{}')¶
- property layer_ind_id¶
- layer_ind(site)¶
- layer_gate_raw(G, where, iso=True, new_sites=None, tags=None, all_site_tags=None)¶
Build out this MERA by placing either a new unitary, isometry or tree tensor, given by
G
, at the sites given bywhere
. This handles propagating the lightcone of tags and marking the correct indices of theIsoTensor
asleft_inds
.- Parameters:
G (array_like) – The raw array to place at the sites. Its shape determines whether it is a unitary or isometry/tree. It should have
k + len(where)
dimensions. For a unitaryk == len(where)
. If it is an isometry/tree,k
will generally be1
, or0
to ‘cap’ the MERA. The rightmost indices are those attached to the current open layer indices.where (sequence of hashable) – The sites to layer the tensor above.
iso (bool, optional) – Whether to declare the tensor as an unitary/isometry by marking the left indices. If
iso = False
(a ‘tree’ tensor) then one should havek <= 1
. Once you have such a ‘tree’ tensor you cannot place isometries or unitaries above it. It will also have the lightcone tags of every site. Technically one could place ‘PEPS’ style tensor withiso = False
andk > 1
but some methods might break.new_sites (sequence of hashable, optional) – Which sites to make new open sites. If not given, defaults to the first
k
sites inwhere
.tags (sequence of str, optional) – Custom tags to add to the new tensor, in addition to the automatically generated site tags.
all_site_tags (sequence of str, optional) – For performance, supply all site tags to avoid recomputing them.
- layer_gate_fill_fn(fill_fn, operation, where, max_bond, new_sites=None, tags=None, all_site_tags=None)¶
Build out this MERA by placing either a new unitary, isometry or tree tensor at sites
where
, generating the data array usingfill_fn
and maximum bond dimensionmax_bond
.- Parameters:
fill_fn (callable) – A function with signature
fill_fn(shape) -> array_like
.operation ({"iso", "uni", "cap", "tree", "treecap"}) – The type of tensor to place.
where (sequence of hashable) – The sites to layer the tensor above.
max_bond (int) – The maximum bond dimension of the tensor. This only applies for isometries and trees and when the product of the lower dimensions is greater than
max_bond
.new_sites (sequence of hashable, optional) – Which sites to make new open sites. If not given, defaults to the first
k
sites inwhere
.tags (sequence of str, optional) – Custom tags to add to the new tensor, in addition to the automatically generated site tags.
all_site_tags (sequence of str, optional) – For performance, supply all site tags to avoid recomputing them.
See also
- partial_trace(keep, optimize='auto-hq', rehearse=False, preserve_tensor=False, **contract_opts)¶
Partial trace out all sites except those in
keep
, making use of the lightcone structure of the MERA.- Parameters:
keep (sequence of hashable) – The sites to keep.
optimize (str or PathOptimzer, optional) – The contraction ordering strategy to use.
rehearse ({False, "tn", "tree"}, optional) –
Whether to rehearse the contraction rather than actually performing it. If:
False
: perform the contraction and return the reduced density matrix,”tn”: just the lightcone tensor network is returned,
”tree”: just the contraction tree that will be used is returned.
contract_opts – Additional options to pass to
tensor_contract()
.
- Returns:
The reduced density matrix on sites
keep
.- Return type:
array_like
- local_expectation(G, where, optimize='auto-hq', rehearse=False, **contract_opts)¶
Compute the expectation value of a local operator
G
at siteswhere
. This is done by contracting the lightcone tensor network to form the reduced density matrix, before taking the trace withG
.- Parameters:
G (array_like) – The local operator to compute the expectation value of.
where (sequence of hashable) – The sites to compute the expectation value at.
optimize (str or PathOptimzer, optional) – The contraction ordering strategy to use.
rehearse ({False, "tn", "tree"}, optional) – Whether to rehearse the contraction rather than actually performing it. See
partial_trace()
for details.contract_opts – Additional options to pass to
tensor_contract()
.
- Returns:
The expectation value of
G
at siteswhere
.- Return type:
See also
- compute_local_expectation(terms, optimize='auto-hq', return_all=False, rehearse=False, executor=None, progbar=False, **contract_opts)¶
Compute the expectation value of a collection of local operators
terms
at siteswhere
. This is done by contracting the lightcone tensor network to form the reduced density matrices, before taking the trace with eachG
interms
.- Parameters:
terms (dict[tuple[hashable], array_like]) – The local operators to compute the expectation value of, keyed by the sites they act on.
optimize (str or PathOptimzer, optional) – The contraction ordering strategy to use.
return_all (bool, optional) – Whether to return all the expectation values, or just the sum.
rehearse ({False, "tn", "tree"}, optional) – Whether to rehearse the contraction rather than actually performing it. See
partial_trace()
for details.executor (Executor, optional) – The executor to use for parallelism.
progbar (bool, optional) – Whether to show a progress bar.
contract_opts – Additional options to pass to
tensor_contract()
.
- expand_bond_dimension(new_bond_dim, rand_strength=0.0, inds_to_expand=None, inplace=False)¶
Expand the maxmimum bond dimension of this isometric tensor network to
new_bond_dim
. Unlikeexpand_bond_dimension()
this proceeds from the physical indices upwards, and only increases a bonds size ifnew_bond_dim
is larger than product of the lower indices dimensions.- Parameters:
new_bond_dim (int) – The new maximum bond dimension to expand to.
rand_strength (float, optional) – The strength of random noise to add to the new array entries, if any.
inds_to_expand (sequence of str, optional) – The indices to expand, if not all.
inplace (bool, optional) – Whether to expand this tensor network in place, or return a new one.
- Return type:
- expand_bond_dimension_¶
- quimb.experimental.merabuilder.calc_1d_unis_isos(sites, block_size, cyclic, group_from_right)¶
Given
sites
, assumed to be in a 1D order, though not neccessarily contiguous, calculate unitary and isometry groupings:│ │ <- new grouped site ┐ ┌─────┐ ┌─────┐ ┌ │ │ ISO │ │ ISO │ │ ┘ └─────┘ └─────┘ └ │ │..│..│ │..│..│ │ ┌───┐ │ ┌───┐ │ ┌───┐ │UNI│ │ │UNI│ │ │UNI│ └───┘ │ └───┘ │ └───┘ │ │ ... │ │ ... │ │ ^^^^^^^ <- isometry groupings of size, block_size ^^^^^ ^^^^^ <- unitary groupings of size 2
- Parameters:
sites (sequence of hashable) – The sites to apply a layer to.
block_size (int) – How many sites to group together per isometry block. Note that currently the unitaries will only ever act on blocks of size 2 across isometry block boundaries.
cyclic (bool) – Whether to apply disentangler / unitaries across the boundary. The isometries will never be applied across the boundary, but since they always form a tree such a bipartition is natural.
group_from_right (bool) – Wether to group the sites starting from the left or right. This only matters if
block_size
does not divide the number of sites. Alternating between left and right more evenly tiles the unitaries and isometries, especially at lower layers.
- Returns:
unis (list[tuple]) – The unitary groupings.
isos (list[tuple]) – The isometry groupings.
- class quimb.experimental.merabuilder.MERA(*args, **kwargs)¶
Bases:
quimb.tensor.tensor_1d.TensorNetwork1DVector
,TensorNetworkGenIso
Replacement class for
MERA
which uses the new infrastructure and thus has methods likecompute_local_expectation
.- _EXTRA_PROPS¶
- _CONTRACT_STRUCTURED = False¶
- _num_layers = None¶
- classmethod from_fill_fn(fill_fn, L, D, phys_dim=2, block_size=2, cyclic=True, uni_fill_fn=None, iso_fill_fn=None, cap_fill_fn=None, **kwargs)¶
Create a 1D MERA using
fill_fn(shape) -> array_like
to fill the tensors.- Parameters:
fill_fn (callable) – A function which takes a shape and returns an array_like of that shape. You can override this specfically for the unitaries, isometries and cap tensors using the kwargs
uni_fill_fn
,iso_fill_fn
andcap_fill_fn
.L (int) – The number of sites.
D (int) – The maximum bond dimension.
phys_dim (int, optional) – The dimension of the physical indices.
block_size (int, optional) – The size of the isometry blocks. Binary MERA is the default, ternary MERA is
block_size=3
.cyclic (bool, optional) – Whether to apply disentangler / unitaries across the boundary. The isometries will never be applied across the boundary, but since they always form a tree such a bipartition is natural.
uni_fill_fn (callable, optional) – A function which takes a shape and returns an array_like of that shape. This is used to fill the unitary tensors. If
None
thenfill_fn
is used.iso_fill_fn (callable, optional) – A function which takes a shape and returns an array_like of that shape. This is used to fill the isometry tensors. If
None
thenfill_fn
is used.cap_fill_fn (callable, optional) – A function which takes a shape and returns an array_like of that shape. This is used to fill the cap tensors. If
None
thenfill_fn
is used.kwargs – Supplied to
TensorNetworkGenIso.__init__
.
- classmethod rand(L, D, seed=None, block_size=2, phys_dim=2, cyclic=True, isometrize_method='svd', **kwargs)¶
Return a random (optionally isometrized) MERA.
- Parameters:
L (int) – The number of sites.
D (int) – The maximum bond dimension.
seed (int, optional) – A random seed.
block_size (int, optional) – The size of the isometry blocks. Binary MERA is the default, ternary MERA is
block_size=3
.phys_dim (int, optional) – The dimension of the physical indices.
cyclic (bool, optional) – Whether to apply disentangler / unitaries across the boundary. The isometries will never be applied across the boundary, but since they always form a tree such a bipartition is natural.
isometrize_method (str or None, optional) – If given, the method to use to isometrize the MERA. If
None
then the MERA is not isometrized.
- property num_layers¶
- quimb.experimental.merabuilder.TTN_randtree_rand(sites, D, phys_dim=2, group_size=2, iso=False, seed=None, **kwargs)¶
Return a randomly constructed tree tensor network.
- Parameters:
sites (list of hashable) – The sites of the tensor network.
D (int) – The maximum bond dimension.
phys_dim (int, optional) – The dimension of the physical indices.
group_size (int, optional) – How many sites to group together in each tensor.
iso (bool, optional) – Whether to build the tree with an isometric flow towards the top.
seed (int, optional) – A random seed.
kwargs – Supplied to
TensorNetworkGenIso.empty
.
- Returns:
ttn – The tree tensor network.
- Return type: