quimb.tensor.tensor_arbgeom_tebd

Tools for performing TEBD like algorithms on arbitrary lattices.

Module Contents

Classes

LocalHamGen

Representation of a local hamiltonian defined on a general graph. This

TEBDGen

Generic class for performing time evolving block decimation on an

SimpleUpdateGen

Simple update for arbitrary geometry hamiltonians.

class quimb.tensor.tensor_arbgeom_tebd.LocalHamGen(H2, H1=None)[source]

Representation of a local hamiltonian defined on a general graph. This combines all two site and one site terms into a single interaction per lattice pair, and caches operations on the terms such as getting their exponential. The sites (nodes) should be hashable and comparable.

Parameters:
  • H2 (dict[tuple[node], array_like]) – The interaction terms, with each key being an tuple of nodes defining an edge and each value the local hamilotonian term for those two nodes.

  • H1 (array_like or dict[node, array_like], optional) – The one site term(s). If a single array is given, assume to be the default onsite term for all terms. If a dict is supplied, the keys should represent specific coordinates like (i, j) with the values the array representing the local term for that site. A default term for all remaining sites can still be supplied with the key None.

terms

The total effective local term for each interaction (with single site terms appropriately absorbed). Each key is a pair of coordinates site_a, site_b with site_a < site_b.

Type:

dict[tuple, array_like]

property nsites

The number of sites in the system.

graph[source]
items()[source]

Iterate over all terms in the hamiltonian. This is mostly for convenient compatibility with compute_local_expectation.

_convert_from_qarray_cached(x)[source]
_flip_cached(x)[source]
_add_cached(x, y)[source]
_div_cached(x, y)[source]
_op_id_cached(x)[source]
_id_op_cached(x)[source]
_expm_cached(x, y)[source]
get_gate(where)[source]

Get the local term for pair where, cached.

get_gate_expm(where, x)[source]

Get the local term for pair where, matrix exponentiated by x, and cached.

apply_to_arrays(fn)[source]

Apply the function fn to all the arrays representing terms.

_nx_color_ordering(strategy='smallest_first', interchange=True)[source]

Generate a term ordering based on a coloring on the line graph.

get_auto_ordering(order='sort', **kwargs)[source]

Get an ordering of the terms to use with TEBD, for example. The default is to sort the coordinates then greedily group them into commuting sets.

Parameters:

order ({'sort', None, 'random', str}) –

How to order the terms before greedily grouping them into commuting (non-coordinate overlapping) sets:

  • 'sort' will sort the coordinate pairs first.

  • None will use the current order of terms which should match the order they were supplied to this LocalHam2D instance.

  • 'random' will randomly shuffle the coordinate pairs before grouping them - not the same as returning a completely random order.

  • 'random-ungrouped' will randomly shuffle the coordinate pairs but not group them at all with respect to commutation.

Any other option will be passed as a strategy to networkx.coloring.greedy_color to generate the ordering.

Returns:

Sequence of coordinate pairs.

Return type:

list[tuple[node]]

__repr__()[source]

Return repr(self).

draw(ordering='sort', show_norm=True, figsize=None, fontsize=8, legend=True, ax=None, **kwargs)[source]

Plot this Hamiltonian as a network.

Parameters:
  • ordering ({'sort', None, 'random'}, optional) – An ordering of the termns, or an argument to be supplied to quimb.tensor.tensor_arbgeom_tebd.LocalHamGen.get_auto_ordering() to generate this automatically.

  • show_norm (bool, optional) – Show the norm of each term as edge labels.

  • figsize (None or tuple[int], optional) – Size of the figure, defaults to size of Hamiltonian.

  • fontsize (int, optional) – Font size for norm labels.

  • legend (bool, optional) – Whether to show the legend of which terms are in which group.

  • ax (None or matplotlib.Axes, optional) – Add to a existing set of axes.

class quimb.tensor.tensor_arbgeom_tebd.TEBDGen(psi0, ham, tau=0.01, D=None, imag=True, gate_opts=None, ordering=None, second_order_reflect=False, compute_energy_every=None, compute_energy_final=True, compute_energy_opts=None, compute_energy_fn=None, compute_energy_per_site=False, callback=None, keep_best=False, progbar=True)[source]

Generic class for performing time evolving block decimation on an arbitrary graph, i.e. applying the exponential of a Hamiltonian using a product formula that involves applying local exponentiated gates only.

property state

Return a copy of the current state.

property n

The number of sweeps performed.

property D

The maximum bond dimension.

property energy

Return the energy of current state, computing it only if necessary.

sweep(tau)[source]

Perform a full sweep of gates at every pair.

\[\psi \rightarrow \prod_{\{ij\}} \exp(-\tau H_{ij}) \psi\]
_update_progbar(pbar)[source]
evolve(steps, tau=None, progbar=None)[source]

Evolve the state with the local Hamiltonian for steps steps with time step tau.

_check_energy()[source]

Logic for maybe computing the energy if needed.

get_state()[source]

The default method for retrieving the current state - simply a copy. Subclasses can override this to perform additional transformations.

set_state(psi)[source]

The default method for setting the current state - simply a copy. Subclasses can override this to perform additional transformations.

presweep(i)[source]

Perform any computations required before the sweep (and energy computation). For the basic TEBD this is nothing.

gate(U, where)[source]

Perform single gate U at coordinate pair where. This is the the most common method to override.

compute_energy()[source]

Compute and return the energy of the current state. Subclasses can override this with a custom method to compute the energy.

__repr__()[source]

Return repr(self).

class quimb.tensor.tensor_arbgeom_tebd.SimpleUpdateGen(psi0, ham, tau=0.01, D=None, imag=True, gate_opts=None, ordering=None, second_order_reflect=False, compute_energy_every=None, compute_energy_final=True, compute_energy_opts=None, compute_energy_fn=None, compute_energy_per_site=False, callback=None, keep_best=False, progbar=True)[source]

Bases: TEBDGen

Simple update for arbitrary geometry hamiltonians.

gate(U, where)[source]

Perform single gate U at coordinate pair where. This is the the most common method to override.

compute_energy()[source]

Compute and return the energy of the current state. Subclasses can override this with a custom method to compute the energy.

get_state(absorb_gauges=True)[source]

The default method for retrieving the current state - simply a copy. Subclasses can override this to perform additional transformations.

set_state(psi, gauges=None)[source]

The default method for setting the current state - simply a copy. Subclasses can override this to perform additional transformations.