MultiSpeciesLattice

  • full name: tenpy.models.lattice.MultiSpeciesLattice

  • parent module: tenpy.models.lattice

  • type: class

Inheritance Diagram

Inheritance diagram of tenpy.models.lattice.MultiSpeciesLattice

Methods

MultiSpeciesLattice.__init__(simple_lattice, ...)

MultiSpeciesLattice.copy()

Shallow copy of self.

MultiSpeciesLattice.count_neighbors([u, key])

Count e.g. the number of nearest neighbors for a site in the bulk.

MultiSpeciesLattice.coupling_shape(dx)

Calculate correct shape of the strengths for a coupling.

MultiSpeciesLattice.distance(u1, u2, dx)

Get the distance for a given coupling between two sites in the lattice.

MultiSpeciesLattice.enlarge_mps_unit_cell([...])

Repeat the unit cell for infinite MPS boundary conditions; in place.

MultiSpeciesLattice.extract_segment([first, ...])

Extract a finite segment from an infinite/large system.

MultiSpeciesLattice.find_coupling_pairs([...])

Automatically find coupling pairs grouped by distances.

MultiSpeciesLattice.from_hdf5(hdf5_loader, ...)

Load instance from a HDF5 file.

MultiSpeciesLattice.lat2mps_idx(lat_idx)

Translate lattice indices (x_0, ..., x_{D-1}, u) to MPS index i.

MultiSpeciesLattice.mps2lat_idx(i)

Translate MPS index i to lattice indices (x_0, ..., x_{dim-1}, u).

MultiSpeciesLattice.mps2lat_values(A[, axes, u])

Reshape/reorder A to replace an MPS index by lattice indices.

MultiSpeciesLattice.mps2lat_values_masked(A)

Reshape/reorder an array A to replace an MPS index by lattice indices.

MultiSpeciesLattice.mps_idx_fix_u([u])

return an index array of MPS indices for which the site within the unit cell is u.

MultiSpeciesLattice.mps_lat_idx_fix_u([u])

Similar as mps_idx_fix_u(), but return also the corresponding lattice indices.

MultiSpeciesLattice.mps_sites()

Return a list of sites for all MPS indices.

MultiSpeciesLattice.multi_coupling_shape(dx)

Calculate correct shape of the strengths for a multi_coupling.

MultiSpeciesLattice.number_nearest_neighbors([u])

Deprecated.

MultiSpeciesLattice.number_next_nearest_neighbors([u])

Deprecated.

MultiSpeciesLattice.ordering(order)

Define orderings as for the simple_lattice with priority for within the unit cell.

MultiSpeciesLattice.plot_basis(ax[, origin, ...])

Plot arrows indicating the basis vectors of the lattice.

MultiSpeciesLattice.plot_bc_identified(ax[, ...])

Mark two sites identified by periodic boundary conditions.

MultiSpeciesLattice.plot_coupling(ax[, ...])

Plot lines connecting nearest neighbors of the lattice.

MultiSpeciesLattice.plot_order(ax[, order, ...])

Plot a line connecting sites in the specified "order" and text labels enumerating them.

MultiSpeciesLattice.plot_sites(ax[, ...])

Plot the sites of the lattice with markers.

MultiSpeciesLattice.position(lat_idx)

return 'space' position of one or multiple sites.

MultiSpeciesLattice.possible_couplings(u1, ...)

Find possible MPS indices for two-site couplings.

MultiSpeciesLattice.possible_multi_couplings(ops)

Generalization of possible_couplings() to couplings with more than 2 sites.

MultiSpeciesLattice.save_hdf5(hdf5_saver, ...)

Export self into a HDF5 file.

MultiSpeciesLattice.self_u_to_simple_u(self_u)

Get index u of the simple_lattice from index u in self.

MultiSpeciesLattice.self_u_to_species_idx(self_u)

Get species index for unit cell index.

MultiSpeciesLattice.simple_u_to_species_u(...)

Get index u in self from the u in the simple_lattice and the species index.

MultiSpeciesLattice.site(i)

return Site instance corresponding to an MPS index i

MultiSpeciesLattice.test_sanity()

Sanity check.

Class Attributes and Properties

MultiSpeciesLattice.Lu

unknown number of sites in the unit cell

MultiSpeciesLattice.boundary_conditions

Human-readable list of boundary conditions from bc and bc_shift.

MultiSpeciesLattice.cylinder_axis

Direction of the cylinder axis.

MultiSpeciesLattice.dim

The dimension of the lattice.

MultiSpeciesLattice.nearest_neighbors

MultiSpeciesLattice.next_nearest_neighbors

MultiSpeciesLattice.next_next_nearest_neighbors

MultiSpeciesLattice.order

Defines an ordering of the lattice sites, thus mapping the lattice to a 1D chain.

class tenpy.models.lattice.MultiSpeciesLattice(simple_lattice, species_sites, species_names=None)[source]

Bases: Lattice

A variant of a SimpleLattice replacing the elementary site with a set of sites.

An initialized MultiSpeciesLattice replaces each site in the given simple_lattice with the species_sites. This is useful e.g. if you want to place spin-full fermions (or bosons) on a regular lattice, without falling back to the GroupedSite causing a larger, local dimension.

This class defines pairs with appended '_all' and the species_names for all existing pairs in the simple_lattice; see the example below.

Parameters:
  • simple_lattice (Lattice) – A regular lattice (not necessarily a SimpleLattice). Typically one of the predefined ones, e.g., Square, Triangular or Honeycomb. You can initialize the simple_lattice with any local sites, e.g. sites=None, since the sites of the MultiSpeciesLattice are defined below.

  • species_sites (list of Site) – A collection of sites representing different species.

  • species_names (list of str | None) – Names for the species. Defaults to ['0', '1', ...]. Used to define separate pairs by appending '_{name}' to the lattice.

N_species

Number of species

Type:

int

species_names

Names for the species.

Type:

list of str

simple_lattice

The regular lattice on which self is based.

Type:

Lattice

simple_Lu

Length of the unit cell of the simple_lattice.

Type:

int

Examples

When defining the sites, you should probably call set_common_charges() (see examples there!) to adjust the charges, especially if you have different sites. Consider a combination of Fermions and Spins:

>>> simple_lat = Square(2, 3, None)
>>> f = tenpy.networks.site.FermionSite(conserve='N')
>>> s = tenpy.networks.site.SpinHalfSite(conserve='Sz')
>>> tenpy.networks.site.set_common_charges([f, s], 'independent')
[array([0, 1]), array([1, 0])]
>>> fs_lat = MultiSpeciesLattice(simple_lat, [f, s], ['f', 's'])

There are corresponding coupling pairs defined:

>>> for key in fs_lat.pairs.keys():
...     if key.startswith('nearest'):
...          print(key)
nearest_neighbors_f-f
nearest_neighbors_f-s
nearest_neighbors_s-f
nearest_neighbors_s-s
nearest_neighbors_all-all
nearest_neighbors_diag
>>> print(fs_lat.pairs['nearest_neighbors_diag'])
[(0, 0, array([1, 0])), (0, 0, array([0, 1])), (1, 1, array([1, 0])), (1, 1, array([0, 1]))]

We further have “onsite” terms for couplings between the species defined. Here, there is no 'onsite_s-f', as it would be the same as 'onsite_f-s'.

>>> for key in fs_lat.pairs.keys():
...     if key.startswith('onsite'):
...          print(key)
onsite_f-s

Note that the “simple lattice” can also have a non-trivial unit cell itself, e.g. the Honeycomb already has two sites in its unit cell:

>>> simple_lat = Honeycomb(2, 3, None)
>>> f = tenpy.networks.site.FermionSite(conserve='N')
>>> sites = [f, copy(f)]
>>> tenpy.networks.site.set_common_charges(sites, 'same')  # same = total N conserved
[array([0, 1]), array([0, 1])]
>>> spinful_fermion_Honeycomb = MultiSpeciesLattice(simple_lat, sites, ['up', 'down'])

In this case, you could also call tenpy.networks.site.spin_half_species().

Lu = None

unknown number of sites in the unit cell

ordering(order)[source]

Define orderings as for the simple_lattice with priority for within the unit cell.

See Lattice.ordering() for arguments.

self_u_to_simple_u(self_u)[source]

Get index u of the simple_lattice from index u in self.

Parameters:

self_u (int) – Unit cell index u in self.

Returns:

simple_u – Unit cell index u in the simple_lattice.

Return type:

int

self_u_to_species_idx(self_u)[source]

Get species index for unit cell index.

Parameters:

self_u (int) – Unit cell index u in self.

Returns:

simple_u – Unit cell index u in the simple_lattice.

Return type:

int

simple_u_to_species_u(simple_u, species_idx)[source]

Get index u in self from the u in the simple_lattice and the species index.

Parameters:
  • simple_u (int) – Unit cell index u in the simple_lattice.

  • species_idx (int) – Index of the species.

Returns:

self_u – Unit cell index u in self.

Return type:

int

property boundary_conditions

Human-readable list of boundary conditions from bc and bc_shift.

Returns:

boundary_conditions – List of "open" or "periodic", one entry for each direction of the lattice.

Return type:

list of str

copy()[source]

Shallow copy of self.

count_neighbors(u=0, key='nearest_neighbors')[source]

Count e.g. the number of nearest neighbors for a site in the bulk.

Parameters:
  • u (int) – Specifies the site in the unit cell, for which we should count the number of neighbors (or whatever key specifies).

  • key (str) – Key of pairs to select what to count.

Returns:

number – Number of nearest neighbors (or whatever key specified) for the u-th site in the unit cell, somewhere in the bulk of the lattice. Note that it might not be the correct value at the edges of a lattice with open boundary conditions.

Return type:

int

coupling_shape(dx)[source]

Calculate correct shape of the strengths for a coupling.

Parameters:

dx (tuple of int) – Translation vector in the lattice for a coupling of two operators. Corresponds to dx argument of tenpy.models.model.CouplingModel.add_multi_coupling().

Returns:

  • coupling_shape (tuple of int) – Len dim. The correct shape for an array specifying the coupling strength. lat_indices has only rows within this shape.

  • shift_lat_indices (array) – Translation vector from origin to the lower left corner of box spanned by dx.

property cylinder_axis

Direction of the cylinder axis.

For an infinite cylinder (bc_MPS=’infinite’ and ``boundary_conditions[1] == ‘open’`), this property gives the direction of the cylinder axis, in the same coordinates as the basis, as a normalized vector. For a 1D lattice or for open boundary conditions along y, it’s just along basis[0].

property dim

The dimension of the lattice.

distance(u1, u2, dx)[source]

Get the distance for a given coupling between two sites in the lattice.

The u1, u2, dx parameters are defined in analogy with add_coupling(), i.e., this function calculates the distance between a pair of operators added with add_coupling (using the basis and unit_cell_positions of the lattice).

Warning

This function ignores “wrapping” around the cylinder in the case of periodic boundary conditions.

Parameters:
  • u1 (int) – Indices within the unit cell; the u1 and u2 of add_coupling()

  • u2 (int) – Indices within the unit cell; the u1 and u2 of add_coupling()

  • dx (array) – Length dim. The translation in terms of basis vectors for the coupling.

Returns:

distance – The distance between site at lattice indices [x, y, u1] and [x + dx[0], y + dx[1], u2], ignoring any boundary effects. In case of non-trivial position_disorder, an array is returned. This array is compatible with the shape/indexing required for add_coupling(). For example to add a Z-Z interaction of strength J/r with r the distance, you can do something like this in init_terms():

for u1, u2, dx in self.lat.pairs[‘nearest_neighbors’]:

dist = self.lat.distance(u1, u2, dx) self.add_coupling(J/dist, u1, ‘Sz’, u2, ‘Sz’, dx)

Return type:

float | ndarray

enlarge_mps_unit_cell(factor=2)[source]

Repeat the unit cell for infinite MPS boundary conditions; in place.

Parameters:

factor (int) – The new number of sites in the MPS unit cell will be increased from N_sites to factor*N_sites_per_ring. Since MPS unit cells are repeated in the x-direction in our convention, the lattice shape goes from (Lx, Ly, ..., Lu) to (Lx*factor, Ly, ..., Lu).

extract_segment(first=0, last=None, enlarge=None)[source]

Extract a finite segment from an infinite/large system.

Parameters:
  • first (int) – The first and last site to include into the segment. last defaults to L - 1, i.e., the MPS unit cell for infinite MPS.

  • last (int) – The first and last site to include into the segment. last defaults to L - 1, i.e., the MPS unit cell for infinite MPS.

  • enlarge (int) – Instead of specifying the first and last site, you can specify this factor by how much the MPS unit cell should be enlarged.

Returns:

copy – A copy of self with “segment” bc_MPS and segment_first_last set.

Return type:

Lattice

find_coupling_pairs(max_dx=3, cutoff=None, eps=1e-10)[source]

Automatically find coupling pairs grouped by distances.

Given the unit_cell_positions and basis, the coupling pairs of nearest_neighbors, next_nearest_neighbors etc at a given distance are basically fixed (although not uniquely, since we take out half of them to avoid double-counting couplings in both directions A_i B_j + B_i A_i). This function iterates through all possible couplings up to a given cutoff distance and then determines the possible pairs at fixed distances (up to round-off errors).

Parameters:
  • max_dx (int) – Maximal index for each index of dx to iterate over. You need large enough values to include every possible coupling up to the desired distance, but choosing it too large might make this function run for a long time.

  • cutoff (float) – Maximal distance (in the units in which basis and unit_cell_positions is given).

  • eps (float) – Tolerance up to which to distances are considered the same.

Returns:

coupling_pairs – Keys are distances of nearest-neighbors, next-nearest-neighbors etc. Values are [(u1, u2, dx), ...] as in pairs.

Return type:

dict

classmethod from_hdf5(hdf5_loader, h5gr, subpath)[source]

Load instance from a HDF5 file.

This method reconstructs a class instance from the data saved with save_hdf5().

Parameters:
  • hdf5_loader (Hdf5Loader) – Instance of the loading engine.

  • h5gr (Group) – HDF5 group which is representing the object to be constructed.

  • subpath (str) – The name of h5gr with a '/' in the end.

Returns:

obj – Newly generated class instance containing the required data.

Return type:

cls

lat2mps_idx(lat_idx)[source]

Translate lattice indices (x_0, ..., x_{D-1}, u) to MPS index i.

Parameters:

lat_idx (array_like [..., dim+1]) – The last dimension corresponds to lattice indices (x_0, ..., x_{D-1}, u). All lattice indices should be positive and smaller than the corresponding entry in self.shape. Exception: for “infinite” or “segment” bc_MPS, an x_0 outside indicates shifts across the boundary.

Returns:

i – MPS index/indices corresponding to lat_idx. Has the same shape as lat_idx without the last dimension.

Return type:

array_like

mps2lat_idx(i)[source]

Translate MPS index i to lattice indices (x_0, ..., x_{dim-1}, u).

Parameters:

i (int | array_like of int) – MPS index/indices.

Returns:

lat_idx – First dimensions like i, last dimension has len dim`+1 and contains the lattice indices ``(x_0, …, x_{dim-1}, u)` corresponding to i. For i across the MPS unit cell and “infinite” or “segment” bc_MPS, we shift x_0 accordingly.

Return type:

array

mps2lat_values(A, axes=0, u=None)[source]

Reshape/reorder A to replace an MPS index by lattice indices.

Parameters:
  • A (ndarray) – Some values. Must have A.shape[axes] = self.N_sites if u is None, or A.shape[axes] = self.N_cells if u is an int.

  • axes ((iterable of) int) – chooses the axis which should be replaced.

  • u (None | int) – Optionally choose a subset of MPS indices present in the axes of A, namely the indices corresponding to self.unit_cell[u], as returned by mps_idx_fix_u(). The resulting array will not have the additional dimension(s) of u.

Returns:

res_A – Reshaped and reordered version of A. Such that MPS indices along the specified axes are replaced by lattice indices, i.e., if MPS index j maps to lattice site (x0, x1, x2), then res_A[..., x0, x1, x2, ...] = A[..., j, ...].

Return type:

ndarray

Examples

Say you measure expectation values of an onsite term for an MPS, which gives you an 1D array A, where A[i] is the expectation value of the site given by self.mps2lat_idx(i). Then this function gives you the expectation values ordered by the lattice:

>>> print(lat.shape, A.shape)
(10, 3, 2) (60,)
>>> A_res = lat.mps2lat_values(A)
>>> A_res.shape
(10, 3, 2)
>>> A_res[tuple(lat.mps2lat_idx(5))] == A[5]
True

If you have a correlation function C[i, j], it gets just slightly more complicated:

>>> print(lat.shape, C.shape)
(10, 3, 2) (60, 60)
>>> lat.mps2lat_values(C, axes=[0, 1]).shape
(10, 3, 2, 10, 3, 2)

If the unit cell consists of different physical sites, an onsite operator might be defined only on one of the sites in the unit cell. Then you can use mps_idx_fix_u() to get the indices of sites it is defined on, measure the operator on these sites, and use the argument u of this function.

>>> u = 0
>>> idx_subset = lat.mps_idx_fix_u(u)
>>> A_u = A[idx_subset]
>>> A_u_res = lat.mps2lat_values(A_u, u=u)
>>> A_u_res.shape
(10, 3)
>>> np.all(A_res[:, :, u] == A_u_res[:, :])
True
mps2lat_values_masked(A, axes=-1, mps_inds=None, include_u=None)[source]

Reshape/reorder an array A to replace an MPS index by lattice indices.

This is a generalization of mps2lat_values() allowing for the case of an arbitrary set of MPS indices present in each axis of A.

Parameters:
  • A (ndarray) – Some values.

  • axes ((iterable of) int) – Chooses the axis of A which should be replaced. If multiple axes are given, you also need to give multiple index arrays as mps_inds.

  • mps_inds ((list of) 1D ndarray) – Specifies for each axis in axes, for which MPS indices we have values in the corresponding axis of A. Defaults to [np.arange(A.shape[ax]) for ax in axes]. For indices across the MPS unit cell and “infinite” bc_MPS, we shift x_0 accordingly.

  • include_u ((list of) bool) – Specifies for each axis in axes, whether the u index of the lattice should be included into the output array res_A. Defaults to len(self.unit_cell) > 1.

Returns:

res_A – Reshaped and reordered copy of A. Such that MPS indices along the specified axes are replaced by lattice indices, i.e., if MPS index j maps to lattice site (x0, x1, x2), then res_A[..., x0, x1, x2, ...] = A[..., mps_inds[j], ...].

Return type:

np.ma.MaskedArray

mps_idx_fix_u(u=None)[source]

return an index array of MPS indices for which the site within the unit cell is u.

If you have multiple sites in your unit-cell, an onsite operator is in general not defined for all sites. This functions returns an index array of the mps indices which belong to sites given by self.unit_cell[u].

Parameters:

u (None | int) – Selects a site of the unit cell. None (default) means all sites.

Returns:

mps_idx – MPS indices for which self.site(i) is self.unit_cell[u]. Ordered ascending.

Return type:

array

mps_lat_idx_fix_u(u=None)[source]

Similar as mps_idx_fix_u(), but return also the corresponding lattice indices.

Parameters:

u (None | int) – Selects a site of the unit cell. None (default) means all sites.

Returns:

  • mps_idx (array) – MPS indices i for which self.site(i) is self.unit_cell[u].

  • lat_idx (2D array) – The row j contains the lattice index (without u) corresponding to mps_idx[j].

mps_sites()[source]

Return a list of sites for all MPS indices.

Equivalent to [self.site(i) for i in range(self.N_sites)].

This should be used for sites of 1D tensor networks (MPS, MPO,…).

multi_coupling_shape(dx)[source]

Calculate correct shape of the strengths for a multi_coupling.

Parameters:

dx (2D array, shape (N_ops, dim)) – dx[i, :] is the translation vector in the lattice for the i-th operator. Corresponds to the dx of each operator given in the argument ops of tenpy.models.model.CouplingModel.add_multi_coupling().

Returns:

  • coupling_shape (tuple of int) – Len dim. The correct shape for an array specifying the coupling strength. lat_indices has only rows within this shape.

  • shift_lat_indices (array) – Translation vector from origin to the lower left corner of box spanned by dx. (Unlike for coupling_shape() it can also contain entries > 0)

number_nearest_neighbors(u=0)[source]

Deprecated.

Deprecated since version 0.5.0: Use count_neighbors() instead.

number_next_nearest_neighbors(u=0)[source]

Deprecated.

Deprecated since version 0.5.0: Use count_neighbors() instead.

property order

Defines an ordering of the lattice sites, thus mapping the lattice to a 1D chain.

Each row of the array contains the lattice indices for one site, the order of the rows thus specifies a path through the lattice, along which an MPS will wind through the lattice.

You can visualize the order with plot_order().

plot_basis(ax, origin=(0.0, 0.0), shade=None, **kwargs)[source]

Plot arrows indicating the basis vectors of the lattice.

Parameters:
  • ax (matplotlib.axes.Axes) – The axes on which we should plot.

  • **kwargs – Keyword arguments for ax.arrow.

plot_bc_identified(ax, direction=-1, origin=None, cylinder_axis=False, **kwargs)[source]

Mark two sites identified by periodic boundary conditions.

Works only for lattice with a 2-dimensional basis.

Parameters:
  • ax (matplotlib.axes.Axes) – The axes on which we should plot.

  • direction (int) – The direction of the lattice along which we should mark the identified sites. If None, mark it along all directions with periodic boundary conditions.

  • cylinder_axis (bool) – Whether to plot the cylinder axis as well.

  • origin (None | np.ndarray) – The origin starting from where we mark the identified sites. Defaults to the first entry of unit_cell_positions.

  • **kwargs – Keyword arguments for the used ax.plot.

plot_coupling(ax, coupling=None, wrap=False, **kwargs)[source]

Plot lines connecting nearest neighbors of the lattice.

Parameters:
  • ax (matplotlib.axes.Axes) – The axes on which we should plot.

  • coupling (list of (u1, u2, dx)) – By default (None), use self.pairs['nearest_neighbors']. Specifies the connections to be plotted; iterating over lattice indices (i0, i1, …), we plot a connection from the site (i0, i1, ..., u1) to the site (i0+dx[0], i1+dx[1], ..., u2), taking into account the boundary conditions.

  • wrap (bool) – If True, plot couplings going around the boundary by directly connecting the sites it connects. This might be hard to see, as this puts lines from one end of the lattice to the other. If False, plot the couplings as dangling lines.

  • **kwargs – Further keyword arguments given to ax.plot().

plot_order(ax, order=None, textkwargs={'color': 'r'}, **kwargs)[source]

Plot a line connecting sites in the specified “order” and text labels enumerating them.

Parameters:
  • ax (matplotlib.axes.Axes) – The axes on which we should plot.

  • order (None | 2D array (self.N_sites, self.dim+1)) – The order as returned by ordering(); by default (None) use order.

  • textkwargs (None | dict) – If not None, we add text labels enumerating the sites in the plot. The dictionary can contain keyword arguments for ax.text().

  • **kwargs – Further keyword arguments given to ax.plot().

plot_sites(ax, markers=['o', '^', 's', 'p', 'h', 'D'], labels=None, **kwargs)[source]

Plot the sites of the lattice with markers.

Parameters:
  • ax (matplotlib.axes.Axes) – The axes on which we should plot.

  • markers (list) – List of values for the keyword marker of ax.plot() to distinguish the different sites in the unit cell, a site u in the unit cell is plotted with a marker markers[u % len(markers)].

  • labels (list of str) – Labels for the different sites in the unit cell.

  • **kwargs – Further keyword arguments given to ax.plot().

position(lat_idx)[source]

return ‘space’ position of one or multiple sites.

Parameters:

lat_idx (ndarray, (... , dim+1)) – Lattice indices.

Returns:

pos – The position of the lattice sites specified by lat_idx in real-space. If position_disorder is non-trivial, it can shift the positions!

Return type:

ndarray, (..., Dim)

possible_couplings(u1, u2, dx, strength=None)[source]

Find possible MPS indices for two-site couplings.

For periodic boundary conditions (bc[a] == False) the index x_a is taken modulo Ls[a] and runs through range(Ls[a]). For open boundary conditions, x_a is limited to 0 <= x_a < Ls[a] and 0 <= x_a+dx[a] < lat.Ls[a].

Parameters:
  • u1 (int) – Indices within the unit cell; the u1 and u2 of add_coupling()

  • u2 (int) – Indices within the unit cell; the u1 and u2 of add_coupling()

  • dx (array) – Length dim. The translation in terms of basis vectors for the coupling.

  • strength (array_like | None) – If given, instead of returning lat_indices and coupling_shape directly return the correct strength_12.

Returns:

  • mps1, mps2 (1D array) – For each possible two-site coupling the MPS indices for the u1 and u2.

  • strength_vals (1D array) – (Only returned if strength is not None.) Such that for (i, j, s) in zip(mps1, mps2, strength_vals): iterates over all possible couplings with s being the strength of that coupling. Couplings where strength_vals == 0. are filtered out.

  • lat_indices (2D int array) – (Only returned if strength is None.) Rows of lat_indices correspond to entries of mps1 and mps2 and contain the lattice indices of the “lower left corner” of the box containing the coupling.

  • coupling_shape (tuple of int) – (Only returned if strength is None.) Len dim. The correct shape for an array specifying the coupling strength. lat_indices has only rows within this shape.

possible_multi_couplings(ops, strength=None)[source]

Generalization of possible_couplings() to couplings with more than 2 sites.

Parameters:

ops (list of (opname, dx, u)) – Same as the argument ops of add_multi_coupling().

Returns:

  • mps_ijkl (2D int array) – Each row contains MPS indices i,j,k,l,…` for each of the operators positions. The positions are defined by dx (j,k,l,… relative to i) and boundary conditions of self (how much the box for given dx can be shifted around without hitting a boundary - these are the different rows).

  • strength_vals (1D array) – (Only returned if strength is not None.) Such that for  (ijkl, s) in zip(mps_ijkl, strength_vals): iterates over all possible couplings with s being the strength of that coupling. Couplings where strength_vals == 0. are filtered out.

  • lat_indices (2D int array) – (Only returned if strength is None.) Rows of lat_indices correspond to rows of mps_ijkl and contain the lattice indices of the “lower left corner” of the box containing the coupling.

  • coupling_shape (tuple of int) – (Only returned if strength is None.) Len dim. The correct shape for an array specifying the coupling strength. lat_indices has only rows within this shape.

save_hdf5(hdf5_saver, h5gr, subpath)[source]

Export self into a HDF5 file.

This method saves all the data it needs to reconstruct self with from_hdf5().

Specifically, it saves unit_cell, Ls, unit_cell_positions, basis, boundary_conditions, pairs under their name, bc_MPS as "boundary_conditions_MPS", and order as "order_for_MPS". Moreover, it saves dim and N_sites as HDF5 attributes.

Parameters:
  • hdf5_saver (Hdf5Saver) – Instance of the saving engine.

  • h5gr (:class`Group`) – HDF5 group which is supposed to represent self.

  • subpath (str) – The name of h5gr with a '/' in the end.

site(i)[source]

return Site instance corresponding to an MPS index i

test_sanity()[source]

Sanity check.

Raises ValueErrors, if something is wrong.