quimb.experimental.operatorbuilder.operatorbuilder

Tools for defining and constructing sparse operators with:

  • arbitrary geometries,

  • numba acceleration,

  • support for reduced bases,

  • efficient parallelization,

and optionally producing:

  • sparse matrix form

  • matrix product operators,

  • dict of local gates form

  • VMC ‘coupled configs’ form

Currently only supports composing operators which are sums of products of diagonal or anti-diagonal real dimension 2 operators.

TODO:

- [ ] fix sparse matrix being built in opposite direction
- [ ] product of operators generator (e.g. for PEPS DMRG)
- [ ] complex and single precision support (lower priority)
- [ ] support for non-diagonal and qudit operators (lower priority)

DONE:

- [x] use compact bitbasis
- [x] design interface for HilbertSpace / OperatorBuilder interaction
- [x] automatic symbolic jordan wigner transformation
- [x] numba accelerated coupled config
- [x] general definition and automatic 'symbolic' jordan wigner
- [x] multithreaded sparse matrix construction
- [x] LocalHam generator (e.g. for simple update, normal PEPS algs)
- [x] automatic MPO generator

Module Contents

Classes

HilbertSpace

Take a set of 'sites' (any sequence of sortable, hashable objects), and

SparseOperatorBuilder

Object for building operators with sparse structure. Specifically,

Functions

get_local_size(n, rank, world_size)

Given global size n, and a rank in [0, world_size), return the size of

get_local_range(n, rank, world_size)

Given global size n, and a rank in [0, world_size), return the range of

parse_edges_to_unique(edges)

Given a list of edges, return a sorted list of unique sites and edges.

get_mat(op[, dtype])

simplify_single_site_ops(coeff, ops)

Simplify a sequence of operators acting on the same site.

get_nth_bit(val, n)

Get the nth bit of val.

flip_nth_bit(val, n)

Flip the nth bit of val.

comb(n, k)

Compute the binomial coefficient n choose k.

get_all_equal_weight_bits(n, k[, dtype])

Get an array of all 'bits' (integers), with n bits, and k of them set.

bit_to_rank(b, n, k)

Given a bitstring b, return the rank of the bitstring in the

rank_to_bit(r, n, k)

Given a rank r, return the bitstring of length n with k bits set

_recursively_fill_flatconfigs(flatconfigs, n, k, c, r)

get_all_equal_weight_flatconfigs(n, k)

Get every flat configuration of length n with k bits set.

flatconfig_to_bit(flatconfig)

Given a flat configuration, return the corresponding bitstring.

flatconfig_to_rank(flatconfig, n, k)

Given a flat configuration flatconfig, return the rank of the

rank_to_flatconfig(r, n, k)

Given a rank r, return the flat configuration of length n

product_of_bits(b1, b2, n2[, dtype])

Get the outer product of two bit arrays.

get_number_bitbasis(*nk_pairs[, dtype])

Create a bit basis with number conservation.

build_bitmap(configs)

Build of map of bits to linear indices, suitable for use with numba.

build_coupling_numba(term_store, site_to_reg)

Create a sparse nested dictionary of how each term couples each

build_coupling(term_store, site_to_reg)

coupled_flatconfigs_numba(flatconfig, coupling_map)

Get the coupled flat configurations for a given flat configuration

coupled_bits_numba(bi, coupling_map)

_build_coo_numba_core(bits, coupling_map[, bitmap])

build_coo_numba(bits, coupling_map[, parallel])

Build an operator in COO form, using the basis bits and the

fermi_hubbard_from_edges(edges[, t, U, mu])

fermi_hubbard_spinless_from_edges(edges[, t, mu])

heisenberg_from_edges(edges[, j, b, hilbert_space])

Create a Heisenberg Hamiltonian on the graph defined by edges.

Attributes

quimb.experimental.operatorbuilder.operatorbuilder.get_local_size(n, rank, world_size)

Given global size n, and a rank in [0, world_size), return the size of the portion assigned to this rank.

quimb.experimental.operatorbuilder.operatorbuilder.get_local_range(n, rank, world_size)

Given global size n, and a rank in [0, world_size), return the range of indices assigned to this rank.

quimb.experimental.operatorbuilder.operatorbuilder.parse_edges_to_unique(edges)

Given a list of edges, return a sorted list of unique sites and edges.

Parameters:

edges (Iterable[tuple[hashable, hashable]]]) – The edges to parse.

Returns:

  • sites (list of hashable) – The unique sites in the edges, sorted.

  • edges (list of (hashable, hashable)) – The unique edges, sorted.

class quimb.experimental.operatorbuilder.operatorbuilder.HilbertSpace(sites, order=None)

Take a set of ‘sites’ (any sequence of sortable, hashable objects), and map this into a ‘register’ or linearly indexed range, optionally using a particular ordering. One can then calculate the size of the Hilbert space, including number conserving subspaces, and get a compact ‘bitbasis’ with which to construct sparse operators.

Parameters:
  • sites (int or sequence of hashable objects) – The sites to map into a linear register. If an integer, simply use range(sites).

  • order (callable or sequence of hashable objects, optional) – If provided, use this to order the sites. If a callable, it should be a sorting key. If a sequence, it should be a permutation of the sites, and key=order.index will be used.

property sites

The ordered tuple of all sites in the Hilbert space.

property nsites

The total number of sites in the Hilbert space.

set_ordering(order=None)
classmethod from_edges(edges, order=None)

Construct a HilbertSpace from a set of edges, which are pairs of sites.

site_to_reg(site)

Convert a site to a linear register index.

reg_to_site(reg)

Convert a linear register index back to a site.

has_site(site)

Check if this HilbertSpace contains a given site.

config_to_bit(config)

Encode a ‘configuration’ as a bit.

Parameters:

config (dict[hashable, int]) – A dictionary mapping sites to their occupation number / spin.

Returns:

bit – The bit corresponding to this configuration.

Return type:

int

config_to_flatconfig(config)

Turn a configuration into a flat configuration, assuming the order given by this HilbertSpace.

Parameters:

config (dict[hashable, int]) – A dictionary mapping sites to their occupation number / spin.

Returns:

flatconfig – A flat configuration, with the occupation number of each site in the order given by this HilbertSpace.

Return type:

ndarray[uint8]

bit_to_config(bit)

Decode a bit to a configuration.

Parameters:

bit (int) – The bit to decode.

Returns:

config – A dictionary mapping sites to their occupation number / spin.

Return type:

dict[hashable, int]

flatconfig_to_config(flatconfig)

Turn a flat configuration into a configuration, assuming the order given by this HilbertSpace.

Parameters:

flatconfig (ndarray[uint8]) – A flat configuration, with the occupation number of each site in the order given by this HilbertSpace.

Returns:

config – A dictionary mapping sites to their occupation number / spin.

Return type:

dict[hashable, int]

rand_config(k=None)

Get a random configuration, optionally requiring it has k bits set.

get_size(*k)

Compute the size of this Hilbert space, optionally taking into account number / z-spin conservation.

Parameters:

k (int or tuple of (int, int)) –

If provided, compute the size of number conserving subspace(s):

- If a single int, compute the size of the subspace with
  ``k`` particles / up states: ``comb(nsites, k)``.
- If a tuple of (int, int), compute the size of the subspace
  of the product of spaces where each pair (n, k) corresponds
  to n sites with k particles / up states. The sum of every n
  should equal ``nsites``.

get_bitbasis(*k, dtype=np.int64)

Get a basis for the Hilbert space, in terms of an integer bitarray, optionally taking into account number / z-spin conservation.

Parameters:
  • k (int or tuple of (int, int)) –

    If provided, get the basis for a number conserving subspace(s):

    - If a single int, compute the size of the subspace with
      ``k`` particles / up states: ``comb(nsites, k)``.
    - If a tuple of (int, int), compute the size of the subspace
      of the product of spaces where each pair (n, k) corresponds
      to n sites with k particles / up states. The sum of every n
      should equal ``nsites``.
    

  • dtype (np.dtype, optional) – The dtype of the bitarray, should be an integer type with at least nsites bits.

Returns:

bits – The basis, each integer element being a binary representation of a configuration.

Return type:

numpy.ndarray

__repr__()

Return repr(self).

quimb.experimental.operatorbuilder.operatorbuilder._OPMAP
quimb.experimental.operatorbuilder.operatorbuilder.get_mat(op, dtype=None)
quimb.experimental.operatorbuilder.operatorbuilder.simplify_single_site_ops(coeff, ops)

Simplify a sequence of operators acting on the same site.

Parameters:
  • coeff (float or complex) – The coefficient of the operator sequence.

  • ops (tuple of str) – The operator sequence.

Returns:

  • new_coeff (float or complex) – The simplified coefficient.

  • new_op (str) – The single, simplified operator that the sequence maps to, up to scalar multiplication.

Examples

>>> simplify_single_site_ops(1.0, ('+', 'z', 'z', 'z', 'z', '-'))
(1.0, 'n')
>>> simplify_single_site_ops(1.0, ("x", "y", "z"))
(-1j, 'I')
class quimb.experimental.operatorbuilder.operatorbuilder.SparseOperatorBuilder(terms=(), hilbert_space=None)

Object for building operators with sparse structure. Specifically, a sum of terms, where each term is a product of operators, where each of these local operators acts on a single site and has at most one entry per row.

Parameters:
  • terms (sequence, optional) – The terms to initialize the builder with. add_term is simply called on each of these.

  • hilbert_space (HilbertSpace) – The Hilbert space to build the operator in. If this is not supplied then a minimal Hilbert space will be constructed from the sites used, when required.

property sites_used

A tuple of the sorted coordinates/sites seen so far.

property nsites

The total number of coordinates/sites seen so far.

property terms

A tuple of the simplified terms seen so far.

property nterms

The total number of terms seen so far.

property locality

The locality of the operator, the maximum support of any term.

property hilbert_space

The Hilbert space of the operator. Created from the sites seen so far if not supplied at construction.

property coupling_map
site_to_reg(site)

Get the register / linear index of coordinate site.

reg_to_site(reg)

Get the site of register / linear index reg.

add_term(*coeff_ops)

Add a term to the operator.

Parameters:
  • coeff (float, optional) – The overall coefficient of the term.

  • ops (sequence of tuple[str, hashable]) –

    The operators of the term, together with the sites they act on. Each term should be a pair of (operator, site), where operator can be:

    • 'x', 'y', 'z': Pauli matrices

    • 'sx', 'sy', 'sz': spin operators (i.e. scaled Pauli matrices)

    • '+', '-': creation/annihilation operators

    • 'n', 'sn', or 'h': number, symmetric number (n - 1/2) and hole (1 - n) operators

    And site is a hashable object that represents the site that the operator acts on.

__iadd__(term)
jordan_wigner_transform()

Transform the terms in this operator by pre-prending pauli Z strings to all creation and annihilation operators, and then simplifying the resulting terms.

build_coo_data(*k, parallel=False)

Build the raw data for a sparse matrix in COO format. Optionally in a reduced k basis and in parallel.

Parameters:
  • k (int or tuple of (int, int)) –

    If provided, get the basis for a number conserving subspace(s):

    - If a single int, compute the size of the subspace with
      ``k`` particles / up states: ``comb(nsites, k)``.
    - If a tuple of (int, int), compute the size of the subspace
      of the product of spaces where each pair (n, k) corresponds
      to n sites with k particles / up states. The sum of every n
      should equal ``nsites``.
    

  • parallel (bool, optional) – Whether to build the matrix in parallel (multi-threaded).

Returns:

  • coo (array) – The data entries for the sparse matrix in COO format.

  • cis (array) – The row indices for the sparse matrix in COO format.

  • cjs (array) – The column indices for the sparse matrix in COO format.

  • N (int) – The total number of basis states.

build_sparse_matrix(*k, stype='csr', parallel=False)

Build a sparse matrix in the given format. Optionally in a reduced k basis and in parallel.

Parameters:
  • k (int or tuple of (int, int)) –

    If provided, get the basis for a number conserving subspace(s):

    - If a single int, compute the size of the subspace with
      ``k`` particles / up states: ``comb(nsites, k)``.
    - If a tuple of (int, int), compute the size of the subspace
      of the product of spaces where each pair (n, k) corresponds
      to n sites with k particles / up states. The sum of every n
      should equal ``nsites``.
    

  • parallel (bool, optional) – Whether to build the matrix in parallel (multi-threaded).

Return type:

scipy.sparse matrix

build_dense()

Get the dense (numpy.ndarray) matrix representation of this operator.

build_local_terms()

Get a dictionary of local terms, where each key is a sorted tuple of sites, and each value is the local matrix representation of the operator on those sites. For use with e.g. tensor network algorithms.

Note terms acting on the same sites are summed together and the size of each local matrix is exponential in the locality of that term.

Returns:

Hk – The local terms.

Return type:

dict[tuple[hashable], numpy.ndarray]

config_coupling(config)

Get a list of other configurations coupled to config by this operator, and the corresponding coupling coefficients. This is for use with VMC for example.

Parameters:

config (dict[site, int]) – The configuration to get the coupling for.

Returns:

  • coupled_configs (list[dict[site, int]]) – Each distinct configuration coupled to config.

  • coeffs (list[float]) – The corresponding coupling coefficients.

show(filler='.')

Print an ascii representation of the terms in this operator.

build_state_machine_greedy()
draw_state_machine(method='greedy', figsize='auto')
build_mpo(method='greedy', **mpo_opts)
__repr__()

Return repr(self).

quimb.experimental.operatorbuilder.operatorbuilder.get_nth_bit(val, n)

Get the nth bit of val.

Examples

>>> get_nth_bit(0b101, 1)
0
quimb.experimental.operatorbuilder.operatorbuilder.flip_nth_bit(val, n)

Flip the nth bit of val.

Examples

>>> bin(flip_nth_bit(0b101, 1))
0b111
quimb.experimental.operatorbuilder.operatorbuilder.comb(n, k)

Compute the binomial coefficient n choose k.

quimb.experimental.operatorbuilder.operatorbuilder.get_all_equal_weight_bits(n, k, dtype=np.int64)

Get an array of all ‘bits’ (integers), with n bits, and k of them set.

quimb.experimental.operatorbuilder.operatorbuilder.bit_to_rank(b, n, k)

Given a bitstring b, return the rank of the bitstring in the basis of all bitstrings of length n with k bits set. Adapted from https://dlbeer.co.nz/articles/kwbs.html.

quimb.experimental.operatorbuilder.operatorbuilder.rank_to_bit(r, n, k)

Given a rank r, return the bitstring of length n with k bits set that has rank r in the basis of all bitstrings of length n with k bits set. Adapted from https://dlbeer.co.nz/articles/kwbs.html.

quimb.experimental.operatorbuilder.operatorbuilder._recursively_fill_flatconfigs(flatconfigs, n, k, c, r)
quimb.experimental.operatorbuilder.operatorbuilder.get_all_equal_weight_flatconfigs(n, k)

Get every flat configuration of length n with k bits set.

quimb.experimental.operatorbuilder.operatorbuilder.flatconfig_to_bit(flatconfig)

Given a flat configuration, return the corresponding bitstring.

quimb.experimental.operatorbuilder.operatorbuilder.flatconfig_to_rank(flatconfig, n, k)

Given a flat configuration flatconfig, return the rank of the bitstring in the basis of all bitstrings of length n with k bits set. Adapted from https://dlbeer.co.nz/articles/kwbs.html.

quimb.experimental.operatorbuilder.operatorbuilder.rank_to_flatconfig(r, n, k)

Given a rank r, return the flat configuration of length n with k bits set that has rank r in the basis of all bitstrings of length n with k bits set. Adapted from https://dlbeer.co.nz/articles/kwbs.html.

quimb.experimental.operatorbuilder.operatorbuilder.product_of_bits(b1, b2, n2, dtype=np.int64)

Get the outer product of two bit arrays.

Parameters:
  • b1 (array) – The bit arrays to take the outer product of.

  • b2 (array) – The bit arrays to take the outer product of.

  • n2 (int) – The number of bits in b2

quimb.experimental.operatorbuilder.operatorbuilder.get_number_bitbasis(*nk_pairs, dtype=np.int64)

Create a bit basis with number conservation.

Parameters:

nk_pairs (sequence of (int, int)) – Each element is a pair (n, k) where n is the number of bits, and k is the number of bits that are set. The basis will be the product of all supplied pairs.

Returns:

basis – An array of integers, each representing a bit string. The size will be prod(comb(n, k) for n, k in nk_pairs).

Return type:

ndarray

Examples

A single number conserving basis:

>>> for i, b in enumerate(get_number_bitbasis((4, 2))):
>>>    print(f"{i}: {b:0>4b}")
0: 0011
1: 0101
2: 0110
3: 1001
4: 1010
5: 1100

A product of two number conserving bases, e.g. n_up and n_down:

>>> for b in get_number_bitbasis((3, 2), (3, 1)):
>>>     print(f"{b:0>6b}")
011001
011010
011100
101001
101010
101100
110001
110010
110100
quimb.experimental.operatorbuilder.operatorbuilder.build_bitmap(configs)

Build of map of bits to linear indices, suitable for use with numba.

quimb.experimental.operatorbuilder.operatorbuilder.build_coupling_numba(term_store, site_to_reg)

Create a sparse nested dictionary of how each term couples each local site configuration to which other local site configuration, and with what coefficient, suitable for use with numba.

Parameters:
  • term_store (dict[term, coeff]) – The terms of the operator.

  • site_to_reg (callable) – A function that maps a site to a linear register index.

Returns:

coupling_map – A nested numba dictionary of the form {term: {reg: {bit_in: (bit_out, coeff), ...}, ...}, ...}.

Return type:

numba.typed.Dict

quimb.experimental.operatorbuilder.operatorbuilder.build_coupling(term_store, site_to_reg)
quimb.experimental.operatorbuilder.operatorbuilder.coupled_flatconfigs_numba(flatconfig, coupling_map)

Get the coupled flat configurations for a given flat configuration and coupling map.

Parameters:
  • flatconfig (ndarray[uint8]) – The flat configuration to get the coupled configurations for.

  • coupling_map (numba.typed.Dict) – A nested numba dictionary of the form {term: {reg: {bit_in: (bit_out, coeff), ...}, ...}, ...}.

Returns:

  • coupled_flatconfigs (ndarray[uint8]) – A list of coupled flat configurations, each with the corresponding coefficient.

  • coeffs (ndarray[float64]) – The coefficients for each coupled flat configuration.

quimb.experimental.operatorbuilder.operatorbuilder.coupled_bits_numba(bi, coupling_map)
quimb.experimental.operatorbuilder.operatorbuilder._build_coo_numba_core(bits, coupling_map, bitmap=None)
quimb.experimental.operatorbuilder.operatorbuilder.build_coo_numba(bits, coupling_map, parallel=False)
Build an operator in COO form, using the basis bits and the

coupling_map, optionally multithreaded.

Parameters:
  • bits (array) – An array of integers, each representing a bit string.

  • coupling_map (Dict[int, Dict[int, Dict[int, Tuple[int, float]]]]) – A nested numba dictionary of couplings. The outermost key is the term index, the next key is the register, and the innermost key is the bit index. The value is a tuple of (coupled bit index, coupling coefficient).

  • parallel (bool or int, optional) – Whether to parallelize the computation. If an integer is given, it specifies the number of threads to use.

Returns:

  • data (array) – The non-zero elements of the operator.

  • cis (array) – The row indices of the non-zero elements.

  • cjs (array) – The column indices of the non-zero elements.

quimb.experimental.operatorbuilder.operatorbuilder.fermi_hubbard_from_edges(edges, t=1.0, U=1.0, mu=0.0)
quimb.experimental.operatorbuilder.operatorbuilder.fermi_hubbard_spinless_from_edges(edges, t=1.0, mu=0.0)
quimb.experimental.operatorbuilder.operatorbuilder.heisenberg_from_edges(edges, j=1.0, b=0.0, hilbert_space=None)

Create a Heisenberg Hamiltonian on the graph defined by edges.

Parameters:
  • edges (Iterable[tuple[hashable, hashable]]) – The edges, as pairs of hashable ‘sites’, that define the graph. Multiple edges are allowed, and will be treated as a single edge.

  • j (float or tuple[float, float, float], optional) – The Heisenberg exchange coupling constant(s). If a single float is given, it is used for all three terms. If a tuple of three floats is given, they are used for the xx, yy, and zz terms respectively. Note that positive values of j correspond to antiferromagnetic coupling.

  • b (float or tuple[float, float, float], optional) – The magnetic field strength(s). If a single float is given, it is used taken as a z-field. If a tuple of three floats is given, they are used for the x, y, and z fields respectively.

  • hilbert_space (HilbertSpace, optional) – The Hilbert space to use. If not given, one will be constructed automatically from the edges.

Returns:

H

Return type:

SparseOperatorBuilder