quimb.tensor.tensor_builder¶
Build specific tensor networks, including states and operators.
Attributes¶
Classes¶
Simple class to allow |
|
Class for easily building custom spin hamiltonians in MPO or LocalHam1D |
Functions¶
|
Get the delta symbol AKA COPY-tensor as an ndarray. |
|
Generate a random tensor with specified shape and inds. |
|
Generate a random tensor with specified shape and inds, and randomly |
|
Get a random symmetric array, i.e. one that is invariant under |
|
Generate a random symmetric tensor with specified local dimension and |
|
A product state in general tensor network form. |
|
A computational basis state in general tensor network form. |
|
|
|
Compute a dictionary of edge frequencies for a list of strings, |
|
Create a tensor network from a sequence of edges defining a graph, |
|
Create a tensor network from a sequence of edges defining a graph, |
|
Create a tensor network from a sequence of edges defining a graph, |
|
Create a random tensor network with geometry defined from a sequence |
|
Create a random regular tensor network. |
|
Create a random tree tensor network. |
|
|
|
Create a random, possibly hyper, tensor network, with a mix of normal |
|
reate a CP-decomposition structured hyper tensor network from a |
|
Create a CP-decomposition structured hyper tensor network state from a |
|
Construct a CP form hyper tensor network of the sum of matrix strings: |
|
Create a hyper tensor network with a tensor on each bond and a hyper |
|
Convert |
|
A scalar 2D lattice tensor network with tensors filled by a function. |
|
A scalar 2D lattice tensor network initialized with empty tensors. |
|
A scalar 2D lattice tensor network with every element set to |
|
A random scalar 2D lattice tensor network. |
|
Create a random 2D lattice tensor network where every tensor is |
|
Build a 2D 'corner double line' (CDL) tensor network. Each plaquette |
|
|
|
Convert |
|
A scalar 3D lattice tensor network with tensors filled by a function. |
|
A scalar 3D lattice tensor network initialized with empty tensors. |
|
A scalar 2D lattice tensor network with every element set to |
|
A random scalar 3D lattice tensor network. |
|
|
|
|
|
The interaction term for the classical ising model. |
|
The magnetic field term for the classical ising model. |
|
The sqrt factorized interaction term for the classical ising model. |
|
The single effective TN site for the classical ising model. |
Parse the |
|
|
Hyper tensor network representation of the 2D classical ising model |
|
Hyper tensor network representation of the 3D classical ising model |
|
The tensor network representation of the 2D classical ising model |
|
Tensor network representation of the 3D classical ising model |
|
Build a hyper tensor network representation of a classical ising model |
|
Build a regular tensor network representation of a classical ising model |
|
Take a coupling matrix or dict of pairwise strengths in either upper or |
Construct a (triangular) '2D' tensor network representation of the |
|
|
|
|
Make a tensor network from sequence of graph edges that counts the |
|
Encode a clause as a single integer |
|
Get the array representing satisfiability of |
|
Get the tensor representing satisfiability of |
|
Get the set of MPS tensors representing satisfiability of |
|
Get the set of PARAFAC arrays representing satisfiability of |
|
Get the set of PARAFAC tensors representing satisfiability of |
|
Given a list of clauses, create a hyper tensor network, with a single |
|
Parse a DIMACS style 'cnf' file into a list of clauses, and possibly a |
|
Create a hyper tensor network from a '.cnf' or '.wcnf' file - i.e. a |
|
Create a random k-SAT instance. |
|
Create a random k-SAT instance encoded as a hyper tensor network. |
|
Create a tensor network with the same outer indices as |
|
Generate a random matrix product state. |
|
Generate a product state in MatrixProductState form, i,e, |
|
A computational basis state in Matrix Product State form. |
|
Generate the neel state in Matrix Product State form. |
|
Build a matrix product state representation of the COPY tensor. |
|
Build the chi=2 OBC MPS representation of the GHZ state. |
|
Build the chi=2 OBC MPS representation of the W state. |
|
Generate a random computation basis state, like '01101001010'. |
|
The all-zeros MPS state, of given bond-dimension. |
|
A product state for sampling tensor network traces. Seen as a vector it |
|
Generate an identity MPO of size |
|
Return an identity matrix operator with the same physical index and |
|
Generate a zeros MPO of size |
|
Return a zeros matrix product operator with the same physical index and |
|
Return an MPO of bond dimension 1 representing the product of raw |
|
Generate a random matrix product state. |
|
Generate a random hermitian matrix product operator. |
Check if |
|
|
Generate tensor(s) for a spin hamiltonian MPO. |
|
|
|
Ising Hamiltonian in MPO form. |
|
Ising Hamiltonian in |
|
|
|
XY-Hamiltonian in MPO form. |
|
XY-Hamiltonian in |
|
|
|
Heisenberg Hamiltonian in MPO form. |
|
Heisenberg Hamiltonian in |
|
XXZ-Hamiltonian in MPO form. |
|
XXZ-Hamiltonian in |
|
|
|
Hamiltonian of one-dimensional bilinear biquadratic chain in MPO form, |
|
Hamiltonian of one-dimensional bilinear biquadratic chain in LocalHam1D |
|
|
|
The many-body-localized spin hamiltonian in MPO form. |
|
The many-body-localized spin hamiltonian in |
|
Ising Hamiltonian in |
|
Heisenberg Hamiltonian in |
|
Heisenberg Hamiltonian in |
|
Heisenberg Hamiltonian in |
Module Contents¶
- quimb.tensor.tensor_builder.delta_array(shape, dtype='float64')¶
Get the delta symbol AKA COPY-tensor as an ndarray.
- Parameters:
- Return type:
- quimb.tensor.tensor_builder.rand_tensor(shape, inds, tags=None, dtype='float64', dist='normal', scale=1.0, loc=0.0, left_inds=None, **randn_opts)¶
Generate a random tensor with specified shape and inds.
- Parameters:
shape (sequence of int) – Size of each dimension.
inds (sequence of str) – Names of each dimension.
tags (sequence of str) – Labels to tag this tensor with.
dtype ({'float64', 'complex128', 'float32', 'complex64'}, optional) – The underlying data type.
dist ({'normal', 'uniform', 'rademacher', 'exp'}, optional) – Type of random number to generate, defaults to ‘normal’.
scale (float, optional) – A multiplicative factor to scale the random numbers by.
loc (float, optional) – An additive offset to add to the random numbers.
left_inds (sequence of str, optional) – Which, if any, indices to group as ‘left’ indices of an effective matrix. This can be useful, for example, when automatically applying unitary constraints to impose a certain flow on a tensor network but at the atomistic (Tensor) level.
- Return type:
- quimb.tensor.tensor_builder.rand_phased(shape, inds, tags=None, dtype=complex)¶
Generate a random tensor with specified shape and inds, and randomly ‘phased’ (distributed on the unit circle) data, such that
T.H @ T == T.norm()**2 == T.size
.
- quimb.tensor.tensor_builder.rand_symmetric_array(d, ndim, dist='normal', loc=0.0, scale=1.0, seed=None, dtype='float64', fill_fn=None)¶
Get a random symmetric array, i.e. one that is invariant under permutation of its indices. It has (ndim + d - 1) choose (d - 1) unique elements.
- Parameters:
d (int) – The size of each dimension.
ndim (int) – The number of dimensions.
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.
scale (float, optional) – A multiplicative factor to scale the random numbers by.
seed (int, optional) – A random seed.
dtype ({'float64', 'complex128', 'float32', 'complex64'}, optional) – The underlying data type.
fill_fn (callable, optional) – A function that takes a shape (here always ()) and returns an array, rather that using the random number generator.
- Return type:
- quimb.tensor.tensor_builder.rand_tensor_symmetric(d, inds, tags=None, dist='normal', loc=0.0, scale=1.0, seed=None)¶
Generate a random symmetric tensor with specified local dimension and inds.
- Parameters:
d (int) – Size of each dimension.
inds (sequence of str) – Names of each dimension.
tags (sequence of str) – Labels to tag this tensor with.
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.
scale (float, optional) – A multiplicative factor to scale the random numbers by.
seed (int, optional) – A random seed.
- Return type:
- quimb.tensor.tensor_builder.TN_from_sites_product_state(site_map, site_tag_id='I{}', site_ind_id='k{}')¶
A product state in general tensor network form.
- Parameters:
- Return type:
- quimb.tensor.tensor_builder.TN_from_sites_computational_state(site_map, site_tag_id='I{}', site_ind_id='k{}', dtype='float64')¶
A computational basis state in general tensor network form.
- Parameters:
site_map (dict[hashable, str]) – Mapping of site to computational state, which should be one of
('0', '1', '+', '-')
.site_tag_id (str, optional) – Format string for site tag labels.
site_ind_id (str, optional) – Format string for site index labels.
dtype ({'float64', 'complex128', 'float32', 'complex64'}, optional) – The data type to use for the array representation.
- Return type:
- quimb.tensor.tensor_builder.gen_unique_edges(edges)¶
- quimb.tensor.tensor_builder.compute_string_edge_frequencies(strings)¶
Compute a dictionary of edge frequencies for a list of strings, including plaquettes.
- quimb.tensor.tensor_builder.TN_from_edges_and_fill_fn(fill_fn, edges, D, phys_dim=None, site_tag_id='I{}', site_ind_id='k{}')¶
Create a tensor network from a sequence of edges defining a graph, and a ‘fill’ function that maps shapes to data.
- Parameters:
fill_fn (callable) – A function with signature
fill_fn(shape) -> array
, used to fill each tensor.edges (sequence of tuple[hashable, hashable]) – The graph edges, as a sequence of pairs of hashable objects, for example integers, representing the nodes. You can redundantly specify
(u, v)
and(v, u)
and only one edge will be added.D (int) – The bond dimension connecting tensors.
phys_dim (int, optional) – If not
None
, give each tensor a ‘physical’, free index of this size at each node.site_tag_id (str, optional) – String with formatter to tag sites.
site_ind_id (str or (str, str), optional) – String with formatter to tag indices (if
phys_dim
specified). If a single str is supplied, the tensor network will have a single index at each site, representing a vector. If a pair of strings is supplied, the tensor network will have two indices at each site, representing an operator with upper and lower indices.
- Return type:
TensorNetworkGen, TensorNetworkGenVector or TensorNetworkGenOperator
- quimb.tensor.tensor_builder.TN_from_edges_empty(edges, D, phys_dim=None, site_tag_id='I{}', site_ind_id='k{}', dtype='float64')¶
Create a tensor network from a sequence of edges defining a graph, initialized with empty tensors.
- Parameters:
edges (sequence of tuple[hashable, hashable]) – The graph edges, as a sequence of pairs of hashable objects, for example integers, representing the nodes. You can redundantly specify
(u, v)
and(v, u)
and only one edge will be added.D (int) – The bond dimension connecting tensors.
phys_dim (int, optional) – If not
None
, give each tensor a ‘physical’, free index of this size at each node.site_tag_id (str, optional) – String with formatter to tag sites.
site_ind_id (str, optional) – String with formatter to tag indices (if
phys_dim
specified).dtype (str, optional) – The data type of the tensors.
- Return type:
TensorNetworkGen, TensorNetworkGenVector or TensorNetworkGenOperator
- quimb.tensor.tensor_builder.TN_from_edges_with_value(value, edges, D, phys_dim=None, site_tag_id='I{}', site_ind_id='k{}', dtype=None)¶
Create a tensor network from a sequence of edges defining a graph, initialized with a constant value. This uses
numpy.broadcast_to
and therefore essentially no memory.- Parameters:
value (scalar) – The value to fill the tensors with.
edges (sequence of tuple[hashable, hashable]) – The graph edges, as a sequence of pairs of hashable objects, for example integers, representing the nodes. You can redundantly specify
(u, v)
and(v, u)
and only one edge will be added.D (int) – The bond dimension connecting tensors.
phys_dim (int, optional) – If not
None
, give each tensor a ‘physical’, free index of this size at each node.site_tag_id (str, optional) – String with formatter to tag sites.
site_ind_id (str, optional) – String with formatter to tag indices (if
phys_dim
specified).dtype (str, optional) – The data type of the tensors.
- Return type:
TensorNetworkGen, TensorNetworkGenVector or TensorNetworkGenOperator
- quimb.tensor.tensor_builder.TN_from_edges_rand(edges, D, phys_dim=None, seed=None, dtype='float64', site_tag_id='I{}', site_ind_id='k{}', **randn_opts)¶
Create a random tensor network with geometry defined from a sequence of edges defining a graph.
- Parameters:
G (sequence of tuple[node, node]) – The edges defining a graph, each element should be a pair of nodes described by hashable objects.
D (int) – The bond dimension connecting tensors.
phys_dim (int, optional) – If not
None
, give each tensor a ‘physical’, free index of this size to mimic a wavefunction oflen(G)
sites.seed (int, optional) – A random seed.
dtype (str, optional) – The data type of the tensors.
site_tag_id (str, optional) – String with formatter to tag sites.
site_ind_id (str, optional) – String with formatter to tag indices (if
phys_dim
specified).randn_opts – Supplied to
randn()
.
- Return type:
TensorNetworkGen, TensorNetworkGenVector or TensorNetworkGenOperator
- quimb.tensor.tensor_builder.TN_rand_from_edges¶
- quimb.tensor.tensor_builder.TN_rand_reg(n, reg, D, phys_dim=None, seed=None, dtype='float64', site_tag_id='I{}', site_ind_id='k{}', **randn_opts)¶
Create a random regular tensor network.
- Parameters:
n (int) – The number of tensors.
reg (int) – The degree of the tensor network (how many tensors each tensor connects to).
D (int) – The bond dimension connecting tensors.
phys_dim (int, optional) – If not
None
, give each tensor a ‘physical’, free index of this size to mimic a wavefunction ofn
sites.seed (int, optional) – A random seed.
site_tag_id (str, optional) – String with formatter to tag sites.
site_ind_id (str, optional) – String with formatter to tag indices (if
phys_dim
specified).
- Return type:
TensorNetworkGen, TensorNetworkGenVector or TensorNetworkGenOperator
- quimb.tensor.tensor_builder.TN_rand_tree(n, D, phys_dim=None, max_degree=None, seed=None, dtype='float64', site_tag_id='I{}', site_ind_id='k{}', **randn_opts)¶
Create a random tree tensor network.
- Parameters:
n (int) – The number of tensors.
D (int) – The bond dimension connecting tensors.
phys_dim (int, optional) – If not
None
, give each tensor a ‘physical’, free index of this size to mimic a wavefunction ofn
sites.max_degree (int, optional) – The maximum degree of any node in the tree, for example 3 means that each tensor can connect to at most 3 other tensors - creating a binary tree (1 parent and two children).
seed (int, optional) – A random seed.
site_tag_id (str, optional) – String with formatter to tag sites.
site_ind_id (str, optional) – String with formatter to tag indices (if
phys_dim
specified).
- Return type:
TensorNetworkGen, TensorNetworkGenVector or TensorNetworkGenOperator
- quimb.tensor.tensor_builder.TN_from_strings(strings, fill_fn=None, line_dim=2, allow_plaquettes=True, site_tag_id='I{}', random_rewire=False, random_rewire_seed=None, join=False, join_avoid_self_loops=True, normalize=False, contract_sites=True, fuse_multibonds=True, **contract_opts)¶
- quimb.tensor.tensor_builder.HTN_rand(n, reg, n_out=0, n_hyper_in=0, n_hyper_out=0, d_min=2, d_max=3, seed=None, dtype='float64', dist='normal', scale=1.0, loc=0.0, site_ind_id='k{}')¶
Create a random, possibly hyper, tensor network, with a mix of normal and hyper inner and outer indices. Useful for testing edges cases.
- Parameters:
n (int) – The number of tensors.
reg (int) – The average degree (number of dimensions per tensor) of the tensor network (prior to placing extra indices).
n_out (int, optional) – The number of normal outer indices to add (i.e. appear exactly once).
n_hyper_in (int, optional) – The number of hyper inner indices to add (i.e. appear on three or more tensors).
n_hyper_out (int, optional) – The number of hyper outer indices to add (i.e. appear in output and two or more tensors).
d_min (int, optional) – The minimum size of any dimension.
d_max (int, optional) – The maximum size of any dimension.
seed (int, optional) – A random seed.
- Return type:
- quimb.tensor.tensor_builder.HTN_CP_from_inds_and_fill_fn(fill_fn, inds, sizes, D, tags=None, bond_ind=None)¶
reate a CP-decomposition structured hyper tensor network from a sequence of indices and a fill function.
- Parameters:
fill_fn (callable) – A function that takes a shape and returns an array.
inds (sequence of str) – The outer indices of the network.
sizes (sequence of int) – The outer sizes of the network.
D (int) – The bond dimension of the inner hyper index.
tags (sequence of str, optional) – A tag for each tensor if supplied.
bond_ind (str, optional) – If given, a specific name for the inner hyper index.
- quimb.tensor.tensor_builder.HTN_CP_from_sites_and_fill_fn(fill_fn, sites, D, phys_dim=2, site_tag_id='I{}', site_ind_id='k{}', bond_ind=None)¶
Create a CP-decomposition structured hyper tensor network state from a sequence of sites and a fill function.
- Parameters:
fill_fn (callable) – A function that takes a shape and returns an array.
sites (sequence of hashable) – The sites of the tensor network.
D (int) – The hyper bond dimension connecting tensors.
phys_dim (int, optional) – The size of the outer, physical indices.
site_tag_id (str, optional) – String with formatter to tag sites.
site_ind_id (str, optional) – String with formatter to tag indices (if
phys_dim
specified).
- Return type:
- quimb.tensor.tensor_builder.HTN_CP_operator_from_products(array_seqs, upper_inds, lower_inds, tags_each=None, tags_all=None, bond_ind=None)¶
Construct a CP form hyper tensor network of the sum of matrix strings:
\[\sum_i A_i B_i C_i \ldots\]using a single hyper index to enumerate each sequence of operators.
- Parameters:
array_seqs (sequence[arrays_like]) – The arrays to use in each product, each of which should be 2D.
upper_inds (sequence[sequence[str]]) – The upper indices to use for each array in each product.
lower_inds (sequence[sequence[str]]) – The lower indices to use for each array in each product.
tags_each (sequence[sequence[str]], optional) – The tags to use for each array in each product.
tags_all (sequence[str], optional) – The tags to use for the whole hyper tensor network.
bond_ind (str, optional) – The name of the hyper index to use, by default randomly generated.
- Returns:
htn
- Return type:
- quimb.tensor.tensor_builder.HTN_dual_from_edges_and_fill_fn(fill_fn, edges, D, phys_dim=None, site_tag_id='I{}', site_ind_id='k{}')¶
Create a hyper tensor network with a tensor on each bond and a hyper index on each node.
- quimb.tensor.tensor_builder.convert_to_2d(tn, Lx=None, Ly=None, site_tag_id='I{},{}', x_tag_id='X{}', y_tag_id='Y{}', inplace=False)¶
Convert
tn
to aTensorNetwork2D
, assuming that is has a generic geometry with sites labelled by (i, j) coordinates already. Useful for constructing 2D tensor networks from functions that only require a list of edges etc.
- quimb.tensor.tensor_builder.TN2D_from_fill_fn(fill_fn, Lx, Ly, D, cyclic=False, site_tag_id='I{},{}', x_tag_id='X{}', y_tag_id='Y{}')¶
A scalar 2D lattice tensor network with tensors filled by a function.
- Parameters:
fill_fn (callable) – A function with signature
fill_fn(shape) -> array
, used to fill each tensor.Lx (int) – Length of side x.
Ly (int) – Length of side y.
D (int) – The bond dimension connecting sites.
cyclic (bool or (bool, bool), optional) – Whether to use periodic boundary conditions. X and Y can be specified separately using a tuple.
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 tags.
y_tag_id (str, optional) – String specifier for naming convention of column tags.
- Return type:
- quimb.tensor.tensor_builder.TN2D_empty(Lx, Ly, D, cyclic=False, site_tag_id='I{},{}', x_tag_id='X{}', y_tag_id='Y{}', dtype='float64')¶
A scalar 2D lattice tensor network initialized with empty tensors.
- Parameters:
Lx (int) – Length of side x.
Ly (int) – Length of side y.
D (int) – The bond dimension connecting sites.
cyclic (bool or (bool, bool), optional) – Whether to use periodic boundary conditions. X and Y can be specified separately using a tuple.
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 tags.
y_tag_id (str, optional) – String specifier for naming convention of column tags.
dtype (str, optional) – The data type of the tensors.
- Return type:
- quimb.tensor.tensor_builder.TN2D_with_value(value, Lx, Ly, D, cyclic=False, site_tag_id='I{},{}', x_tag_id='X{}', y_tag_id='Y{}', dtype=None)¶
A scalar 2D lattice tensor network with every element set to
value
. This usesnumpy.broadcast_to
and therefore essentially no memory.- Parameters:
value (scalar) – The value to fill the tensors with.
Lx (int) – Length of side x.
Ly (int) – Length of side y.
D (int) – The bond dimension connecting sites.
cyclic (bool or (bool, bool), optional) – Whether to use periodic boundary conditions. X and Y can be specified separately using a tuple.
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 tags.
y_tag_id (str, optional) – String specifier for naming convention of column tags.
dtype (str, optional) – The data type of the tensors.
- Return type:
- quimb.tensor.tensor_builder.TN2D_rand(Lx, Ly, D, cyclic=False, site_tag_id='I{},{}', x_tag_id='X{}', y_tag_id='Y{}', dist='normal', loc=0, scale=1, seed=None, dtype='float64')¶
A random scalar 2D lattice tensor network.
- Parameters:
Lx (int) – Length of side x.
Ly (int) – Length of side y.
D (int) – The bond dimension connecting sites.
cyclic (bool or (bool, bool), optional) – Whether to use periodic boundary conditions. X and Y can be specified separately using a tuple.
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 tags.
y_tag_id (str, optional) – String specifier for naming convention of column tags.
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.
scale (float, optional) – A multiplicative factor to scale the random numbers by.
seed (int, optional) – A random seed.
dtype (dtype, optional) – Data type of the random arrays.
- Return type:
- quimb.tensor.tensor_builder.TN2D_rand_symmetric(Lx, Ly, D, cyclic=False, site_tag_id='I{},{}', x_tag_id='X{}', y_tag_id='Y{}', dist='normal', loc=0, scale=1, seed=None, dtype='float64')¶
Create a random 2D lattice tensor network where every tensor is symmetric up to index permutations.
- Parameters:
Lx (int) – Length of side x.
Ly (int) – Length of side y.
D (int) – The bond dimension connecting sites.
cyclic (bool or (bool, bool), optional) – Whether to use periodic boundary conditions. X and Y can be specified separately using a tuple.
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 tags.
y_tag_id (str, optional) – String specifier for naming convention of column tags.
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.
scale (float, optional) – A multiplicative factor to scale the random numbers by.
seed (int, optional) – A random seed.
dtype (dtype, optional) – Data type of the random arrays.
- Return type:
- quimb.tensor.tensor_builder.TN2D_corner_double_line(Lx, Ly, line_dim=2, tiling=2, fill_missing_edges=True, site_tag_id='I{},{}', x_tag_id='X{}', y_tag_id='Y{}', **kwargs)¶
Build a 2D ‘corner double line’ (CDL) tensor network. Each plaquette contributes a matrix (by default the identity) at each corner, connected in a loop. The corners for each site are then grouped and optionally contracted. Such a tensor network has strong local correlations. See https://arxiv.org/abs/1412.0732. If the sites are not contracted, the resulting network is a product of loops that can be easily and exactly contracted.
Note that if identity matrices are used, the contracted value of the tensor network is
line_dim**num_plaquettes
.- Parameters:
Lx (int) – Length of side x.
Ly (int) – Length of side y.
line_dim (int, optional) – The dimension of the matrices at each corner. If contract is True, then the resulting bonds with have dimension line_dim**tiling.
tiling ({1, 2}, optional) – How to tile the plaquettes. If
1
, the plaquettes are tiled in a checkerboard pattern resulting in a single line per edge. If2
, the plaquettes are tiled in a dense pattern resulting in two lines per edge.fill_missing_edges (bool, optional) – Whether to fill in the missing edges around the border with open strings, ensuring every bond exists and has the same dimension.
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 tags.
y_tag_id (str, optional) – String specifier for naming convention of column tags.
kwargs – Additional keyword arguments are passed to
TN_from_strings()
.
- Return type:
See also
- quimb.tensor.tensor_builder.convert_to_3d(tn, Lx=None, Ly=None, Lz=None, site_tag_id='I{},{},{}', x_tag_id='X{}', y_tag_id='Y{}', z_tag_id='Z{}', inplace=False)¶
Convert
tn
to aTensorNetwork3D
, assuming that is has a generic geometry with sites labelled by (i, j, k) coordinates already. Useful for constructing 3D tensor networks from functions that only require a list of edges etc.
- quimb.tensor.tensor_builder.TN3D_from_fill_fn(fill_fn, Lx, Ly, Lz, D, cyclic=False, site_tag_id='I{},{},{}', x_tag_id='X{}', y_tag_id='Y{}', z_tag_id='Z{}')¶
A scalar 3D lattice tensor network with tensors filled by a function.
- Parameters:
fill_fn (callable) – A function with signature
fill_fn(shape) -> array
, used to fill each tensor.Lx (int) – Length of side x.
Ly (int) – Length of side y.
Lz (int) – Length of side z.
D (int) – The bond dimension connecting sites.
cyclic (bool or (bool, bool, bool), optional) – Whether to use periodic boundary conditions. X, Y and Z can be specified separately using a tuple.
site_tag_id (str, optional) – String formatter specifying how to label each site.
dtype (dtype, optional) – Data type of the random arrays.
- Return type:
- quimb.tensor.tensor_builder.TN3D_empty(Lx, Ly, Lz, D, cyclic=False, site_tag_id='I{},{},{}', x_tag_id='X{}', y_tag_id='Y{}', z_tag_id='Z{}', dtype='float64')¶
A scalar 3D lattice tensor network initialized with empty tensors.
- Parameters:
Lx (int) – Length of side x.
Ly (int) – Length of side y.
Lz (int) – Length of side z.
D (int) – The bond dimension connecting sites.
cyclic (bool or (bool, bool, bool), optional) – Whether to use periodic boundary conditions. X, Y and Z can be specified separately using a tuple.
site_tag_id (str, optional) – String formatter specifying how to label each site.
dtype (dtype, optional) – Data type of the random arrays.
seed (int, optional) – Random seed.
- Return type:
- quimb.tensor.tensor_builder.TN3D_with_value(value, Lx, Ly, Lz, D, cyclic=False, site_tag_id='I{},{},{}', x_tag_id='X{}', y_tag_id='Y{}', z_tag_id='Z{}', dtype=None)¶
A scalar 2D lattice tensor network with every element set to
value
. This usesnumpy.broadcast_to
and therefore essentially no memory.- Parameters:
value (scalar) – The value to fill the tensors with.
Lx (int) – Length of side x.
Ly (int) – Length of side y.
Lz (int) – Length of side z.
D (int) – The bond dimension connecting sites.
cyclic (bool or (bool, bool, bool), optional) – Whether to use periodic boundary conditions. X, Y and Z can be specified separately using a tuple.
site_tag_id (str, optional) – String formatter specifying how to label each site.
dtype (dtype, optional) – Data type of the random arrays.
seed (int, optional) – Random seed.
- Return type:
- quimb.tensor.tensor_builder.TN3D_rand(Lx, Ly, Lz, D, cyclic=False, site_tag_id='I{},{},{}', x_tag_id='X{}', y_tag_id='Y{}', z_tag_id='Z{}', dist='normal', loc=0.0, scale=1.0, seed=None, dtype='float64')¶
A random scalar 3D lattice tensor network.
- Parameters:
Lx (int) – Length of side x.
Ly (int) – Length of side y.
Lz (int) – Length of side z.
D (int) – The bond dimension connecting sites.
cyclic (bool or (bool, bool, bool), optional) – Whether to use periodic boundary conditions. X, Y and Z can be specified separately using a tuple.
site_tag_id (str, optional) – String formatter specifying how to label each site.
dtype (dtype, optional) – Data type of the random arrays.
seed (int, optional) – Random seed.
- Return type:
- quimb.tensor.tensor_builder.TN3D_corner_double_line(Lx, Ly, Lz, line_dim=2, tiling=2, fill_missing_edges=True, site_tag_id='I{},{},{}', x_tag_id='X{}', y_tag_id='Y{}', z_tag_id='Z{}', **kwargs)¶
- quimb.tensor.tensor_builder.classical_ising_S_matrix(beta, j=1.0)¶
The interaction term for the classical ising model.
- quimb.tensor.tensor_builder.classical_ising_H_matrix(beta, h=0.0)¶
The magnetic field term for the classical ising model.
- quimb.tensor.tensor_builder.classical_ising_sqrtS_matrix(beta, j=1.0, asymm=None)¶
The sqrt factorized interaction term for the classical ising model. If
j
is negative you can supplyasymm='l'
or'r'
to keep the matrix real, but it must be paired with the opposite in a tensor network.
- quimb.tensor.tensor_builder.classical_ising_T_matrix(beta, j=1.0, h=0.0, directions='lrud', output=False, asymm=None)¶
The single effective TN site for the classical ising model.
- quimb.tensor.tensor_builder.parse_j_coupling_to_function(j)¶
Parse the
j
argument to a function that can be called to get the coupling strength between two nodes. The input can be a constant, dict or function.
- quimb.tensor.tensor_builder.HTN2D_classical_ising_partition_function(Lx, Ly, beta, h=0.0, j=1.0, cyclic=False, ind_id='s{},{}', site_tag_id='I{},{}')¶
Hyper tensor network representation of the 2D classical ising model partition function. The indices will be shared by 4 or 5 tensors depending on whether
h
is non-zero. As opposed to the ‘normal’ tensor network, here each classical spin is still a single index, which is easier to contract exactly.- Parameters:
Lx (int) – Length of side x.
Ly (int) – Length of side y.
beta (float) – The inverse temperature.
h (float, optional) – The magnetic field strength.
j (float, dict, or callable, optional) – The interaction strength, positive being ferromagnetic. If a dict should contain entries for each edge, keyed by the edge. If a callable should have the signature
j(node_a, node_b)
and return a float.cyclic (bool or (bool, bool), optional) – Whether to use periodic boundary conditions. X and Y can be specified separately using a tuple.
ind_id (str, optional) – How to label the indices i.e.
ind_id.format(i, j)
, each of which corresponds to a single classical spin.site_tag_id (str, optional) – How to label the site tags, note that in the hyper tensor network representation each tensor will have two site tags for the two sites it connects.
- Return type:
- quimb.tensor.tensor_builder.HTN3D_classical_ising_partition_function(Lx, Ly, Lz, beta, j=1.0, h=0.0, cyclic=False, ind_id='s{},{},{}', site_tag_id='I{},{},{}')¶
Hyper tensor network representation of the 3D classical ising model partition function. The indices will be shared by 6 or 7 tensors depending on whether
h
is non-zero. As opposed to the ‘normal’ tensor network, here each classical spin is still a single index, which is easier to contract exactly.- Parameters:
Lx (int) – Length of side x.
Ly (int) – Length of side y.
Lz (int) – Length of side z.
beta (float) – The inverse temperature.
j (float, dict, or callable, optional) – The interaction strength, positive being ferromagnetic. If a dict should contain entries for each edge, keyed by the edge. If a callable should have the signature
j(node_a, node_b)
and return a float.h (float, optional) – The magnetic field strength.
cyclic (bool or (bool, bool, bool), optional) – Whether to use periodic boundary conditions. X, Y and Z can be specified separately using a tuple.
ind_id (str, optional) – How to label the indices i.e.
ind_id.format(i, j, k)
, each of which corresponds to a single classical spin.site_tag_id (str, optional) – How to label the site tags, note that in the hyper tensor network representation each tensor will have two site tags for the two sites it connects.
- Return type:
- quimb.tensor.tensor_builder.TN2D_classical_ising_partition_function(Lx, Ly, beta, j=1.0, h=0.0, cyclic=False, site_tag_id='I{},{}', x_tag_id='X{}', y_tag_id='Y{}', outputs=(), ind_id='s{},{}')¶
The tensor network representation of the 2D classical ising model partition function.
- Parameters:
Lx (int) – Length of side x.
Ly (int) – Length of side y.
beta (float) – The inverse temperature.
j (float, dict, or callable, optional) – The interaction strength, positive being ferromagnetic. If a dict should contain entries for each edge, keyed by the edge. If a callable should have the signature
j(node_a, node_b)
and return a float.h (float, optional) – The magnetic field strength.
cyclic (bool or (bool, bool), optional) – Whether to use periodic boundary conditions. X and Y can be specified separately using a tuple.
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 tags.
y_tag_id (str, optional) – String specifier for naming convention of column tags.
outputs (sequence of tuple[int, int], optional) – Which sites to generate output indices (i.e. dangling legs) for. The index is named according to
ind_id
.ind_id (str, optional) – How to label the indices i.e.
ind_id.format(i, j)
, each of which corresponds to a single classical spin inoutputs
.
- Return type:
- quimb.tensor.tensor_builder.TN3D_classical_ising_partition_function(Lx, Ly, Lz, beta, j=1.0, h=0.0, cyclic=False, site_tag_id='I{},{},{}', x_tag_id='X{}', y_tag_id='Y{}', z_tag_id='Z{}', outputs=(), ind_id='s{},{},{}')¶
Tensor network representation of the 3D classical ising model partition function.
- Parameters:
Lx (int) – Length of side x.
Ly (int) – Length of side y.
Lz (int) – Length of side z.
beta (float) – The inverse temperature.
j (float, dict, or callable, optional) – The interaction strength, positive being ferromagnetic. If a dict should contain entries for each edge, keyed by the edge. If a callable should have the signature
j(node_a, node_b)
and return a float.h (float, optional) – The magnetic field strength.
cyclic (bool or (bool, bool, bool), optional) – Whether to use periodic boundary conditions. X, Y and Z can be specified separately using a tuple.
site_tag_id (str, optional) – String formatter specifying how to label each site.
x_tag_id (str, optional) – String formatter specifying how to label each x-plane.
y_tag_id (str, optional) – String formatter specifying how to label each y-plane.
z_tag_id (str, optional) – String formatter specifying how to label each z-plane.
outputs (sequence of tuple[int, int, int], optional) – Which sites to generate output indices (i.e. dangling legs) for. The index is named according to
ind_id
.ind_id (str, optional) – How to label the indices i.e.
ind_id.format(i, j, k)
, each of which corresponds to a single classical spin inoutputs
.
- Return type:
- quimb.tensor.tensor_builder.HTN_classical_partition_function_from_edges(edges, beta, j=1.0, h=0.0, site_ind_id='s{}', site_tag_id='I{}', bond_tag_id='B{},{}')¶
Build a hyper tensor network representation of a classical ising model partition function by specifying graph edges. There will be a single tensor per interaction rather than per site, as well as a single tensor for each site, if
h != 0.0
.- Parameters:
edges (sequence of tuple[hashable, hashable]) – The graph edges, as a sequence of pairs of hashable objects, for example integers, representing the nodes. You can redundantly specify
(u, v)
and(v, u)
and only one edge will be added.beta (float, optional) – The inverse temperature.
j (float, dict, or callable, optional) – The interaction strength, positive being ferromagnetic. If a dict should contain entries for each edge, keyed by the edge. If a callable should have the signature
j(node_a, node_b)
and return a float.h (float, or callable, optional) – The magnetic field strength. If a callable should have the signature
h(node)
and return a float.site_ind_id (str, optional) – A string formatter for naming tensor indices like
site_ind_id.format(node)
.site_tag_id (str, optional) – A string formatter for naming tensor tags like
site_tag_id.format(node)
.bond_tag_id (str, optional) – A string formatter for naming tensor tags like
bond_tag_id.format(node_a, node_b)
.
- Return type:
- quimb.tensor.tensor_builder.TN_classical_partition_function_from_edges(edges, beta, j=1.0, h=0.0, site_tag_id='I{}', bond_ind_id='b{},{}', outputs=(), ind_id='s{}')¶
Build a regular tensor network representation of a classical ising model partition function by specifying graph edges. There will be a single tensor per site.
- Parameters:
edges (sequence of tuple[hashable, hashable]) – The graph edges, as a sequence of pairs of hashable objects, for example integers, representing the nodes. You can redundantly specify
(u, v)
and(v, u)
and only one edge will be added.beta (float, optional) – The inverse temperature.
j (float, dict, or callable, optional) – The interaction strength, positive being ferromagnetic. If a dict should contain entries for each edge, keyed by the edge. If a callable should have the signature
j(node_a, node_b)
and return a float.h (float, or callable, optional) – The magnetic field strength. If a callable should have the signature
h(node)
and return a float.site_tag_id (str, optional) – A string formatter for naming tensor tags like
site_ind_id.format(node)
.bond_ind_id (str, optional) – A string formatter for naming the indices bewteen tensors like
bond_ind_id.format(node_a, node_b)
.
- Return type:
- quimb.tensor.tensor_builder.make_couplings_matrix_symmetric(J, UPLO='auto')¶
Take a coupling matrix or dict of pairwise strengths in either upper or lower triangular form and return a symmetric matrix.
- Parameters:
J (array_like or dict[(int, int), float]) – The couplings, as a square matrix or explicit dict. If a matrix, the upper or lower triangle (corresponding to
UPLO
) should the pairwise interaction strengths. If a dict, then each non-zero interaction can be specified explicitly.UPLO ({"auto", "L", "U"}, optional) – Whether to read the couplings in from
J
in lower (‘L’:i >= j
) or upper (‘U’:i <= j
) form. IfJ
is a dict then this only matters if both are present and in that case acts as a preference. If"auto"
then either the lower or upper, or symmetric distances can be supplied, but an error is raised if both are.
- Returns:
J – The couplings as a symmetric matrix.
- Return type:
array_like
- quimb.tensor.tensor_builder.TN2D_embedded_classical_ising_partition_function(Jij, beta, outputs=(), ordering=None, sites_location='side', UPLO='auto', contract_sites=True, site_tag_id='I{},{}', x_tag_id='X{}', y_tag_id='Y{}', ind_id='s{}')¶
Construct a (triangular) ‘2D’ tensor network representation of the classical ising model partition function with all-to-all interactions, by embedding it in a 2D lattice. For N spins this will have N(N-1)/2 tensors, present at sites (i, j) where i > j. If outputs is supplied, then dangling indices will be added for each output, to give a ‘standard’ TN representation of the unnormalized marginal distribution over those outputs. Since each spin is effectively “delocalized” into a COPY-tensor MPS spanning the lattice, there is a choice of where to put the outputs, specified by sites_location.
sites_location="side"
:... . . . ╭── │ j=2 . . ╭──O──... │ │ j=1 . ╭──O──O──... │ │ │ j=0 ╭──O──O──O──... │ │ │ │ s0 s1 s2 s3 <- site indices i=0, 1, 2, 3, ...
sites_location="diag"
:... . . s3 ╭── ╲│ j=2 . s2 ╭──O──... ╲│ │ j=1 s1 ╭──O──O──... ╲│ │ │ j=0 s0──O──O──O──... i=0, 1, 2, 3, ...
- Parameters:
Jij (array_like or dict[(int, int), float]) – The couplings, as a square matrix or explicit dict. If a matrix, the upper or lower triangle (corresponding to
UPLO
) should the pairwise interaction strengths. If a dict, then each non-zero interaction can be specified explicitly.beta (float) – The inverse temperature.
outputs (sequence of int, optional) – Which sites to generate output indices (i.e. dangling legs) for. The index is named according to
ind_id
. The location of the output depends onsites_location
.ordering (sequence of int, optional) – If supplied, the ordering of the embedded spins into the lattice.
sites_location ({"side", "diag"}, optional) – Whether to place the output indices on the side or diagonally. If
"side"
then the output indices will be placed along the ‘y-min’ border, i.e.[(j, 0) for j in range(1, N)]
, with two on the first tensor. If"diag"
then the output indices will be placed along the central diagonal fold of the lattice, i.e.[(j + 1, j) for j in range(1, N)]
, with two on the first tensor.UPLO ({"auto", "L", "U"}, optional) – Whether to read the couplings in from
J
in lower (‘L’:i >= j
) or upper (‘U’:i <= j
) form. IfJ
is a dict then this only matters if both are present and in that case acts as a preference. If"auto"
then either the lower or upper, or symmetric distances can be supplied, but an error is raised if both are.contract_sites (bool, optional) – Whether to contract the sites of the TN, which after initial construction will have one horizontol spin-rail, one vertical spin-rail and one interaction tensor per site.
site_tag_id (str, optional) – String formatter specifying how to label each site.
x_tag_id (str, optional) – String formatter specifying how to label each x-plane.
y_tag_id (str, optional) – String formatter specifying how to label each y-plane.
ind_id (str, optional) – How to label the indices i.e.
ind_id.format(i)
, each of which corresponds to a single classical spin inoutputs
.
- Return type:
- quimb.tensor.tensor_builder.dimer_data(d, cover_count=1, dtype=float)¶
- quimb.tensor.tensor_builder.TN_dimer_covering_from_edges(edges, cover_count=1, site_tag_id='I{}', bond_ind_id='b{},{}', dtype=float)¶
Make a tensor network from sequence of graph edges that counts the number of ways to cover the graph exactly with dimers. See https://arxiv.org/abs/1805.10598 for the construction.
- Parameters:
edges (sequence of tuple) – The edges, each item should be a pair of hashable objects describing nodes linked.
cover_count (int, optional) – The exact number of times each node must be ‘covered’. For example 1 for a standard dimer covering or 2 for ‘ice rules’.
site_tag_id (str, optional) – A string formatter for naming tensor tags like
site_ind_id.format(node)
.bond_ind_id (str, optional) – A string formatter for naming the indices bewteen tensors like
bond_ind_id.format(node_a, node_b)
.
- Return type:
- quimb.tensor.tensor_builder.clause_negmask(clause)¶
Encode a clause as a single integer
m
.
- quimb.tensor.tensor_builder.or_clause_data(ndim, m=0, dtype=float, q=2)¶
Get the array representing satisfiability of
ndim
clauses with unsatisfied condition encoded inm
.
- quimb.tensor.tensor_builder.or_clause_tensor(ndim, m, inds, tags=None, dtype='float64')¶
Get the tensor representing satisfiability of
ndim
clauses with unsatisfied condition encoded inm
labelled byinds
andtags
.
- quimb.tensor.tensor_builder.or_clause_mps_tensors(ndim, m, inds, tags=None, dtype='float64')¶
Get the set of MPS tensors representing satisfiability of
ndim
clauses with unsatisfied condition encoded inm
labelled byinds
andtags
.
- quimb.tensor.tensor_builder.or_clause_parafac_data(ndim, m, dtype='float64')¶
Get the set of PARAFAC arrays representing satisfiability of
ndim
clauses with unsatisfied condition encoded inm
.
- quimb.tensor.tensor_builder.clause_parafac_tensors(ndim, m, inds, tags=None, dtype='float64')¶
Get the set of PARAFAC tensors representing satisfiability of
ndim
clauses with unsatisfied condition encoded inm
labelled byinds
andtags
.
- quimb.tensor.tensor_builder.HTN_from_clauses(clauses, weights=None, mode='parafac', dtype='float64', clause_tag_id='CLAUSE{}', var_ind_id='var{}', weight_tag_id='WEIGHT{}')¶
Given a list of clauses, create a hyper tensor network, with a single hyper index for each variable, and single tensor or tensor decomposition for each clause. If weights are given, there will also be a single tensor for each non-trivially weighted variable.
- Parameters:
clauses (sequence of tuple[int]) – The clauses as a sequence of tuples of integers. Each integer represents a variable, and the sign indicates whether it is negated. The variables thus must be non-zero integers.
weights (dict[int, float], optional) – The weights for each variable. Each key should be a signed variable integer, such that relative weights for a variable
v
are(weights[-v], weights[v])
. If only one is given of this pair, the other is assumed to sum to 1. If a variable is not supplied, orweights=None
, then both weights are assumed to be 1 and no tensor is created for the variable.mode ({'parafac', 'mps', 'dense', int}, optional) –
How to represent the clauses:
’parafac’ - N rank-2 tensors connected by a single hyper index. You could further call
hyperinds_resolve()
for more options to convert the hyper index into a (decomposed) COPY-tensor.’mps’ - N rank-3 tensors connected along a 1D line.
’dense’ - contract the hyper index.
int - use the ‘parafac’ mode, but only if the length of a clause is larger than this threshold.
Note that variables are always represented by a single (hyper) index, which is like an implicit PARAFAC decomposition.
dtype (str) – The data type of the tensors.
clause_tag_id (str) – The tag to use for the clause tensors. The tag will be formatted with the clause index.
var_ind_id (str) – The index to use for the variable tensors. The index will be formatted with the variable index.
weight_tag_id (str) – The tag to use for the weight tensors. The tag will be formatted with the variable index.
- Returns:
htn
- Return type:
- quimb.tensor.tensor_builder.cnf_file_parse(fname)¶
Parse a DIMACS style ‘cnf’ file into a list of clauses, and possibly a dictionary of weights. The weights, if present, can be specified either as:
(CACHET format): a line per weight like
w {signed_var} {weight}
, wheresigned_var
is an integer whose sign specifies the sign of the weight being set.(MC2021 competition format): the same as above, but with each line specified as
c p weight {signed_var} {weight}
.(MINIC2D format): a single line of the form
c weights {wp_1} {wm_1} {wp_2} {wm_2}... ``, where ``wp_i
andwn_i
are the positive and negative weights for variablei
. Weights specified this way are overriden by the previous two formats.
- quimb.tensor.tensor_builder.HTN_from_cnf(fname, mode='parafac', dtype='float64', clause_tag_id='CLAUSE{}', var_ind_id='var{}', weight_tag_id='WEIGHT{}', **kwargs)¶
Create a hyper tensor network from a ‘.cnf’ or ‘.wcnf’ file - i.e. a model counting or weighted model counting instance specification.
- Parameters:
fname (str) – Path to a ‘.cnf’ or ‘.wcnf’ file.
mode ({'parafac', 'mps', 'dense', int}, optional) –
How to represent the clauses:
’parafac’ - N rank-2 tensors connected by a single hyper index. You could further call
hyperinds_resolve()
for more options to convert the hyper index into a (decomposed) COPY-tensor.’mps’ - N rank-3 tensors connected along a 1D line.
’dense’ - contract the hyper index.
int - use the ‘parafac’ mode, but only if the length of a clause is larger than this threshold.
dtype (str or dtype, optional) – Data type of the tensors.
clause_tag_id (str, optional) – Format string for clause tags.
var_ind_id (str, optional) – Format string for variable indices.
weight_tag_id (str, optional) – Format string for weight tags.
kwargs – Additional keyword arguments passed to
HTN_from_clauses()
.
- Returns:
htn
- Return type:
- quimb.tensor.tensor_builder.random_ksat_instance(k, num_variables, num_clauses=None, alpha=None, seed=None, allow_repeat_variables=False)¶
Create a random k-SAT instance.
- Parameters:
k (int) – Number of variables per clause.
num_variables (int) – Number of variables in the instance.
num_clauses (int, optional) – Number of clauses in the instance. If not specified, will be determined from alpha.
alpha (float, optional) – If num_clauses is not directly specified then the average number of clauses per variable. Taken as a Poisson parameter. Either this or num_clauses must be specified.
seed (int, optional) – Random seed.
allow_repeat_variables (bool, optional) – Whether to allow the same variable to appear multiple times in a single clause.
- Returns:
instance – Dictionary with keys ‘num_variables’, ‘num_clauses’, ‘clauses’. The ‘clauses’ key contains a list of tuples, each tuple representing a clause. Within each tuple, each element is an integer representing a variable, with the sign of the integer representing the sign of the variable in the clause.
- Return type:
- quimb.tensor.tensor_builder.HTN_random_ksat(k, num_variables, num_clauses=None, alpha=None, seed=None, allow_repeat_variables=False, mode='parafac', dtype='float64', clause_tag_id='CLAUSE{}', variable_ind_id='var{}')¶
Create a random k-SAT instance encoded as a hyper tensor network.
- Parameters:
k (int) – Number of variables per clause.
num_variables (int) – Number of variables in the instance.
num_clauses (int, optional) – Number of clauses in the instance. If not specified, will be determined from alpha.
alpha (float, optional) – If num_clauses is not directly specified then the average number of clauses per variable. Taken as a Poisson parameter. Either this or num_clauses must be specified.
seed (int, optional) – Random seed.
allow_repeat_variables (bool, optional) – Whether to allow the same variable to appear multiple times in a single clause.
mode ({'parafac', 'mps', 'dense', int}, optional) –
How to represent the clauses:
’parafac’ - N rank-2 tensors connected by a single hyper index. You could further call
hyperinds_resolve()
for more options to convert the hyper index into a (decomposed) COPY-tensor.’mps’ - N rank-3 tensors connected along a 1D line.
’dense’ - contract the hyper index.
int - use the ‘parafac’ mode, but only if the length of a clause is larger than this threshold.
Note that variables are always represented by a single (hyper) index, which is like an implicit PARAFAC decomposition.
dtype (str, optional) – Data type of the tensors.
clause_tag_id (str, optional) – Format string for clause tags. Should contain a single {} which will be replaced by the clause number.
variable_ind_id (str, optional) – Format string for variable indices. Should contain a single {} which will be replaced by the variable number.
- Return type:
- quimb.tensor.tensor_builder.TN_matching(tn, max_bond, site_tags=None, fill_fn=None, dtype=None, **randn_opts)¶
Create a tensor network with the same outer indices as
tn
but with a single tensor per site with bond dimensionmax_bond
between each connected site. Generally to be used as an initial guess for fitting.- Parameters:
tn (TensorNetwork) – The tensor network to match, it can have arbitrary local structure and output indices, as long as
site_tags
effectively partitions it.max_bond (int) – The bond dimension to use between each site.
site_tags (sequence of str, optional) – The tags to use to select the tensors from
tn
. If not given, usestn.site_tags
. The tensor network built will have one tensor per site, in the order given bysite_tags
.dtype (dtype, optional) – The data type to use for the new tensors, if not given uses the same as the original tensors.
randn_opts – Supplied to
randn()
.
- Return type:
- quimb.tensor.tensor_builder.MPS_rand_state(L, bond_dim, phys_dim=2, normalize=True, cyclic=False, dist='normal', loc=0.0, scale=1.0, dtype='float64', trans_invar=False, **mps_opts)¶
Generate a random matrix product state.
- Parameters:
L (int) – The number of sites.
bond_dim (int) – The bond dimension.
phys_dim (int, optional) – The physical (site) dimensions, defaults to 2.
normalize (bool, optional) – Whether to normalize the state.
cyclic (bool, optional) – Generate a MPS with periodic boundary conditions or not, default is open boundary conditions.
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.
scale (float, optional) – A multiplicative factor to scale the random numbers by.
dtype ({float, complex} or numpy dtype, optional) – Data type of the tensor network.
trans_invar (bool (optional)) – Whether to generate a translationally invariant state, requires cyclic=True.
mps_opts – Supplied to
MatrixProductState
.
- quimb.tensor.tensor_builder.MPS_product_state(arrays, cyclic=False, **mps_opts)¶
Generate a product state in MatrixProductState form, i,e, with bond dimension 1, from single site vectors described by
arrays
.
- quimb.tensor.tensor_builder.MPS_computational_state(binary, dtype='float64', cyclic=False, **mps_opts)¶
A computational basis state in Matrix Product State form.
- Parameters:
binary (str or sequence of int) – String specifying the state, e.g.
'00101010111'
or[0, 0, 1]
.dtype ({'float64', 'complex128', 'float32', 'complex64'}, optional) – The data type to use for the array representation.
cyclic (bool, optional) – Generate a MPS with periodic boundary conditions or not, default open boundary conditions.
mps_opts – Supplied to MatrixProductState constructor.
- quimb.tensor.tensor_builder.MPS_neel_state(L, down_first=False, dtype='float64', **mps_opts)¶
Generate the neel state in Matrix Product State form.
- quimb.tensor.tensor_builder.MPS_COPY(L, phys_dim=2, dtype='float64', **mps_opts)¶
Build a matrix product state representation of the COPY tensor.
- Parameters:
- Return type:
- quimb.tensor.tensor_builder.MPS_ghz_state(L, dtype='float64', **mps_opts)¶
Build the chi=2 OBC MPS representation of the GHZ state.
- Parameters:
L (int) – Number of qubits.
dtype ({'float64', 'complex128', 'float32', 'complex64'}, optional) – The underlying data type.
mps_opts – Supplied to
MatrixProductState
.
- quimb.tensor.tensor_builder.MPS_w_state(L, dtype='float64', **mps_opts)¶
Build the chi=2 OBC MPS representation of the W state.
- Parameters:
L (int) – Number of qubits.
dtype ({'float64', 'complex128', 'float32', 'complex64'}, optional) – The underlying data type.
mps_opts – Supplied to
MatrixProductState
.
- quimb.tensor.tensor_builder.MPS_rand_computational_state(L, dtype='float64', **mps_opts)¶
Generate a random computation basis state, like ‘01101001010’.
- Parameters:
L (int) – The number of qubits.
seed (int, optional) – The seed to use.
dtype ({float, complex} or numpy dtype, optional) – Data type of the tensor network.
mps_opts – Supplied to
MatrixProductState
.
- quimb.tensor.tensor_builder.MPS_zero_state(L, bond_dim=1, phys_dim=2, cyclic=False, dtype='float64', **mps_opts)¶
The all-zeros MPS state, of given bond-dimension.
- Parameters:
L (int) – The number of sites.
bond_dim (int, optional) – The bond dimension, defaults to 1.
phys_dim (int, optional) – The physical (site) dimensions, defaults to 2.
cyclic (bool, optional) – Generate a MPS with periodic boundary conditions or not, default is open boundary conditions.
dtype ({float, complex} or numpy dtype, optional) – Data type of the tensor network.
mps_opts – Supplied to
MatrixProductState
.
- quimb.tensor.tensor_builder.MPS_sampler(L, dtype=complex, squeeze=True, **mps_opts)¶
A product state for sampling tensor network traces. Seen as a vector it has the required property that
psi.H @ psi == d
always for hilbert space sized
.
- quimb.tensor.tensor_builder.MPO_identity(L, sites=None, phys_dim=2, dtype='float64', cyclic=False, **mpo_opts)¶
Generate an identity MPO of size
L
.- Parameters:
L (int) – The number of sites.
phys_dim (int, optional) – The physical (site) dimensions, defaults to 2.
dtype ({float, complex} or numpy dtype, optional) – Data type of the tensor network.
cyclic (bool, optional) – Generate a MPO with periodic boundary conditions or not, default is open boundary conditions.
mpo_opts – Supplied to
MatrixProductOperator
.
- quimb.tensor.tensor_builder.MPO_identity_like(mpo, **mpo_opts)¶
Return an identity matrix operator with the same physical index and inds/tags as
mpo
.
- quimb.tensor.tensor_builder.MPO_zeros(L, phys_dim=2, dtype='float64', cyclic=False, **mpo_opts)¶
Generate a zeros MPO of size
L
.- Parameters:
L (int) – The number of sites.
phys_dim (int, optional) – The physical (site) dimensions, defaults to 2.
dtype ({float, complex} or numpy dtype, optional) – Data type of the tensor network.
cyclic (bool, optional) – Generate a MPO with periodic boundary conditions or not, default is open boundary conditions.
mpo_opts – Supplied to
MatrixProductOperator
.
- Return type:
- quimb.tensor.tensor_builder.MPO_zeros_like(mpo, **mpo_opts)¶
Return a zeros matrix product operator with the same physical index and inds/tags as
mpo
.- Parameters:
mpo (MatrixProductOperator) – The MPO to copy the shape of.
- Return type:
- quimb.tensor.tensor_builder.MPO_product_operator(arrays, cyclic=False, **mpo_opts)¶
Return an MPO of bond dimension 1 representing the product of raw operators given in
arrays
.- Parameters:
arrays (sequence of 2D array_like) – The operators to form a tensor product of.
cyclic (bool, optional) – Whether to generate a cyclic MPO or not.
mpo_opts – Supplied to
MatrixProductOperator
.
- Return type:
- quimb.tensor.tensor_builder.MPO_rand(L, bond_dim, phys_dim=2, normalize=True, cyclic=False, herm=False, dtype='float64', dist='normal', loc=0.0, scale=1.0, **mpo_opts)¶
Generate a random matrix product state.
- Parameters:
L (int) – The number of sites.
bond_dim (int) – The bond dimension.
phys_dim (int, optional) – The physical (site) dimensions, defaults to 2.
normalize (bool, optional) – Whether to normalize the operator such that
trace(A.H @ A) == 1
.cyclic (bool, optional) – Generate a MPO with periodic boundary conditions or not, default is open boundary conditions.
dtype ({float, complex} or numpy dtype, optional) – Data type of the tensor network.
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.
scale (float, optional) – A multiplicative factor to scale the random numbers by.
herm (bool, optional) – Whether to make the matrix hermitian (or symmetric if real) or not.
mpo_opts – Supplied to
MatrixProductOperator
.
- quimb.tensor.tensor_builder.MPO_rand_herm(L, bond_dim, phys_dim=2, normalize=True, dtype='float64', **mpo_opts)¶
Generate a random hermitian matrix product operator. See
MPO_rand
.
- quimb.tensor.tensor_builder.maybe_make_real(X)¶
Check if
X
is real, if so, convert to contiguous array.
- quimb.tensor.tensor_builder.spin_ham_mpo_tensor(one_site_terms, two_site_terms, S=1 / 2, left_two_site_terms=None, which=None, cyclic=False)¶
Generate tensor(s) for a spin hamiltonian MPO.
- Parameters:
one_site_terms (sequence of (scalar, operator)) – The terms that act on a single site, each
operator
can be a string suitable to be sent tospin_operator()
or an actual 2d-array.two_site_terms (sequence of (scalar, operator operator)) – The terms that act on two neighbouring sites, each
operator
can be a string suitable to be sent tospin_operator()
or an actual 2d-array.S (fraction, optional) – What size spin to use, defaults to spin-1/2.
left_two_site_terms (sequence of (scalar, operator operator), optional) – If the interaction to the left of this site has different spin terms then the equivalent list of terms for that site.
which ({None, 'L', 'R', 'A'}, optional) – If
None
, generate the middle tensor, if ‘L’ a left-end tensor, if ‘R’ a right-end tensor and if ‘A’ all three.cyclic (bool, optional) – Whether to use periodic boundary conditions - default False.
- Returns:
The middle, left, right or all three MPO tensors.
- Return type:
- class quimb.tensor.tensor_builder._TermAdder(terms, nsite)¶
Simple class to allow
SpinHam1D
syntax likebuilder[i, j] += (1/2, 'Z', 'X')
. This object is temporarily created by the getitem call, accumulates the new term, then has its the new combined list of terms extracted in the setitem call.- terms¶
- nsite¶
- __iadd__(new)¶
- class quimb.tensor.tensor_builder.SpinHam1D(S=1 / 2, cyclic=False)¶
Class for easily building custom spin hamiltonians in MPO or LocalHam1D form. Currently limited to nearest neighbour interactions (and single site terms). It is possible to set ‘default’ translationally invariant terms, but also terms acting on specific sites only (which take precedence). It is also possible to build a sparse matrix version of the hamiltonian (obviously for small sizes only).
- Parameters:
Examples
Initialize the spin hamiltonian builder:
>>> builder = SpinHam1D(S=3 / 2)
Add some two-site terms:
>>> builder += 0.5, '+', '-' >>> builder += 0.5, '-', '+' >>> builder += 1.0, 'Z', 'Z'
Add a single site term:
>>> builder -= 0.3, 'Z'
Build a MPO version of the hamiltonian for use with DMRG:
>>> mpo_ham = builder.build_mpo(100) >>> mpo_ham <MatrixProductOperator(tensors=100, L=100, max_bond=5)>
Build a LocalHam1D version of the hamiltonian for use with TEBD:
>>> builder.build_local_ham(100) <LocalHam1D(L=100, cyclic=False)>
You can also set terms for specific sites (this overides any of the ‘default’, translationally invariant terms set as above):
>>> builder[10, 11] += 0.75, '+', '-' >>> builder[10, 11] += 0.75, '-', '+' >>> builder[10, 11] += 1.5, 'Z', 'Z'
Or specific one-site terms (which again overides any default single site terms set above):
>>> builder[10] += 3.7, 'Z' >>> builder[11] += 0.0, 'I' # '0' term turns off field
- S = 0.5¶
- one_site_terms = []¶
- two_site_terms = []¶
- cyclic = False¶
- var_one_site_terms¶
- var_two_site_terms¶
- add_term(factor, *operators)¶
Add another term to the expression to be built.
- Parameters:
factor (scalar) – Scalar factor to multiply this term by.
*operators (str or array) – The operators to use. Can specify one or two for single or two site terms respectively. Can use strings, which are supplied to
spin_operator()
, or actual arrays as long as they have the correct dimension.
- sub_term(factor, *operators)¶
Subtract a term - simple alias that flips sign of
factor
.
- __iadd__(term)¶
- __isub__(term)¶
- __getitem__(sites)¶
Part of the machinery that allows terms to be added to specific sites like:
>>> builder[i] += 1/2, 'X' >>> builder[45, 46] += 1/2, 'Z', 'Z'
- __setitem__(sites, value)¶
Part of the machinery that allows terms to be added to specific sites like:
>>> builder[i] += 1/2, 'X' >>> builder[45, 46] += 1/2, 'Z', 'Z'
Could also be called directly with a list of terms like:
>>> builder[13, 14] = [(1, 'Z', 'Z'), (0.5, 'X', 'Y')]
Which would overide any terms set so far.
- build_mpo(L, upper_ind_id='k{}', lower_ind_id='b{}', site_tag_id='I{}', tags=None)¶
Build an MPO instance of this spin hamiltonian of size
L
. See alsoMatrixProductOperator
.
- build_sparse(L, **ikron_opts)¶
Build a sparse matrix representation of this Hamiltonian.
- _get_spin_op(factor, *ss)¶
- _sum_spin_ops(terms)¶
- build_local_ham(L=None, **local_ham_1d_opts)¶
Build a nearest neighbour interactor instance of this spin hamiltonian of size
L
. See alsoLocalHam1D
.- Parameters:
L (int, optional) – The number of spins, if the hamiltonian only has two-site terms this is optional.
- Return type:
- quimb.tensor.tensor_builder.SpinHam¶
- quimb.tensor.tensor_builder._ham_ising(j=1.0, bx=0.0, *, S=1 / 2, cyclic=False)¶
- quimb.tensor.tensor_builder.MPO_ham_ising(L, j=1.0, bx=0.0, *, S=1 / 2, cyclic=False, **mpo_opts)¶
Ising Hamiltonian in MPO form.
\[H_\mathrm{Ising} = J \sum_{i} S^Z_i S^Z_{i + 1} - B_x \sum_{i} S^X_i\]Note the default convention of antiferromagnetic interactions and spin operators not Pauli matrices.
- Parameters:
L (int) – The number of sites.
j (float, optional) – The ZZ interaction strength. Positive is antiferromagnetic.
bx (float, optional) – The X-magnetic field strength.
S ({1/2, 1, 3/2, ...}, optional) – The underlying spin of the system, defaults to 1/2.
cyclic (bool, optional) – Generate a MPO with periodic boundary conditions or not, default is open boundary conditions.
local_ham_1d_opts (mpo_opts or) – Supplied to
MatrixProductOperator
.
- Return type:
- quimb.tensor.tensor_builder.ham_1d_ising(L=None, j=1.0, bx=0.0, *, S=1 / 2, cyclic=False, **local_ham_1d_opts)¶
Ising Hamiltonian in
LocalHam1D
form.\[H_\mathrm{Ising} = J \sum_{i} S^Z_i S^Z_{i + 1} - B_x \sum_{i} S^X_i\]Note the default convention of antiferromagnetic interactions and spin operators not Pauli matrices.
- Parameters:
L (int) – The number of sites.
j (float, optional) – The ZZ interaction strength. Positive is antiferromagnetic.
bx (float, optional) – The X-magnetic field strength.
S ({1/2, 1, 3/2, ...}, optional) – The underlying spin of the system, defaults to 1/2.
cyclic (bool, optional) – Generate a hamiltonian with periodic boundary conditions or not, default is open boundary conditions.
local_ham_1d_opts (mpo_opts or) – Supplied to
LocalHam1D
.
- Return type:
- quimb.tensor.tensor_builder.NNI_ham_ising¶
- quimb.tensor.tensor_builder._ham_XY(j=1.0, bz=0.0, *, S=1 / 2, cyclic=False)¶
- quimb.tensor.tensor_builder.MPO_ham_XY(L, j=1.0, bz=0.0, *, S=1 / 2, cyclic=False, **mpo_opts)¶
XY-Hamiltonian in MPO form.
\[H_\mathrm{XY} = \sum_{i} ( J_X S^X_i S^X_{i + 1} + J_Y S^Y_i S^Y_{i + 1} ) - B_x \sum_{i} S^Z_i\]Note the default convention of antiferromagnetic interactions and spin operators not Pauli matrices.
- Parameters:
L (int) – The number of sites.
j (float or (float, float), optional) – The XX and YY interaction strength. Positive is antiferromagnetic.
bz (float, optional) – The Z-magnetic field strength.
S ({1/2, 1, 3/2, ...}, optional) – The underlying spin of the system, defaults to 1/2.
cyclic (bool, optional) – Generate a MPO with periodic boundary conditions or not, default is open boundary conditions.
local_ham_1d_opts (mpo_opts or) – Supplied to
MatrixProductOperator
.
- Return type:
- quimb.tensor.tensor_builder.ham_1d_XY(L=None, j=1.0, bz=0.0, *, S=1 / 2, cyclic=False, **local_ham_1d_opts)¶
XY-Hamiltonian in
LocalHam1D
form.\[H_\mathrm{XY} = \sum_{i} ( J_X S^X_i S^X_{i + 1} + J_Y S^Y_i S^Y_{i + 1} ) - B_Z \sum_{i} S^Z_i\]Note the default convention of antiferromagnetic interactions and spin operators not Pauli matrices.
- Parameters:
L (int) – The number of sites.
j (float or (float, float), optional) – The XX and YY interaction strength. Positive is antiferromagnetic.
bz (float, optional) – The Z-magnetic field strength.
S ({1/2, 1, 3/2, ...}, optional) – The underlying spin of the system, defaults to 1/2.
cyclic (bool, optional) – Generate a hamiltonian with periodic boundary conditions or not, default is open boundary conditions.
local_ham_1d_opts – Supplied to
LocalHam1D
.
- Return type:
- quimb.tensor.tensor_builder.NNI_ham_XY¶
- quimb.tensor.tensor_builder._ham_heis(j=1.0, bz=0.0, *, S=1 / 2, cyclic=False)¶
- quimb.tensor.tensor_builder.MPO_ham_heis(L, j=1.0, bz=0.0, *, S=1 / 2, cyclic=False, **mpo_opts)¶
Heisenberg Hamiltonian in MPO form.
\[H_\mathrm{Heis} = \sum_{i} ( J_X S^X_i S^X_{i + 1} + J_Y S^Y_i S^Y_{i + 1} + J_Z S^Z_i S^Z_{i + 1} ) - B_Z \sum_{i} S^Z_i\]Note the default convention of antiferromagnetic interactions and spin operators not Pauli matrices.
- Parameters:
L (int) – The number of sites.
j (float or (float, float, float), optional) – The XX, YY and ZZ interaction strength. Positive is antiferromagnetic.
bz (float, optional) – The Z-magnetic field strength.
S ({1/2, 1, 3/2, ...}, optional) – The underlying spin of the system, defaults to 1/2.
cyclic (bool, optional) – Generate a MPO with periodic boundary conditions or not, default is open boundary conditions.
mpo_opts – Supplied to
MatrixProductOperator
.
- Return type:
- quimb.tensor.tensor_builder.ham_1d_heis(L=None, j=1.0, bz=0.0, *, S=1 / 2, cyclic=False, **local_ham_1d_opts)¶
Heisenberg Hamiltonian in
LocalHam1D
form.\[H_\mathrm{Heis} = \sum_{i} ( J_X S^X_i S^X_{i + 1} + J_Y S^Y_i S^Y_{i + 1} + J_Z S^Z_i S^Z_{i + 1} ) - B_Z \sum_{i} S^Z_i\]Note the default convention of antiferromagnetic interactions and spin operators not Pauli matrices.
- Parameters:
L (int) – The number of sites.
j (float or (float, float, float), optional) – The XX, YY and ZZ interaction strength. Positive is antiferromagnetic.
bz (float, optional) – The Z-magnetic field strength.
S ({1/2, 1, 3/2, ...}, optional) – The underlying spin of the system, defaults to 1/2.
cyclic (bool, optional) – Generate a hamiltonian with periodic boundary conditions or not, default is open boundary conditions.
local_ham_1d_opts – Supplied to
LocalHam1D
.
- Return type:
- quimb.tensor.tensor_builder.NNI_ham_heis¶
- quimb.tensor.tensor_builder.MPO_ham_XXZ(L, delta, jxy=1.0, *, S=1 / 2, cyclic=False, **mpo_opts)¶
XXZ-Hamiltonian in MPO form.
\[H_\mathrm{XXZ} = \sum_{i} ( J_X S^X_i S^X_{i + 1} + J_Y S^Y_i S^Y_{i + 1} + \Delta S^Z_i S^Z_{i + 1} ) - B_Z \sum_{i} S^Z_i\]Note the default convention of antiferromagnetic interactions and spin operators not Pauli matrices.
- Parameters:
L (int) – The number of sites.
delta (float) – The ZZ-interaction strength. Positive is antiferromagnetic.
jxy (float, optional) – The X- and Y- interaction strength, defaults to 1. Positive is antiferromagnetic.
S ({1/2, 1, 3/2, ...}, optional) – The underlying spin of the system, defaults to 1/2.
cyclic (bool, optional) – Generate a MPO with periodic boundary conditions or not, default is open boundary conditions.
mpo_opts – Supplied to
MatrixProductOperator
.
- Return type:
- quimb.tensor.tensor_builder.ham_1d_XXZ(L=None, delta=None, jxy=1.0, *, S=1 / 2, cyclic=False, **local_ham_1d_opts)¶
XXZ-Hamiltonian in
LocalHam1D
form.\[H_\mathrm{XXZ} = \sum_{i} ( J_X S^X_i S^X_{i + 1} + J_Y S^Y_i S^Y_{i + 1} + \Delta S^Z_i S^Z_{i + 1} ) - B_Z \sum_{i} S^Z_i\]Note the default convention of antiferromagnetic interactions and spin operators not Pauli matrices.
- Parameters:
L (int) – The number of sites.
delta (float) – The ZZ-interaction strength. Positive is antiferromagnetic.
jxy (float, optional) – The X- and Y- interaction strength, defaults to 1. Positive is antiferromagnetic.
S ({1/2, 1, 3/2, ...}, optional) – The underlying spin of the system, defaults to 1/2.
cyclic (bool, optional) – Generate a hamiltonian with periodic boundary conditions or not, default is open boundary conditions.
local_ham_1d_opts – Supplied to
LocalHam1D
.
- Return type:
- quimb.tensor.tensor_builder.NNI_ham_XXZ¶
- quimb.tensor.tensor_builder._ham_bilinear_biquadratic(theta, *, S=1 / 2, cyclic=False)¶
- quimb.tensor.tensor_builder.MPO_ham_bilinear_biquadratic(L=None, theta=0, *, S=1 / 2, cyclic=False, compress=True, **mpo_opts)¶
Hamiltonian of one-dimensional bilinear biquadratic chain in MPO form, see PhysRevB.93.184428.
- Parameters:
L (int) – The number of sites.
theta (float or (float, float), optional) – The parameter for linear and non-linear term of interaction strength, defaults to 0.
S ({1/2, 1, 3/2, ...}, optional) – The underlying spin of the system, defaults to 1/2.
cyclic (bool, optional) – Generate a hamiltonian with periodic boundary conditions or not, default is open boundary conditions.
mpo_opts – Supplied to
MatrixProductOperator
.
- Return type:
- quimb.tensor.tensor_builder.ham_1d_bilinear_biquadratic(L=None, theta=0, *, S=1 / 2, cyclic=False, **local_ham_1d_opts)¶
Hamiltonian of one-dimensional bilinear biquadratic chain in LocalHam1D form, see PhysRevB.93.184428.
- Parameters:
L (int) – The number of sites.
theta (float or (float, float), optional) – The parameter for linear and non-linear term of interaction strength, defaults to 0.
S ({1/2, 1, 3/2, ...}, optional) – The underlying spin of the system, defaults to 1/2.
cyclic (bool, optional) – Generate a hamiltonian with periodic boundary conditions or not, default is open boundary conditions.
local_ham_1d_opts – Supplied to
LocalHam1D
.
- Return type:
- quimb.tensor.tensor_builder.NNI_ham_bilinear_biquadratic¶
- quimb.tensor.tensor_builder._ham_mbl(L, dh, j=1.0, seed=None, S=1 / 2, *, cyclic=False, dh_dist='s', dh_dim=1, beta=None)¶
- quimb.tensor.tensor_builder.MPO_ham_mbl(L, dh, j=1.0, seed=None, S=1 / 2, *, cyclic=False, dh_dist='s', dh_dim=1, beta=None, **mpo_opts)¶
The many-body-localized spin hamiltonian in MPO form.
\[H_\mathrm{MBL} = \sum_{i} ( J_X S^X_i S^X_{i + 1} + J_Y S^Y_i S^Y_{i + 1} + J_Z S^Z_i S^Z_{i + 1} ) - \sum_{i} h_i S^Z_i\]Note the default convention of antiferromagnetic interactions and spin operators not Pauli matrices.
- Parameters:
L (int) – Number of spins.
dh (float) – Random noise strength.
j (float, or (float, float, float), optional) – Interaction strength(s) e.g. 1 or (1., 1., 0.5). Positive is antiferromagnetic.
seed (int, optional) – Random number to seed the noise with.
S ({1/2, 1, 3/2, ...}, optional) – The underlying spin of the system, defaults to 1/2.
cyclic (bool, optional) – Whether to use periodic boundary conditions - default is False.
dh_dist ({'s', 'g', 'qp'}, optional) – Whether to use sqaure, guassian or quasiperiodic noise.
beta (float, optional) – Frequency of the quasirandom noise, only if
dh_dist='qr'
.mpo_opts – Supplied to
MatrixProductOperator
.
- Return type:
- quimb.tensor.tensor_builder.ham_1d_mbl(L, dh, j=1.0, seed=None, S=1 / 2, *, cyclic=False, dh_dist='s', dh_dim=1, beta=None, **local_ham_1d_opts)¶
The many-body-localized spin hamiltonian in
LocalHam1D
form.\[H_\mathrm{MBL} = \sum_{i} ( J_X S^X_i S^X_{i + 1} + J_Y S^Y_i S^Y_{i + 1} + J_Z S^Z_i S^Z_{i + 1} ) - \sum_{i} h_i S^Z_i\]Note the default convention of antiferromagnetic interactions and spin operators not Pauli matrices.
- Parameters:
L (int) – Number of spins.
dh (float) – Random noise strength.
j (float, or (float, float, float), optional) – Interaction strength(s) e.g. 1 or (1., 1., 0.5). Positive is antiferromagnetic.
seed (int, optional) – Random number to seed the noise with.
S ({1/2, 1, 3/2, ...}, optional) – The underlying spin of the system, defaults to 1/2.
cyclic (bool, optional) – Whether to use periodic boundary conditions - default is False.
dh_dist ({'s', 'g', 'qp'}, optional) – Whether to use sqaure, guassian or quasiperiodic noise.
beta (float, optional) – Frequency of the quasirandom noise, only if
dh_dist='qr'
.local_ham_1d_opts – Supplied to
LocalHam1D
.
- Return type:
- quimb.tensor.tensor_builder.NNI_ham_mbl¶
- quimb.tensor.tensor_builder.ham_2d_ising(Lx, Ly, j=1.0, bx=0.0, **local_ham_2d_opts)¶
Ising Hamiltonian in
LocalHam2D
form.\[H_\mathrm{Ising} = J \sum_{<ij>} S^Z_i S^Z_{j} - B_x \sum_{i} S^X_i\]for nearest neighbors \(<ij>\). Note the default convention of antiferromagnetic interactions and spin operators not Pauli matrices.
- Parameters:
- Return type:
- quimb.tensor.tensor_builder.ham_2d_heis(Lx, Ly, j=1.0, bz=0.0, **local_ham_2d_opts)¶
Heisenberg Hamiltonian in
LocalHam2D
. form.\[H_\mathrm{Heis} = \sum_{<ij>} ( J_X S^X_i S^X_{j} + J_Y S^Y_i S^Y_{j} + J_Z S^Z_i S^Z_{j} ) - B_Z \sum_{i} S^Z_{i}\]for nearest neighbors \(<ij>\). Note the default convention of antiferromagnetic interactions and spin operators not Pauli matrices.
- Parameters:
- Return type:
- quimb.tensor.tensor_builder.ham_2d_j1j2(Lx, Ly, j1=1.0, j2=0.5, bz=0.0, **local_ham_2d_opts)¶
Heisenberg Hamiltonian in
LocalHam2D
. form.\[H_\mathrm{Heis} = \sum_{<ij>} ( J_{1,X} S^X_i S^X_{j} + J_{1,Y} S^Y_i S^Y_{j} + J_{1,Z} S^Z_i S^Z_{j} ) + \sum_{<<ij>>} ( J_{2,X} S^X_i S^X_{j} + J_{2,Y} S^Y_i S^Y_{j} + J_{2,Z} S^Z_i S^Z_{j} ) - B_Z \sum_{i} S^Z_{i}\]for nearest neighbors \(<ij>\) and diagonal next nearest neighbors \(<<ij>>\). Note the default convention of antiferromagnetic interactions and spin operators not Pauli matrices.
- Parameters:
Lx (int) – The number of rows.
Ly (int) – The number of columns.
j2 (float or (float, float, float), optional) – The nearest neighbor XX, YY and ZZ interaction strength. Positive is antiferromagnetic.
j2 – The diagonal next nearest nearest XX, YY and ZZ interaction strength. Positive is antiferromagnetic.
bz (float, optional) – The Z-magnetic field strength.
local_ham_2d_opts – Supplied to
LocalHam2D
.
- Return type:
- quimb.tensor.tensor_builder.ham_3d_heis(Lx, Ly, Lz, j=1.0, bz=0.0, **local_ham_3d_opts)¶
Heisenberg Hamiltonian in
LocalHam3D
. form.\[H_\mathrm{Heis} = \sum_{<ij>} ( J_X S^X_i S^X_{j} + J_Y S^Y_i S^Y_{j} + J_Z S^Z_i S^Z_{j} ) - B_Z \sum_{i} S^Z_{i}\]for nearest neighbors \(<ij>\). Note the default convention of antiferromagnetic interactions and spin operators not Pauli matrices.
- Parameters:
Lx (int) – The number of x-planes.
Ly (int) – The number of y-planes.
Ly – The number of z-planes.
j (float or (float, float, float), optional) – The XX, YY and ZZ interaction strength. Positive is antiferromagnetic.
bz (float, optional) – The Z-magnetic field strength.
local_ham_3d_opts – Supplied to
LocalHam3D
.
- Return type: