quimb.tensor.tensor_arbgeom_tebd¶
Tools for performing TEBD like algorithms on arbitrary lattices.
Classes¶
Representation of a local hamiltonian defined on a general graph. This |
|
Generic class for performing time evolving block decimation on an |
|
Simple update for arbitrary geometry hamiltonians. |
Functions¶
|
Generate an edge coloring for the graph given by |
Module Contents¶
- quimb.tensor.tensor_arbgeom_tebd.edge_coloring(edges, strategy='smallest_last', interchange=True, group=True)[source]¶
Generate an edge coloring for the graph given by
edges
, usingnetworkx.coloring.greedy_color
.- Parameters:
edges (sequence[tuple[hashable, hashable]]) – The edges of the graph.
strategy (str or callable, optional) –
The strategy to use for coloring the edges. Can be:
’largest_first’
’smallest_last’
’random_sequential’
…
interchange (bool, optional) – Whether to use the interchange heuristic. Usually generates better colorings but can be slower.
group (bool, optional) – Whether to group the edges by color or return a flat list.
- 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 keyNone
.
- 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
withsite_a < site_b
.
- _op_cache¶
- terms¶
- sites¶
- _sites_to_covering_terms¶
- default_H1¶
- property nsites¶
- The number of sites in the system.
- items()[source]¶
Iterate over all terms in the hamiltonian. This is mostly for convenient compatibility with
compute_local_expectation
.
- get_gate_expm(where, x)[source]¶
Get the local term for pair
where
, matrix exponentiated byx
, and cached.
- 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 thisLocalHam2D
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:
- 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.
- imag¶
- property state¶
- Return a copy of the current state.
- ham¶
- progbar¶
- callback¶
- tau¶
- last_tau = 0.0¶
- gate_opts¶
- compute_energy_opts¶
- compute_energy_every¶
- compute_energy_final¶
- compute_energy_fn¶
- compute_energy_per_site¶
- second_order_reflect¶
- _n = 0¶
- its = []¶
- taus = []¶
- energies = []¶
- keep_best¶
- best¶
- sweep(tau)[source]¶
Perform a full sweep of gates at every pair.
\[\psi \rightarrow \prod_{\{ij\}} \exp(-\tau H_{ij}) \psi\]
- evolve(steps, tau=None, progbar=None)[source]¶
Evolve the state with the local Hamiltonian for
steps
steps with time steptau
.
- 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.
- 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 pairwhere
. This is the the most common method to override.
- 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 pairwhere
. 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.