MPS

Inheritance Diagram

Inheritance diagram of tenpy.networks.mps.MPS

Methods

MPS.__init__(sites, Bs, SVs[, bc, form, norm])

MPS.add(other, alpha, beta[, cutoff])

Return an MPS which represents alpha|self> + beta |others>.

MPS.apply_local_op(i, op[, unitary, ...])

Apply a local (one or multi-site) operator to self.

MPS.apply_product_op(ops[, unitary, renormalize])

Apply a (global) product of local onsite operators to self.

MPS.average_charge([bond])

Return the average charge for the block on the left of a given bond.

MPS.canonical_form(**kwargs)

Bring self into canonical 'B' form, (re-)calculate singular values.

MPS.canonical_form_finite([renormalize, ...])

Bring a finite (or segment) MPS into canonical form (in place).

MPS.canonical_form_infinite(**kwargs)

Deprecated wrapper around canonical_form_infinite1().

MPS.canonical_form_infinite1([renormalize, ...])

Bring an infinite MPS into canonical form (in place).

MPS.canonical_form_infinite2([renormalize, ...])

Convert infinite MPS to canonical form.

MPS.charge_variance([bond])

Return the charge variance on the left of a given bond.

MPS.compress(options)

Compress an MPS.

MPS.compress_svd(trunc_par)

Compress self with a single sweep of SVDs; in place.

MPS.compute_K(perm[, swap_op, trunc_par, ...])

Compute the momentum quantum numbers of the entanglement spectrum for 2D states.

MPS.convert_form([new_form])

Transform self into different canonical form (by scaling the legs with singular values).

MPS.copy()

Returns a copy of self.

MPS.correlation_function(ops1, ops2[, ...])

Correlation function <bra|op1_i op2_j|ket> of single site operators, sandwiched between bra and ket. For examples the contraction for a two-site operator on site i would look like::.

MPS.correlation_length([target, tol_ev0, ...])

Calculate the correlation length by diagonalizing the transfer matrix.

MPS.correlation_length_charge_sectors([...])

Return possible charge_sector argument for correlation_length().

MPS.enlarge_chi(extra_legs[, random_fct])

Artificially enlarge the bond dimension by the specified extra legs/charges.

MPS.enlarge_mps_unit_cell([factor])

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

MPS.entanglement_entropy([n, bonds, ...])

Calculate the (half-chain) entanglement entropy for all nontrivial bonds.

MPS.entanglement_entropy_segment([segment, ...])

Calculate entanglement entropy for general geometry of the bipartition.

MPS.entanglement_entropy_segment2(segment[, n])

Calculate entanglement entropy for general geometry of the bipartition.

MPS.entanglement_spectrum([by_charge])

return entanglement energy spectrum.

MPS.expectation_value(ops[, sites, axes])

Expectation value <bra|ops|ket> of (n-site) operator(s).

MPS.expectation_value_multi_sites(operators, i0)

Expectation value <bra|op0_{i0}op1_{i0+1}...opN_{i0+N}|ket>.

MPS.expectation_value_term(term[, autoJW])

Expectation value <bra|op_{i0}op_{i1}...op_{iN}|ket>.

MPS.expectation_value_terms_sum(term_list[, ...])

Calculate expectation values for a bunch of terms and sum them up.

MPS.extract_segment(first, last)

Extract an segment from a finite or infinite MPS.

MPS.from_Bflat(sites, Bflat[, SVs, bc, ...])

Construct a matrix product state from a set of numpy arrays Bflat and singular vals.

MPS.from_full(sites, psi[, form, cutoff, ...])

Construct an MPS from a single tensor psi with one leg per physical site.

MPS.from_hdf5(hdf5_loader, h5gr, subpath)

Load instance from a HDF5 file.

MPS.from_lat_product_state(lat, p_state[, ...])

Construct an MPS from a product state given in lattice coordinates.

MPS.from_product_mps_covering(mps_covering, ...)

Create an MPS as a product of (many) local mps covering all sites to be created.

MPS.from_product_state(sites, p_state[, bc, ...])

Construct a matrix product state from a given product state.

MPS.from_singlets(site, L, pairs[, up, ...])

Create an MPS of entangled singlets.

MPS.gauge_total_charge([qtotal, vL_leg, vR_leg])

Gauge the legcharges of the virtual bonds such that the MPS has a total qtotal.

MPS.get_B(i[, form, copy, cutoff, label_p])

Return (view of) B at site i in canonical form.

MPS.get_SL(i)

Return singular values on the left of site i

MPS.get_SR(i)

Return singular values on the right of site i

MPS.get_grouped_mps(blocklen)

Like group_sites(), but make a copy.

MPS.get_op(op_list, i)

Given a list of operators, select the one corresponding to site i.

MPS.get_rho_segment(segment)

Return reduced density matrix for a segment.

MPS.get_theta(i[, n, cutoff, formL, formR])

Calculates the n-site wavefunction on sites[i:i+n].

MPS.get_total_charge([only_physical_legs])

Calculate and return the qtotal of the whole MPS (when contracted).

MPS.group_sites([n, grouped_sites])

Modify self inplace to group sites.

MPS.group_split([trunc_par])

Modify self inplace to split previously grouped sites.

MPS.increase_L([new_L])

Modify self inplace to enlarge the MPS unit cell; in place.

MPS.mutinf_two_site([max_range, n])

Calculate the two-site mutual information \(I(i:j)\).

MPS.norm_test()

Check that self is in canonical form.

MPS.overlap(other[, charge_sector, ...])

Compute overlap <self|other>.

MPS.overlap_translate_finite(psi[, shift])

Contract <self|T^N|psi> for translation T with finite, periodic boundaries.

MPS.permute_sites(perm[, swap_op, ...])

Applies the permutation perm to the state (inplace).

MPS.perturb([randomize_params, close_1, ...])

Locally perturb the state a little bit; in place.

MPS.probability_per_charge([bond])

Return probabilities of charge value on the left of a given bond.

MPS.roll_mps_unit_cell([shift])

Shift the section we define as unit cell of an infinite MPS; in place.

MPS.sample_measurements([first_site, ...])

Sample measurement results in the computational basis.

MPS.save_hdf5(hdf5_saver, h5gr, subpath)

Export self into a HDF5 file.

MPS.set_B(i, B[, form])

Set B at site i.

MPS.set_SL(i, S)

Set singular values on the left of site i

MPS.set_SR(i, S)

Set singular values on the right of site i

MPS.set_svd_theta(i, theta[, trunc_par, ...])

SVD a two-site wave function theta and save it in self.

MPS.spatial_inversion()

Perform a spatial inversion along the MPS.

MPS.swap_sites(i[, swap_op, trunc_par])

Swap the two neighboring sites i and i+1 (inplace).

MPS.term_correlation_function_left(term_L, ...)

Correlation function between (multi-site) terms, moving the left term, fix right term.

MPS.term_correlation_function_right(term_L, ...)

Correlation function between (multi-site) terms, moving the right term, fix left term.

MPS.term_list_correlation_function_right(...)

Correlation function between sums of multi-site terms, moving the right sum of term.

MPS.test_sanity()

Sanity check, raises ValueErrors, if something is wrong.

Class Attributes and Properties

MPS.L

Number of physical sites; for an iMPS the len of the MPS unit cell.

MPS.chi

Dimensions of the (nontrivial) virtual bonds.

MPS.dim

List of local physical dimensions.

MPS.finite

Distinguish MPS vs iMPS.

MPS.nontrivial_bonds

Slice of the non-trivial bond indices, depending on self.bc.

class tenpy.networks.mps.MPS(sites, Bs, SVs, bc='finite', form='B', norm=1.0)[source]

Bases: BaseMPSExpectationValue

A Matrix Product State, finite (MPS) or infinite (iMPS).

Parameters:
  • sites (list of Site) – Defines the local Hilbert space for each site.

  • Bs (list of Array) – The ‘matrices’ of the MPS. Labels are vL, vR, p (in any order).

  • SVs (list of 1D array) – The singular values on each bond. Should always have length L+1. Entries out of nontrivial_bonds are ignored.

  • bc ('finite' | 'segment' | 'infinite') – Boundary conditions as described in the table of the module doc-string.

  • form ((list of) {'B' | 'A' | 'C' | 'G' | 'Th' | None | tuple(float, float)}) – The form of the stored ‘matrices’, see table in module doc-string. A single choice holds for all of the entries.

sites

Defines the local Hilbert space for each site.

Type:

list of Site

bc

Boundary conditions as described in above table.

Type:

{‘finite’, ‘segment’, ‘infinite’}

form

Describes the canonical form on each site. None means non-canonical form. For form = (nuL, nuR), the stored _B[i] are s**form[0] -- Gamma -- s**form[1] (in Vidal’s notation).

Type:

list of {None, tuple(float, float)}

chinfo

The nature of the charge.

Type:

ChargeInfo

dtype

The data type of the _B.

Type:

type

norm

The overall norm of the state, i.e. sqrt(<psi|psi>) - the tensors are kept normalized. Ignored for (normalized) expectation_value() and co, but important for overlap() and expectation value methods of MPSEnvironment.

Type:

float

grouped

Number of sites grouped together, see group_sites().

Type:

int

segment_boundaries

Only defined for ‘segment’ bc if canonical_form_finite() has been called. If defined, it contains the U_L and V_R that would be returned by that function.

Type:

tuple of Array | (None, None)

_B

The ‘matrices’ of the MPS. Labels are vL, vR, p (in any order). We recommend using get_B() and set_B(), which will take care of the different canonical forms.

Type:

list of npc.Array

_S

The singular values on each virtual bond, length L+1. May be None if the MPS is not in canonical form. Otherwise, _S[i] is to the left of _B[i]. We recommend using get_SL(), get_SR(), set_SL(), set_SR(), which takes proper care of the boundary conditions. Sometimes (e.g. during DMRG with an enabled mixer), entries may temporarily be a non-diagonal tenpy.linalg.np_conserved.Array to be inserted between the left canonical ‘A’ tensors on the left and right-canonical _B[i] on the right.

Type:

list of {None, 1D array, Array}

_valid_forms

Class attribute. Mapping for canonical forms to a tuple (nuL, nuR) indicating that self._B[i] = s[i]**nuL -- Gamma[i] -- s[i]**nuR is saved.

Type:

dict

_valid_bc

Class attribute. Possible valid boundary conditions.

Type:

tuple of str

_transfermatrix_keep

How many states to keep at least when diagonalizing a TransferMatrix. Important if the state develops a near-degeneracy.

Type:

int

_p_label, _B_labels

Class attribute. _p_label defines the physical legs of the B-tensors, _B_labels lists all the labels of the B tensors. Used by methods like get_theta() to avoid the necessity of re-implementations for derived classes like the Purification_MPS if just the number of physical legs changed.

Type:

list of str

test_sanity()[source]

Sanity check, raises ValueErrors, if something is wrong.

copy()[source]

Returns a copy of self.

The copy still shares the sites, chinfo, and LegCharges of the B tensors, but the values of B and S are deeply copied.

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 sites, chinfo (under these names), _B as "tensors", _S as "singular_values", bc as "boundary_condition", and form converted to a single array of shape (L, 2) as "canonical_form", Moreover, it saves norm, L, grouped and _transfermatrix_keep (as “transfermatrix_keep”) as HDF5 attributes, as well as the maximum of chi under the name “max_bond_dimension”.

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.

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 represent 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

classmethod from_lat_product_state(lat, p_state, allow_incommensurate=False, **kwargs)[source]

Construct an MPS from a product state given in lattice coordinates.

This is a wrapper around from_product_state(). The purpose is to make the p_state argument independent of the order of the Lattice, and specify it in terms of lattice indices instead.

Parameters:
  • lat (Lattice) – The underlying lattice defining the geometry and Hilbert Space.

  • p_state (array_like of {int | str | 1D array}) – Defines the product state to be represented. Should be of dimension lat.dim`+1, entries are indexed by lattice indices. Entries of the array as for the `p_state argument of from_product_state(). It gets tiled to the shape lat.shape, if it is smaller.

  • allow_incommensurate (bool) – Allow an incommensurate tiling of p_state to the full lattice. For example, if you pass p_state=[['up'], ['down']] with a tenpy.models.lattice.Chain, this function raises an error for an odd number of sites in the Chain, but if you set allow_incommensurate=True, it will still work and give you a state with total Sz = +1/2 for odd sites (since total Sz=0 doesn’t fit).

  • **kwargs – Other keyword arguments as defined in from_product_state(). bc is set by default from lat.bc_MPS.

Returns:

product_mps – An MPS representing the specified product state.

Return type:

MPS

Examples

Let’s first consider a Ladder composed of a SpinHalfSite and a FermionSite.

To initialize a state of up-spins on the spin sites and half-filled fermions, you can use:

Note that the same p_state works for a finite lattice of even length, say L=10, as well. We then just “tile” in x-direction, i.e., repeat the specified state 5 times:

You can also easily half-fill a Honeycomb, for example with only the A sites occupied, or as stripe parallel to the x-direction (stripe_x, alternating along y axis), or as stripes parallel to the y-direction (stripe_y, alternating along x axis).

classmethod from_product_state(sites, p_state, bc='finite', dtype=<class 'numpy.float64'>, permute=True, form='B', chargeL=None)[source]

Construct a matrix product state from a given product state.

Parameters:
  • sites (list of Site) – The sites defining the local Hilbert space.

  • p_state (list of {int | str | 1D array}) – Defines the product state to be represented; one entry for each site of the MPS. An entry of str type is translated to an int with the help of state_labels(). An entry of int type represents the physical index of the state to be used. An entry which is a 1D array defines the complete wavefunction on that site; this allows to make a (local) superposition.

  • bc ({'infinite', 'finite', 'segment'}) – MPS boundary conditions. See docstring of MPS.

  • dtype (type or string) – The data type of the array entries.

  • permute (bool) – The Site might permute the local basis states if charge conservation gets enabled. If permute is True (default), we permute the given p_state locally according to each site’s perm. The p_state entries should then always be given as if conserve=None in the Site.

  • form ((list of) {'B' | 'A' | 'C' | 'G' | None | tuple(float, float)}) – Defines the canonical form. See module doc-string. A single choice holds for all of the entries.

  • chargeL (charges) – Leg charges at bond 0, which are purely conventional.

Returns:

product_mps – An MPS representing the specified product state.

Return type:

MPS

Examples

Example to get a Neel state for a TFIChain:

>>> from tenpy.networks.mps import MPS
>>> L = 10
>>> M = tenpy.models.tf_ising.TFIChain({'L': L})
>>> p_state = ["up", "down"] * (L//2)  # repeats entries L/2 times
>>> psi = MPS.from_product_state(M.lat.mps_sites(), p_state, bc=M.lat.bc_MPS)

The meaning of the labels "up","down" is defined by the Site, in this example a SpinHalfSite.

Extending the example, we can replace the spin in the center with one with arbitrary angles theta, phi in the bloch sphere. However, note that you can not write this bloch state (for theta != 0, pi) when conserving symmetries, as the two physical basis states correspond to different symmetry sectors.

>>> spin = tenpy.networks.site.SpinHalfSite(conserve=None, sort_charge=False)
>>> p_state = ["up", "down"] * (L//2)  # repeats entries L/2 times
>>> theta, phi = np.pi/4, np.pi/6
>>> bloch_sphere_state = np.array([np.cos(theta/2), np.exp(1.j*phi)*np.sin(theta/2)])
>>> p_state[L//2] = bloch_sphere_state   # replace one spin in center
>>> psi = MPS.from_product_state([spin]*L, p_state, bc=M.lat.bc_MPS, dtype=complex)

Note that for the more general SpinChain, the order of the two entries for the bloch_sphere_state would be exactly the opposite (when we keep the the north-pole of the bloch sphere being the up-state). The reason is that the SpinChain uses the general SpinSite, where the states are ordered ascending from 'down' to 'up'. The SpinHalfSite on the other hand uses the order 'up', 'down' where that the Pauli matrices look as usual.

classmethod from_Bflat(sites, Bflat, SVs=None, bc='finite', dtype=None, permute=True, form='B', legL=None)[source]

Construct a matrix product state from a set of numpy arrays Bflat and singular vals.

Parameters:
  • sites (list of Site) – The sites defining the local Hilbert space.

  • Bflat (iterable of numpy ndarrays) – The matrix defining the MPS on each site, with legs 'p', 'vL', 'vR' (physical, virtual left/right).

  • SVs (list of 1D array | None) – The singular values on each bond. Should always have length L+1. By default (None), set all singular values to the same value. Entries out of nontrivial_bonds are ignored.

  • bc ({'infinite', 'finite', 'segment'}) – MPS boundary conditions. See docstring of MPS.

  • dtype (type or string) – The data type of the array entries. Defaults to the common dtype of Bflat.

  • permute (bool) – The Site might permute the local basis states if charge conservation gets enabled. If permute is True (default), we permute the given Bflat locally according to each site’s perm. The Bflat should then always be given as if conserve=None in the Site.

  • form ((list of) {'B' | 'A' | 'C' | 'G' | None | tuple(float, float)}) – Defines the canonical form of Bflat. See module doc-string. A single choice holds for all of the entries.

  • leg_L (LegCharge | None) – Leg charges at bond 0, which are purely conventional. If None, use trivial charges.

Returns:

mps – An MPS with the matrices Bflat converted to npc arrays.

Return type:

MPS

classmethod from_full(sites, psi, form=None, cutoff=1e-16, normalize=True, bc='finite', outer_S=None)[source]

Construct an MPS from a single tensor psi with one leg per physical site.

Performs a sequence of SVDs of psi to split off the B matrices and obtain the singular values, the result will be in canonical form. Obviously, this is only well-defined for finite or segment boundary conditions.

Parameters:
  • sites (list of Site) – The sites defining the local Hilbert space.

  • psi (Array) – The full wave function to be represented as an MPS. Should have labels 'p0', 'p1', ...,  'p{L-1}'. Additionally, it may have (or must have for ‘segment’ bc) the legs 'vL', 'vR', which are trivial for ‘finite’ bc.

  • form ('B' | 'A' | 'C' | 'G' | None) – The canonical form of the resulting MPS, see module doc-string. None defaults to ‘A’ form on the first site and ‘B’ form on all following sites.

  • cutoff (float) – Cutoff of singular values used in the SVDs.

  • normalize (bool) – Whether the resulting MPS should have ‘norm’ 1.

  • bc ('finite' | 'segment') – Boundary conditions.

  • outer_S (None | (array, array)) – For ‘segment’ bc the singular values on the left and right of the considered segment, None for ‘finite’ boundary conditions.

Returns:

psi_mps – MPS representation of psi, in canonical form and possibly normalized.

Return type:

MPS

classmethod from_singlets(site, L, pairs, up='up', down='down', lonely=[], lonely_state='up', bc='finite')[source]

Create an MPS of entangled singlets.

Parameters:
  • site (Site) – The site defining the local Hilbert space, taken uniformly for all sites.

  • L (int) – The number of sites.

  • pairs (list of (int, int)) – Pairs of sites to be entangled; the returned MPS will have a singlet for each pair in pairs. For bc='infinite' MPS, some indices can be outside the range [0, L) to indicate couplings across infinite MPS unit cells.

  • up (int | str) – A singlet is defined as (|up down> - |down up>)/2**0.5, up and down give state indices or labels defined on the corresponding site.

  • down (int | str) – A singlet is defined as (|up down> - |down up>)/2**0.5, up and down give state indices or labels defined on the corresponding site.

  • lonely (list of int) – Sites which are not included into a singlet pair.

  • lonely_state (int | str) – The state for the lonely sites.

  • bc ({'infinite', 'finite', 'segment'}) – MPS boundary conditions. See docstring of MPS.

Returns:

singlet_mps – An MPS representing singlets on the specified pairs of sites.

Return type:

MPS

classmethod from_product_mps_covering(mps_covering, index_map, bc='finite')[source]

Create an MPS as a product of (many) local mps covering all sites to be created.

This is a generalization of from_singlets() to allow arbitrary local, entangled states over multiple sites. Those local states are represented by MPS in mps_covering, such that each site in the final MPS gets it state from exactly one local MPS.

For example to reproduce from_singlets(), you define L/2 local two-site mps in a singlet state, and use the pairs as index_map. (If you have lonely states not entangled into a singlet, add a trivial one-site MPS for them as well.) Indeed, if you look into the source code of from_singlets(), this is exactly what it does (at least since this method was implemented).

More generally, you can easily initialize any kind of valence bond solid state, if you just initialize the corresponding local states and generate a corresponding index_map of the geometry of pairs, see the example below.

Parameters:
  • mps_covering (list of MPS) – List of local, ‘finite’ MPS.

  • index_map (list of tuple of int) – For each local_psi in mps_covering, add one tuple with local_psi.L entries which sites in the returned psi should

Returns:

psi – An MPS constructed as explained above.

Return type:

MPS

Example

Say you have three local MPS with mps_covering = [psi_A, psi_B, psi_C] with A, B and C tensors on 3, 2, and 2 sites, respectively. Using the index_map=[[0, 1, 3], [2, 5], [4, 6]] would combine them as:

|    A0-A1----A2 C0----C1
|    |  |     |  |     |
|    |  |  B0-|--|--B1 |
|    |  |  |  |  |  |  |
|    0  1  2  3  4  5  6

Using the index_map=[[1, 2, 0], [3, 5], [6, 4]] would rather combine them as:

|    A2----.     C1----C0
|    |     |     |     |
|    |  A0-A1 B0-|--B1 |
|    |  |  |  |  |  |  |
|    0  1  2  3  4  5  6

As another example, let’s generalize from_singlets() to “spin-full” fermions represented by the (spin-less!) FermionSite with an extra index for the spin, here u=0,1 on a 2D square lattice.

>>> ferm = FermionSite(conserve='N')
>>> lat = MultiSpeciesLattice(Square(4, 2, None), [ferm]*2, ['up', 'down'])
>>> ferm_up_down = MPS.from_product_state([ferm]*4, ['full', 'empty', 'empty', 'full'])
>>> ferm_down_up = MPS.from_product_state([ferm]*4, ['empty', 'full', 'full', 'empty'])
>>> ferm_singlet = ferm_up_down.add(ferm_down_up, 0.5**0.5, -0.5**0.5)
>>> index_map = [[(x, y, 0), (x, y, 1), (x+1, y, 0), (x+1, y, 1)]
...    for (x, y) in [(0, 0), (0, 1), (2, 0), (2, 1)]]
>>> index_map = [[lat.lat2mps_idx(x_y_u) for x_y_u in pairs] for pairs in index_map]
>>> psi = MPS.from_product_mps_covering([ferm_singlet]*4, index_map)

This will generate a singlet valence bold solid (VBS) state looking like this:

|    x--x  x--x
|
|    x--x  x--x
property L

Number of physical sites; for an iMPS the len of the MPS unit cell.

property dim

List of local physical dimensions.

property finite

Distinguish MPS vs iMPS.

True for an MPS (bc='finite', 'segment'), False for an iMPS (bc='infinite').

property chi

Dimensions of the (nontrivial) virtual bonds.

property nontrivial_bonds

Slice of the non-trivial bond indices, depending on self.bc.

get_B(i, form='B', copy=False, cutoff=1e-16, label_p=None)[source]

Return (view of) B at site i in canonical form.

Parameters:
  • i (int) – Index choosing the site.

  • form ('B' | 'A' | 'C' | 'G' | 'Th' | None | tuple(float, float)) – The (canonical) form of the returned B. For None, return the matrix in whatever form it is. If any of the tuple entry is None, also don’t scale on the corresponding axis.

  • copy (bool) – Whether to return a copy even if form matches the current form.

  • cutoff (float) – During DMRG with a mixer, S may be a matrix for which we need the inverse. This is calculated as the Penrose pseudo-inverse, which uses a cutoff for the singular values.

  • label_p (None | str) – Ignored by default (None). Otherwise replace the physical label 'p' with 'p'+label_p'. (For derived classes with more than one “physical” leg, replace all the physical leg labels accordingly.)

Returns:

B – The MPS ‘matrix’ B at site i with leg labels 'vL', 'p', 'vR'. May be a view of the matrix (if copy=False), or a copy (if the form changed or copy=True).

Return type:

Array

:raises ValueError : if self is not in canonical form and form is not None.:

set_B(i, B, form='B')[source]

Set B at site i.

Parameters:
  • i (int) – Index choosing the site.

  • B (Array) – The ‘matrix’ at site i. No copy is made! Should have leg labels 'vL', 'p', 'vR' (not necessarily in that order).

  • form ('B' | 'A' | 'C' | 'G' | 'Th' | None | tuple(float, float)) – The (canonical) form of the B to set. None stands for non-canonical form.

set_svd_theta(i, theta, trunc_par=None, update_norm=False)[source]

SVD a two-site wave function theta and save it in self.

Parameters:
  • i (int) – theta is the wave function on sites i, i + 1.

  • theta (Array) – The two-site wave function with labels combined into "(vL.p0)", "(p1.vR)", ready for svd.

  • trunc_par (None | dict) – Parameters for truncation, see truncation. If None, no truncation is done.

  • update_norm (bool) – If True, multiply the norm of theta into norm.

get_SL(i)[source]

Return singular values on the left of site i

get_SR(i)[source]

Return singular values on the right of site i

set_SL(i, S)[source]

Set singular values on the left of site i

set_SR(i, S)[source]

Set singular values on the right of site i

get_theta(i, n=2, cutoff=1e-16, formL=1.0, formR=1.0)[source]

Calculates the n-site wavefunction on sites[i:i+n].

Parameters:
  • i (int) – Site index.

  • n (int) – Number of sites. The result lives on sites[i:i+n].

  • cutoff (float) – During DMRG with a mixer, S may be a matrix for which we need the inverse. This is calculated as the Penrose pseudo-inverse, which uses a cutoff for the singular values.

  • formL (float) – Exponent for the singular values to the left.

  • formR (float) – Exponent for the singular values to the right.

Returns:

theta – The n-site wave function with leg labels vL, p0, p1, .... p{n-1}, vR. In Vidal’s notation (with s=lambda, G=Gamma): theta = s**form_L G_i s G_{i+1} s ... G_{i+n-1} s**form_R.

Return type:

Array

convert_form(new_form='B')[source]

Transform self into different canonical form (by scaling the legs with singular values).

Parameters:

new_form ((list of) {'B' | 'A' | 'C' | 'G' | 'Th' | None | tuple(float, float)}) – The form the stored ‘matrices’. The table in module doc-string. A single choice holds for all of the entries.

:raises ValueError : if trying to convert from a None form. Use canonical_form() instead!:

increase_L(new_L=None)[source]

Modify self inplace to enlarge the MPS unit cell; in place.

Deprecated since version 0.5.1: This method will be removed in version 1.0.0. Use the equivalent psi.enlarge_mps_unit_cell(new_L//psi.L) instead of psi.increase_L(new_L).

Parameters:

new_L (int) – New number of sites. Needs to be an integer multiple of L. Defaults to 2*self.L.

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 unit cell will be increased from L to factor*L.

roll_mps_unit_cell(shift=1)[source]

Shift the section we define as unit cell of an infinite MPS; in place.

Suppose we have a unit cell with tensors [A, B, C, D] (repeated on both sites). With shift = 1, the new unit cell will be [D, A, B, C], whereas shift = -1 will give [B, C, D, A].

Parameters:

shift (int) – By how many sites to move the tensors to the right.

overlap_translate_finite(psi, shift=1)[source]

Contract <self|T^N|psi> for translation T with finite, periodic boundaries.

Looks like this for shift=1, with the open virtual legs contracted in the end:

--B[L-1] Th[0] -- B[1] -- B[2] -- ..... B[L-2] --
   |      |       |       |             |
 Th*[0] --B*[1] --B*[2] --B*[3] --..... B*[L-1]

An alternative to calling this method would be to call permute_sites() followed by overlap(). Note that permute_sites uses truncation on the way, though, which could severely affect the precision especially for general (non-local) permutations while this function contracts everything exactly (possibly at higher numerical cost).

Parameters:
  • psi (MPS) – MPS to take overlap with.

  • shift (int) – Translation by how many sites. Note that for large shift, the contraction is \(O(\chi^4)\) compared to DMRG etc scaling as \(O(\chi^3)\).

Returns:

self_Tn_psi – Contraction of <self|T^N|psi>.

Return type:

float

See also

permute_sites

Allows more general permutations of the sites.

overlap

Directly the overlap between two MPS without translation.

roll_mps_unit_cell

Effectively applies T^shift on infinite MPS.

enlarge_chi(extra_legs, random_fct=<built-in method normal of numpy.random.mtrand.RandomState object>)[source]

Artificially enlarge the bond dimension by the specified extra legs/charges. In place.

First converts the MPS in B form. This function then fills up the ‘vR’ leg with zeros, and then groups physical and vR legs to fill up the ‘vL’ leg with orthogonal rows. In than way, we get an MPS with larger bond dimension, still in right-canonical form (assuming no “overcomplete” charge block), representing the same state, with the additional singular values being exactly zero.

Note

You should probably choose the extra charges to be sensible, to expand into the space you are interested in, and not just into a random direction!

Parameters:
  • extra_legs (list of {None, LegCharge, int}) – The extra charges to be added on the virtual legs, with qconj=+1. Length L +1 for finite, length L for infinite, with entry i left of site i. If an int is given, fill up with a single block of charges like the Schmidt state with highest weight. Note that this might force the resulting state to not be in strict “right-canonical” B form if that charge block becomes overcomplete.

  • random_fct – Function generating random entries to choose extra orthogonal rows in the B tensors. Should accept a size keyword for the shape, and return numpy arrays.

Returns:

permutations – Permutation performed on each virtual leg, such that new_S = concatenate(old_S, zeros)[perm].

Return type:

list of 1D array | None

spatial_inversion()[source]

Perform a spatial inversion along the MPS.

Exchanges the first with the last tensor and so on, i.e., exchange site i with site L-1 - i. This is equivalent to a mirror/reflection with the bond left of L/2 (even L) or the site (L-1)/2 (odd L) as a fixpoint. For infinite MPS, the bond between MPS unit cells is another fix point.

group_sites(n=2, grouped_sites=None)[source]

Modify self inplace to group sites.

Group each n sites together using the GroupedSite. This might allow to do TEBD with a Trotter decomposition, or help the convergence of DMRG (in case of too long range interactions).

Parameters:
  • n (int) – Number of sites to be grouped together.

  • grouped_sites (None | list of GroupedSite) – The sites grouped together.

See also

group_split

Reverts the grouping.

group_split(trunc_par=None)[source]

Modify self inplace to split previously grouped sites.

Parameters:

trunc_par (dict) – Parameters for truncation, see truncation. Defaults to {'chi_max': max(self.chi)}.

Returns:

trunc_err – The error introduced by the truncation for the splitting.

Return type:

TruncationError

See also

group_sites

Should have been used before to combine sites.

get_grouped_mps(blocklen)[source]

Like group_sites(), but make a copy.

Parameters:

blocklen (int) – Number of subsequent sites to be combined; n in group_sites().

Returns:

New MPS object with bunched sites.

Return type:

grouped_MPS

extract_segment(first, last)[source]

Extract an segment from a finite or infinite MPS.

Parameters:
  • first (int) – The first and last site to include into the segment.

  • last (int) – The first and last site to include into the segment.

Returns:

psi_segment – Copy of self with ‘segment’ boundary conditions.

Return type:

MPS

get_total_charge(only_physical_legs=False)[source]

Calculate and return the qtotal of the whole MPS (when contracted).

If set, the segment_boundaries are included (unless only_physical_legs is True).

Parameters:

only_physical_legs (bool) – For 'finite' boundary conditions, the total charge can be gauged away by changing the LegCharge of the trivial legs on the left and right of the MPS. This option allows to project out the trivial legs to get the actual “physical” total charge.

Returns:

qtotal – The sum of the qtotal of the individual B tensors.

Return type:

charges

gauge_total_charge(qtotal=None, vL_leg=None, vR_leg=None)[source]

Gauge the legcharges of the virtual bonds such that the MPS has a total qtotal.

Acts in place, i.e. changes the B tensors. Make a (shallow) copy if needed.

Parameters:
  • qtotal ((list of) charges) – If a single set of charges is given, it is the desired total charge of the MPS (which get_total_charge() will return afterwards). By default (None), use 0 charges, unless vL_leg and vR_leg are specified, in which case we adjust the total charge to match these legs.

  • vL_leg (None | LegCharge) – Desired new virtual leg on the very left and right. Needs to have the same block structure as the current legs, but can have shifted charge entries. For infinite MPS, we need vL_leg to be the conjugate leg of vR_leg. For segment MPS, these legs are the outer-most legs, possibly including the segment_boundaries.

  • vR_leg (None | LegCharge) – Desired new virtual leg on the very left and right. Needs to have the same block structure as the current legs, but can have shifted charge entries. For infinite MPS, we need vL_leg to be the conjugate leg of vR_leg. For segment MPS, these legs are the outer-most legs, possibly including the segment_boundaries.

entanglement_entropy(n=1, bonds=None, for_matrix_S=False)[source]

Calculate the (half-chain) entanglement entropy for all nontrivial bonds.

Consider a bipartition of the system into \(A = \{ j: j <= i_b \}\) and \(B = \{ j: j > i_b\}\) and the reduced density matrix \(\rho_A = tr_B(\rho)\). The von-Neumann entanglement entropy is defined as \(S(A, n=1) = -tr(\rho_A \log(\rho_A)) = S(B, n=1)\). The generalization for n != 1, n>0 are the Renyi entropies: \(S(A, n) = \frac{1}{1-n} \log(tr(\rho_A^2)) = S(B, n=1)\)

This function calculates the entropy for a cut at different bonds i, for which the the eigenvalues of the reduced density matrix \(\rho_A\) and \(\rho_B\) is given by the squared schmidt values S of the bond.

Parameters:
  • n (int/float) – Selects which entropy to calculate; n=1 (default) is the usual von-Neumann entanglement entropy.

  • bonds (None | (iterable of) int) – Selects the bonds at which the entropy should be calculated. None defaults to range(0, L+1)[self.nontrivial_bonds], i.e., range(1, L) for ‘finite’ MPS and range(0, L) for ‘infinite’ MPS.

  • for_matrix_S (bool) – Switch calculate the entanglement entropy even if the _S are matrices. Since \(O(\chi^3)\) is expensive compared to the usual \(O(\chi)\), we raise an error by default.

Returns:

entropies – Entanglement entropies for half-cuts. entropies[j] contains the entropy for a cut at bond bonds[j], i.e. between sites bonds[j]-1 and bonds[j]. For infinite systems with default bonds=None, this means that entropies[0] will be a cut left of site 0 and is the one you should look at to e.g. study the scaling of the entanglement with chi or to extract the topological entanglement entropy - don’t take the average over bonds, in particular if you have 2D cylinders or ladders. On the contrary, for finite systems with bonds=None, take the central value of the returned array entropies[len(entropies)//2)] == entropies[(L-1)//2] (and not just entropies[L//2]) to extract the half-chain entanglement entropy.

Return type:

1D ndarray

entanglement_entropy_segment(segment=[0], first_site=None, n=1)[source]

Calculate entanglement entropy for general geometry of the bipartition.

This function is similar as entanglement_entropy(), but for more general geometry of the region A to be a segment of a few sites.

This is achieved by explicitly calculating the reduced density matrix of A and thus works only for small segments. The alternative entanglement_entropy_segment2() might work for larger segments at small enough bond dimensions.

Parameters:
  • segment (list of int) – Given a first site i, the region A_i is defined to be [i+j for j in segment].

  • first_site (None | (iterable of) int) – Calculate the entropy for segments starting at these sites. None defaults to range(L-segment[-1]) for finite or range(L) for infinite boundary conditions.

  • n (int | float) – Selects which entropy to calculate; n=1 (default) is the usual von-Neumann entanglement entropy, otherwise the n-th Renyi entropy.

Returns:

entropiesentropies[i] contains the entropy for the the region A_i defined above.

Return type:

1D ndarray

entanglement_entropy_segment2(segment, n=1)[source]

Calculate entanglement entropy for general geometry of the bipartition.

This function is similar to entanglement_entropy_segment(), but allows more sites in segment. The trick is to exploit that for a pure state (which the MPS represents) and a bipartition into regions A and B, the entropy is the same in both regions, \(S(A) = S(B)\). Hence we can trace out the specified segment and obtain \(\rho_B = tr_A(rho)\), where A is the specified segment. The price is a huge computation cost of \(O(chi^6 d^{3x})\) where x is the number of physical legs not included into segment between min(segment) and max(segment).

Parameters:
  • segment (list of int) – The site indices specifying region A. We calculate and diagonalize the full reduced density matrix of the complement of A.

  • n (int | float) – Selects which entropy to calculate; n=1 (default) is the usual von-Neumann entanglement entropy, otherwise the n-th Renyi entropy.

Returns:

entropy – The entropy for the the region defined by the segment (or equivalently it’s complement).

Return type:

float

entanglement_spectrum(by_charge=False)[source]

return entanglement energy spectrum.

Parameters:

by_charge (bool) – Whether we should sort the spectrum on each bond by the possible charges.

Returns:

ent_spectrum – For each (non-trivial) bond the entanglement spectrum. If by_charge is False, return (for each bond) a sorted 1D ndarray with the convention \(S_i^2 = e^{-\xi_i}\), where \(S_i\) labels a Schmidt value and \(\xi_i\) labels the entanglement ‘energy’ in the returned spectrum. If by_charge is True, return a a list of tuples (charge, sub_spectrum) for each possible charge on that bond.

Return type:

list

get_rho_segment(segment)[source]

Return reduced density matrix for a segment.

Note that the dimension of rho_A scales exponentially in the length of the segment.

Parameters:

segment (iterable of int) – Sites for which the reduced density matrix is to be calculated. Assumed to be sorted.

Returns:

rho – Reduced density matrix of the segment sites. Labels 'p0', 'p1', ..., 'pk', 'p0*', 'p1*', ..., 'pk*' with k=len(segment).

Return type:

Array

probability_per_charge(bond=0)[source]

Return probabilities of charge value on the left of a given bond.

For example for particle number conservation, define \(N_b = sum_{i<b} n_i\) for a given bond b. This function returns the possible values of N_b as rows of charge_values, and for each row the probability that this combination occurs in the given state.

Parameters:

bond (int) – The bond to be considered. The returned charges are summed on the left of this bond.

Returns:

  • charge_values (2D array) – Columns correspond to the different charges in self.chinfo. Rows are the different charge fluctuations at this bond

  • probabilities (1D array) – For each row of charge_values the probability for these values of charge fluctuations.

average_charge(bond=0)[source]

Return the average charge for the block on the left of a given bond.

For example for particle number conservation, define \(N_b = sum_{i<b} n_i\) for a given bond b. Then this function returns \(<\psi| N_b |\psi>\).

Parameters:

bond (int) – The bond to be considered. The returned charges are summed over the sites left of bond.

Returns:

average_charge – For each type of charge in chinfo the average value when summing the charge values over sites left of the given bond.

Return type:

1D array

charge_variance(bond=0)[source]

Return the charge variance on the left of a given bond.

For example for particle number conservation, define \(N_b = sum_{i<b} n_i\) for a given bond b. Then this function returns \(<\psi| N_b^2 |\psi> - (<\psi| N_b |\psi>)^2\).

Parameters:

bond (int) – The bond to be considered. The returned charges are summed over the sites left of bond.

Returns:

average_charge – For each type of charge in chinfo the variance of of the charge values left of the given bond.

Return type:

1D array

mutinf_two_site(max_range=None, n=1)[source]

Calculate the two-site mutual information \(I(i:j)\).

Calculates \(I(i:j) = S(i) + S(j) - S(i,j)\), where \(S(i)\) is the single site entropy on site \(i\) and \(S(i,j)\) the two-site entropy on sites \(i,j\).

Parameters:
  • max_range (int) – Maximal distance |i-j| for which the mutual information should be calculated. None defaults to L-1.

  • n (float) – Selects the entropy to use, see entropy().

Returns:

  • coords (2D array) – Coordinates for the mutinf array.

  • mutinf (1D array) – mutinf[k] is the mutual information \(I(i:j)\) between the sites i, j = coords[k].

overlap(other, charge_sector=None, ignore_form=False, understood_infinite=False, **kwargs)[source]

Compute overlap <self|other>.

Parameters:
  • other (MPS) – An MPS with the same physical sites.

  • charge_sector (None | charges | 0) – Selects the charge sector in which the dominant eigenvector of the TransferMatrix is. None stands for all sectors, 0 stands for the sector of zero charges. If a sector is given, it assumes the dominant eigenvector is in that charge sector.

  • ignore_form (bool) – If False (default), take into account the canonical form form at each site. If True, we ignore the canonical form (i.e., whether the MPS is in left, right, mixed or no canonical form) and just contract all the _B as they are. (This can give different results!)

  • understood_infinite (bool) – Raise a warning to make aware of A warning about infinite MPS. Set understood_infinite=True to suppress the warning.

  • **kwargs – Further keyword arguments given to TransferMatrix.eigenvectors(); only used for infinite boundary conditions.

Returns:

overlap – The contraction <self|other> * self.norm * other.norm (i.e., taking into account the norm of both MPS). For an infinite MPS, <self|other> is the overlap per unit cell, i.e., the largest eigenvalue of the TransferMatrix.

Return type:

dtype.type

expectation_value_terms_sum(term_list, prefactors=None)[source]

Calculate expectation values for a bunch of terms and sum them up.

This is equivalent to the following expression:

sum([self.expectation_value_term(term)*strength for term, strength in term_list])

However, for efficiency, the term_list is converted to an MPO and the expectation value of the MPO is evaluated.

Deprecated since version 0.4.0: prefactor will be removed in version 1.0.0. Instead, directly give just TermList(term_list, prefactors) as argument.

Parameters:
  • term_list (TermList) – The terms and prefactors (strength) to be summed up.

  • prefactors – Instead of specifying a TermList, one can also specify the term_list and strength separately. This is deprecated.

Returns:

  • terms_sum (list of (complex) float) – Equivalent to the expression sum([self.expectation_value_term(term)*strength for term, strength in term_list]).

  • _mpo – Intermediate results: the generated MPO. For a finite MPS, terms_sum = _mpo.expectation_value(self), for an infinite MPS terms_sum = _mpo.expectation_value(self) * self.L

See also

expectation_value_term

evaluates a single term.

tenpy.networks.mpo.MPO.expectation_value

expectation value density of an MPO.

sample_measurements(first_site=0, last_site=None, ops=None, rng=None, norm_tol=1e-12)[source]

Sample measurement results in the computational basis.

This function samples projective measurements on a contiguous range of sites, tracing out the remaining sites.

Note that for infinite boundary conditions, the probability of sampling a set of sigmas is not |psi.overlap(MPS.from_product_state(sigmas, ...))|^2, because the latter would project to the set sigmas on each (translated) MPS unit cell, while this function is only projecting to them in a single MPS unit cell.

Parameters:
  • first_site (int) – Take measurements on the sites in range(first_site, last_site + 1). last_site defaults to L - 1.

  • last_site (int) – Take measurements on the sites in range(first_site, last_site + 1). last_site defaults to L - 1.

  • ops (list of str) – If not None, sample in the eigenbasis of self.sites[i].get_op(ops[(i - first_site) % len(ops)]) and directly return the corresponding eigenvalue in sigmas.

  • rng (numpy.random.Generator) – The random number generator; if None, a new numpy.random.default_rng() is generated.

  • norm_tol (float) – Tolerance

Returns:

  • sigmas (list of int | list of float) – On each site the index of the local basis that was measured, as specified in the corresponding Site in sites. Note that this can change depending on whether/what charges you conserve! Explicitly specifying the measurement operator will avoid that issue.

  • weight (float) – The weight sqrt(trace(|psi><psi|sigmas...><sigmas...|)), i.e., the probability of measuring |sigmas...> is weight**2. For a finite system where we sample all sites (i.e., the trace over the compliment of the sites is trivial), this is the actual overlap <sigmas...|psi> including the phase.

norm_test()[source]

Check that self is in canonical form.

Returns:

norm_error – For each site the norm error to the left and right. The error norm_error[i, 0] is defined as the norm-difference between the following networks:

|   --theta[i]---.       --s[i]--.
|       |        |    vs         |
|   --theta*[i]--.       --s[i]--.

Similarly, norm_error[i, 1] is the norm-difference of:

|   .--theta[i]---         .--s[i+1]--
|   |    |          vs     |
|   .--theta*[i]--         .--s[i+1]--

Return type:

array, shape (L, 2)

canonical_form(**kwargs)[source]

Bring self into canonical ‘B’ form, (re-)calculate singular values.

Simply calls canonical_form_finite() or canonical_form_infinite1(). Keyword arguments are passed on to the corresponding specialized versions.

canonical_form_finite(renormalize=True, cutoff=0.0, envs_to_update=None)[source]

Bring a finite (or segment) MPS into canonical form (in place).

If any site is in form None, it does not use any of the singular values S (for ‘finite’ boundary conditions, or only the very left S for ‘segment’ b.c.). If all sites have a form, it respects the form to ensure that one S is included per bond. The final state is always in right-canonical ‘B’ form.

Performs one sweep left to right doing QR decompositions, and one sweep right to left doing SVDs calculating the singular values.

Parameters:
  • renormalize (bool) – Whether a change in the norm should be discarded or used to update norm. Note that even renormalize=True does not reset the norm to 1. To do that, you would rather have to set psi.norm = 1 explicitly!

  • cutoff (float | None) – Cutoff of singular values used in the SVDs.

  • envs_to_update (None | list of MPSEnvironment) – Clear the environments; for segment also update the left/rightmost LP/RP.

Returns:

U_L, V_R – Only returned for 'segment' boundary conditions. The unitaries defining the new left and right Schmidt states in terms of the old ones, with legs 'vL', 'vR'.

Return type:

Array

canonical_form_infinite(**kwargs)[source]

Deprecated wrapper around canonical_form_infinite1().

canonical_form_infinite1(renormalize=True, tol_xi=1000000.0)[source]

Bring an infinite MPS into canonical form (in place).

If any site is in form None, it does not use any of the singular values S. If all sites have a form, it respects the form to ensure that one S is included per bond. The final state is always in right-canonical ‘B’ form.

Proceeds in three steps, namely 1) diagonalize right and left transfermatrix on a given bond to bring that bond into canonical form, and then 2) sweep right to left, and 3) left to right to bringing other bonds into canonical form.

Parameters:
  • renormalize (bool) – Whether a change in the norm should be discarded or used to update norm. Note that even renormalize=True does not reset the norm to 1. To do that, you would rather have to set psi.norm = 1 explicitly!

  • tol_xi (float) – Raise an error if the correlation length is larger than that (which indicates a degenerate “cat” state, e.g., for spontaneous symmetry breaking).

canonical_form_infinite2(renormalize=True, tol=1e-15, arnoldi_params=None, cutoff=1e-15)[source]

Convert infinite MPS to canonical form.

Implementation following Algorithm 1,2 in [vanderstraeten2019].

Parameters:
  • renormalize (bool) – Whether a change in the norm should be discarded or used to update norm. Note that even renormalize=True does not reset the norm to 1. To do that, you would rather have to set psi.norm = 1 explicitly!

  • tol (float) – Precision down to which the state should be in canonical form.

  • arnoldi_params (dict) – Parameters for Arnoldi.

  • cutoff – Truncation cutoff for small singular values.

correlation_length(target=1, tol_ev0=1e-08, charge_sector=0, return_charges=False)[source]

Calculate the correlation length by diagonalizing the transfer matrix.

Assumes that self is in canonical form.

Works only for infinite MPS, where the transfer matrix is a useful concept. Assuming a single-site unit cell, any correlation function splits into \(C(A_i, B_j) = A'_i T^{j-i-1} B'_j\) with some parts left and right and the \(j-i-1\)-th power of the transfer matrix in between. The largest eigenvalue is 1 (if self is properly normalized) and gives the dominant contribution of \(A'_i E_1 * 1^{j-i-1} * E_1^T B'_j = <A> <B>\), and the second largest one gives a contribution \(\propto \lambda_2^{j-i-1}\). Thus \(\lambda_2 = \exp(-\frac{1}{\xi})\).

More general for a L-site unit cell we get \(\lambda_2 = \exp(-\frac{L}{\xi})\), where the xi is given in units of 1 lattice spacing in the MPS.

Warning

For a higher-dimensional lattice (which the MPS class doesn’t know about), the correct unit is the lattice spacing in x-direction, and the correct formula is \(\lambda_2 = \exp(-\frac{L_x}{\xi})\), where L_x is the number of lattice spacings in the infinite direction within the MPS unit cell, e.g. the number of “rings” of a cylinder in the MPS unit cell. To get to these units, divide the returned xi by the number of sites within a “ring”, for a lattice given in N_sites_per_ring.

Parameters:
  • target (int) – We look for the target + 1 largest eigenvalues.

  • tol_ev0 (float | None) – Print warning if largest eigenvalue deviates from 1 by more than tol_ev0. If None, assume dominant eigenvector in 0-charge-sector to be 1 for non-zero charge_sector.

  • charge_sector (None | list of int | 0 | list of list of int) – Selects the charge sector (=list of int) in which we look for the dominant eigenvalue of the TransferMatrix. None stands for all sectors, 0 stands for the zero-charge sector. Defaults to 0, i.e., assumes the dominant eigenvector is in charge sector 0. If you pass a list of charge sectors (i.e. 2D array of int), this function returns target dominant eigenvalues in each of those sectors.

  • return_charges (bool) – If True, return the charge sectors along with the eigenvalues.

Returns:

  • xi (float | 1D array) – If target = 1, return just the correlation length, otherwise an array of the target largest correlation lengths. It is measured in units of a single spacing between sites in the MPS language, see the warning above.

  • charge_sectors (list of charge sectors) – For each entry in xi the charge sector, i.e., qtotal of the dominant eigenvalue.

See also

correlation_length_charge_sectors

lists possible charge sectors.

correlation_length_charge_sectors(drop_symmetric=True, include_0=True)[source]

Return possible charge_sector argument for correlation_length().

The correlation_length() is calculated from eigenvalues of the transfer matrix. The left/right eigenvectors correspond to contraction of the left/right parts of the transfer-matrix in the network of the correlation_function(). The charge_sector one can pass to the correlation_length() (or TransferMatrix, respectively) is the qtotal that eigenvector.

Since bra and ket are identical for the correlation_length(), one can flip top and bottom and to an overall conjugate, and gets back to the same TransferMatrix, hence eigenvalues of a given eigenvector and it’s dagger (seen as a matrix with legs 'vL', 'vL*') are identical. Since that flips the sign of all charges, we can conclude that the correlation length in a given charge sector and the negative charge sector are identical. The option drop_symmetric hence allows to only return charge sectors where the negative charge sector was not yet returned.

Parameters:

drop_symmetric (bool) – See above.

add(other, alpha, beta, cutoff=1e-15)[source]

Return an MPS which represents alpha|self> + beta |others>.

Works only for ‘finite’, ‘segment’ boundary conditions. For ‘segment’ boundary conditions, the virtual legs on the very left/right are assumed to correspond to each other (i.e. self and other have the same state outside of the considered segment). Takes into account norm.

Parameters:
  • other (MPS) – Another MPS of the same length to be added with self.

  • alpha (complex float) – Prefactors for self and other. We calculate alpha * |self> + beta * |other>

  • beta (complex float) – Prefactors for self and other. We calculate alpha * |self> + beta * |other>

  • cutoff (float | None) – Cutoff of singular values used in the SVDs.

Returns:

sum – An MPS representing alpha|self> + beta |other>. Has same total charge as self.

Return type:

MPS

apply_local_op(i, op, unitary=None, renormalize=False, cutoff=1e-13, understood_infinite=False)[source]

Apply a local (one or multi-site) operator to self.

Note that this destroys the canonical form if the local operator is non-unitary. Therefore, this function calls canonical_form() if necessary.

For infinite MPS, it applies the operator in parallel within each unit site, i.e., really it applies \(\prod_{u \in \mathbb{Z}} (\mathrm{op}_{i + u *L})\) where L is the number of sites in the MPS unit cell.

Parameters:
  • i (int) – (Left-most) index of the site(s) on which the operator should act.

  • op (str | npc.Array) – A physical operator acting on site i, with legs 'p', 'p*' for a single-site operator or with legs ['p0', 'p1', ...], ['p0*', 'p1*', ...] for an operator acting on n>=2 sites. Strings (like 'Id', 'Sz') are translated into single-site operators defined by sites.

  • unitary (None | bool) – Whether op is unitary, i.e., whether the canonical form is preserved (True) or whether we should call canonical_form() (False). None checks whether norm(op dagger(op) - identity) is smaller than cutoff.

  • renormalize (bool) – Whether the final state should keep track of the norm (False, default), or the change of norm should be discarded (True).

  • cutoff (float) – Cutoff for singular values if op acts on more than one site (see from_full()). (And used as cutoff for a unspecified unitary.)

  • understood_infinite (bool) – Raise a warning to make aware of A warning about infinite MPS. Set understood_infinite=True to suppress the warning.

apply_product_op(ops, unitary=None, renormalize=False)[source]

Apply a (global) product of local onsite operators to self.

Note that this destroys the canonical form if any local operator is non-unitary. Therefore, this function calls canonical_form() if necessary.

The result is equivalent to the following loop, but more efficient by avoiding intermediate calls to canonical_form() inside the loop:

for i, op in enumerate(ops):
    self.apply_local_op(i, op, unitary, renormalize, cutoff)
Parameters:
  • ops ((list of) str | npc.Array) – List of onsite operators to apply on each site, with legs 'p', 'p*'. Strings (like 'Id', 'Sz') are translated into single-site operators defined by sites.

  • unitary (None | bool) – Whether op is unitary, i.e., whether the canonical form is preserved (True) or whether we should call canonical_form() (False). None checks whether max(norm(op dagger(op) - identity) for op in ops) < 1.e-14

  • renormalize (bool) – Whether the final state should keep track of the norm (False, default), or the change of norm should be discarded (True).

perturb(randomize_params=None, close_1=True, canonicalize=None)[source]

Locally perturb the state a little bit; in place.

Parameters:
swap_sites(i, swap_op='auto', trunc_par=None)[source]

Swap the two neighboring sites i and i+1 (inplace).

Exchange two neighboring sites: form theta, ‘swap’ the physical legs and split with an svd. While the ‘swap’ is just a transposition/relabeling for bosons, one needs to be careful about the signs from Jordan-Wigner strings for fermions.

Parameters:
  • i (int) – Swap the two sites at positions i and i+1.

  • swap_op (None | 'auto', 'autoInv' | Array) – The operator used to swap the physical legs of the two-site wave function theta. For None, just transpose/relabel the legs. Alternative give an npc Array which represents the full operator used for the swap. Should have legs ['p0', 'p1', 'p0*', 'p1*'] with 'p0', 'p1*' contractible. For 'auto' we try to be smart about fermionic signs, see note below.

  • trunc_par (dict) – Parameters for truncation, see truncation.

Returns:

trunc_err – The error of the represented state introduced by the truncation after the swap.

Return type:

TruncationError

Notes

For fermions, it’s crucial to use the correct swap_op. The swap_op is a two-site operator exchanging ‘p0’ and ‘p1’ legs. For bosons, this is really just a relabeling (done for swap_op=None). Alternatively, you can construct the operator explicitly like this:

siteL, siteR = psi.sites[i], psi.sites[i+1]
dL, dR = siteL.dim, siteR.dim
legL, legR = siteL.leg, siteR.leg
swap_op_dense = np.eye(dL*dR)
swap_op = npc.Array.from_ndarray(swap_op_dense.reshape([dL, dR, dL, dR]),
                                 [legL, legR, legL.conj(), legR.conj()],
                                 labels=['p1', 'p0', 'p0*', 'p1*'])

However, for fermions we need to be very careful about the Jordan-Wigner strings. Let’s derive how the operator should look like.

You can write a state as

\[|\psi> = \sum_{[n_j]} \psi_{[n_j]} \prod_j (c^\dagger_j)^{n_j} |vac>\]

where [n_j] denotes a set of \(n_j \in [0, 1]\) for each physical site j and the product over j is taken in increasing order. Let \(P\) be the operator switching i <-> i+1, with inverse \(P^\dagger\). Then:

\[\begin{split}P |\psi> = \sum_{[n_j]} \psi_{[n_i]} P \prod_j (c^\dagger_i)^{n_j} |vac> \\ = \sum_{[n_j]} \psi_{[n_i]} P \prod_j (c^\dagger_i)^{n_j} |vac>\end{split}\]

When \(P\) acts on the product of \(c^\dagger_{i}\) operators, it commutes \((c^\dagger_i)^{n_i}\) with \((c^\dagger_{i+1})^{n_{i+1}}\). This gives a a sign \((-1)^{n_i * n_{i+1}}\). We must hence include this sign in the swap operator. The n_i in the equations above is given by JW_exponent. This leads to the following swap operator used for fermions with swap_op='auto', suitable to just permute sites:

siteL, siteR = psi.sites[i], psi.sites[i+1]
dL, dR = siteL.dim, siteR.dim
legL, legR = siteL.leg, siteR.leg
n_i_n_j = np.outer(siteL.JW_exponent, siteR.JW_exponent).reshape(dL*dR)
swap_op_dense = np.diag((-1)**n_i_n_j)
swap_op = npc.Array.from_ndarray(swap_op_dense.reshape([dL, dR, dL, dR]),
                                 [legL, legR, legL.conj(), legR.conj()],
                                 labels=['p1', 'p0', 'p0*', 'p1*'])

In some cases you might want to use a more complicated swap operator. As outlined in (the appendix of) [shapourian2017], a typical hamiltonian of the form \(H = -t \sum_i c_i^\dagger c_{i+1} + h.c. + \text{density interaction}\) is invariant under a reflection \(R\) acting as \(R c^e_x R^\dagger = i c^o_{-x}\) and \(R c^o_x R^\dagger = i c^e_{-x}\) for even/odd fermion sites. The following code includes the factor of \(i\), or rather \(-i\) since we have creation operators, into the swap operator and is used with swap_op='autoInv':

siteL, siteR = psi.sites[i], psi.sites[i+1]
dL, dR = siteL.dim, siteR.dim
legL, legR = siteL.leg, siteR.leg
n_i = np.outer(siteL.JW_exponent, np.ones(dR)).reshape(dL*dR)
n_j = np.outer(np.ones(dL), siteR.JW_exponent).reshape(dL*dR)
swap_op_dense = np.diag((-1)**(n_i * n_j) * (-1.j)**n_i * (-1.j)**n_j)
swap_op = npc.Array.from_ndarray(swap_op_dense.reshape([dL, dR, dL, dR]),
                                 [legL, legR, legL.conj(), legR.conj()],
                                 labels=['p1', 'p0', 'p0*', 'p1*'])
permute_sites(perm, swap_op='auto', trunc_par=None, verbose=None)[source]

Applies the permutation perm to the state (inplace).

Deprecated since version 0.8.0: Drop / ignore verbose argument, never print something.

Parameters:
  • perm (ndarray[ndim=1, int]) – The applied permutation, such that psi.permute_sites(perm)[i] = psi[perm[i]] (where [i] indicates the i-th site).

  • swap_op (None | 'auto', 'autoInv' | Array) – The operator used to swap the physical legs of a two-site wave function theta, see swap_sites().

  • trunc_par (dict) – Parameters for truncation, see truncation.

Returns:

trunc_err – The error of the represented state introduced by the truncation after the swaps.

Return type:

TruncationError

compute_K(perm, swap_op='auto', trunc_par=None, canonicalize=1e-06, verbose=None, expected_mean_k=0.0)[source]

Compute the momentum quantum numbers of the entanglement spectrum for 2D states.

Works for an infinite MPS living on a cylinder, infinitely long in x direction and with periodic boundary conditions in y directions. If the state is invariant under ‘rotations’ around the cylinder axis, one can find the momentum quantum numbers of it. (The rotation is nothing more than a translation in y.) This function permutes some sites (on a copy of self) to enact the rotation, and then finds the dominant eigenvector of the mixed transfer matrix to get the quantum numbers, along the lines of [pollmann2012], see also (the appendix and Fig. 11 in the arXiv version of) [cincio2013].

Deprecated since version 0.8.0: Drop / ignore verbose argument, never print something.

Parameters:
  • perm (1D ndarray | Lattice) – Permutation to be applied to the physical indices, see permute_sites(). If a lattice is given, we use it to read out the lattice structure and shift each site by one lattice-vector in y-direction (assuming periodic boundary conditions). (If you have a CouplingModel, give its lat attribute for this argument)

  • swap_op (None | 'auto', 'autoInv' | Array) – The operator used to swap the physical legs of a two-site wave function theta, see swap_sites().

  • trunc_par (dict) – Parameters for truncation, see truncation.

  • canonicalize (float) – Check that self is in canonical form; call canonical_form() if norm_test() yields np.linalg.norm(self.norm_test()) > canonicalize.

  • expected_mean_k (float) – As explained in [cincio2013], the returned W is extracted as eigenvector of a mixed transfer matrix, and hence has an undefined phase. We fix the overall phase such that sum(s[j]**2 exp(iK[j]) == np.sum(W) = np.exp(1.j*expected_mean_k).

Returns:

  • U (Array) – Unitary representation of the applied permutation on left Schmidt states.

  • W (ndarray) – 1D array of the form S**2 exp(i K), where S are the Schmidt values on the left bond. You can use np.abs() and np.angle() to extract the (squared) Schmidt values S and momenta K from W.

  • q (LegCharge) – LegCharge corresponding to W.

  • ov (complex) – The eigenvalue of the mixed transfer matrix <psi|T|psi> per L sites. An absolute value different smaller than 1 indicates that the state is not invariant under the permutation or that the truncation error trunc_err was too large!

  • trunc_err (TruncationError) – The error of the represented state introduced by the truncation after swaps when performing the truncation.

compress(options)[source]

Compress an MPS.

Options

config MPS_compress
option summary

chi_list (from Sweep) in IterativeSweeps.reset_stats

By default (``None``) this feature is disabled. [...]

chi_list_reactivates_mixer (from Sweep) in IterativeSweeps.sweep

If True, the mixer is reset/reactivated each time the bond dimension growth [...]

combine (from Sweep) in Sweep

Whether to combine legs into pipes. This combines the virtual and [...]

compression_method

Mandatory. [...]

init_env_data (from Sweep) in DMRGEngine.init_env

Dictionary as returned by ``self.env.get_initialization_data()`` from [...]

lanczos_params (from Sweep) in Sweep

Lanczos parameters as described in :cfg:config:`Lanczos`.

max_sweeps (from VariationalCompression) in VariationalCompression

Maximum number of sweeps to perform for the compression.

min_sweeps (from VariationalCompression) in VariationalCompression

Minimum number of sweeps to perform for the compression.

mixer (from Sweep) in DMRGEngine.mixer_activate

Specifies which :class:`Mixer` to use, if any. [...]

mixer_params (from Sweep) in DMRGEngine.mixer_activate

Mixer parameters as described in :cfg:config:`Mixer`.

orthogonal_to (from Sweep) in DMRGEngine.init_env

Deprecated in favor of the `orthogonal_to` function argument (forwarded fro [...]

start_env (from Sweep) in DMRGEngine.init_env

Number of sweeps to be performed without optimization to update the environment.

start_env_sites (from VariationalCompression) in VariationalCompression

Number of sites to contract for the initial LP/RP environment in case of in [...]

tol_theta_diff (from VariationalCompression) in VariationalCompression

Stop after less than `max_sweeps` sweeps if the 1-site wave function change [...]

trunc_params

Truncation parameters as described in :cfg:config:`truncation`.

option compression_method: 'SVD' | 'variational'

Mandatory. Selects the method to be used for compression. For the SVD compression, trunc_params is the only other option used.

option trunc_params: dict

Truncation parameters as described in truncation.

compress_svd(trunc_par)[source]

Compress self with a single sweep of SVDs; in place.

Perform a single right-sweep of QR/SVD without truncation, followed by a left-sweep with truncation, very much like canonical_form_finite().

Warning

In case of a strong compression, this does not find the optimal, global solution.

Parameters:

trunc_par (dict) – Parameters for truncation, see truncation.

correlation_function(ops1, ops2, sites1=None, sites2=None, opstr=None, str_on_first=True, hermitian=False, autoJW=True)[source]

Correlation function <bra|op1_i op2_j|ket> of single site operators, sandwiched between bra and ket. For examples the contraction for a two-site operator on site i would look like:

|          .--S--B[i]--B[i+1]--...--B[j]---.
|          |     |     |            |      |
|          |     |     |            op2    |
|          LP[i] |     |            |      RP[j]
|          |     op1   |            |      |
|          |     |     |            |      |
|          .--S--B*[i]-B*[i+1]-...--B*[j]--.

Onsite terms are taken in the order <psi | op1 op2 | psi>. If opstr is given and str_on_first=True, it calculates:

|           for i < j                               for i > j
|
|          .--S--B[i]---B[i+1]--...- B[j]---.     .--S--B[j]---B[j+1]--...- B[i]---.
|          |     |      |            |      |     |     |      |            |      |
|          |     opstr  opstr        op2    |     |     op2    |            |      |
|          LP[i] |      |            |      RP[j] LP[j] |      |            |      RP[i]
|          |     op1    |            |      |     |     opstr  opstr        op1    |
|          |     |      |            |      |     |     |      |            |      |
|          .--S--B*[i]--B*[i+1]-...- B*[j]--.     .--S--B*[j]--B*[j+1]-...- B*[i]--.

For i==j, no opstr is included. For str_on_first=False, the opstr on site min(i, j) is always left out. Strings (like 'Id', 'Sz') in the arguments are translated into single-site operators defined by the Site on which they act. Each operator should have the two legs 'p', 'p*'.

Warning

This function is only evaluating correlation functions by moving right, and hence can be inefficient if you try to vary the left end while fixing the right end. In that case, you might be better off (=faster evaluation) by using term_correlation_function_left() with a small for loop over the right indices.

Parameters:
  • ops1 ((list of) { Array | str }) – First operator of the correlation function (acting after ops2). If a list is given, ops1[i] acts on site i of the MPS. Note that even if a list is given, we still just evaluate two-site correlations! psi.correlation_function(['A','B'], ['C', 'D']) evaluates <A_i C_j> for even i and even j, <B_i C_j> for even i and odd j, <B_i C_j> for odd i and even j, and <B_i D_j> for odd i and odd j.

  • ops2 ((list of) { Array | str }) – Second operator of the correlation function (acting before ops1). If a list is given, ops2[j] acts on site j of the MPS.

  • sites1 (None | int | list of int) – List of site indices i; a single int is translated to range(0, sites1). None defaults to all sites range(0, L). Is sorted before use, i.e. the order is ignored.

  • sites2 (None | int | list of int) – List of site indices; a single int is translated to range(0, sites2). None defaults to all sites range(0, L). Is sorted before use, i.e. the order is ignored.

  • opstr (None | (list of) { Array | str }) – Ignored by default (None). Operator(s) to be inserted between ops1 and ops2. If less than L operators are given, we repeat them periodically. If given as a list, opstr[r] is inserted at site r (independent of sites1 and sites2).

  • str_on_first (bool) – Whether the opstr is included on the site min(i, j). Note the order, which is chosen that way to handle fermionic Jordan-Wigner strings correctly. (In other words: choose str_on_first=True for fermions!)

  • hermitian (bool) – Optimization flag: if sites1 == sites2 and Ops1[i]^\dagger == Ops2[i] (which is not checked explicitly!), the resulting C[x, y] will be hermitian. We can use that to avoid calculations, so hermitian=True will run faster.

  • autoJW (bool) – Ignored if opstr is given. If True, auto-determine if a Jordan-Wigner string is needed. Works only if exclusively strings were used for op1 and op2.

Returns:

C – The correlation function C[x, y] = <bra|ops1[i] ops2[j]|ket>, where ops1[i] acts on site i=sites1[x] and ops2[j] on site j=sites2[y]. If opstr is given, it gives (for str_on_first=True): - For i < j: C[x, y] = <bra|ops1[i] prod_{i <= r < j} opstr[r] ops2[j]|ket>. - For i > j: C[x, y] = <bra|prod_{j <= r < i} opstr[r] ops1[i] ops2[j]|ket>. - For i = j: C[x, y] = <bra|ops1[i] ops2[j]|ket>. The condition <= r is replaced by a strict < r, if str_on_first=False.

Warning

The MPSEnvironment variant of this method takes the accumulated MPS norm into account, which is non-trivial e.g. when you used apply_local_op with non-unitary operators.

In contrast, the MPS variant of this method ignores the norm, i.e. returns the expectation value for the normalized state.

Return type:

2D ndarray

Examples

Let’s prepare a state in alternating |+z>, |+x> states:

>>> spin_half = tenpy.networks.site.SpinHalfSite(conserve=None)
>>> p_state = ['up', [np.sqrt(0.5), -np.sqrt(0.5)]]*3
>>> psi = tenpy.networks.mps.MPS.from_product_state([spin_half]*6, p_state, "infinite")

Default arguments calculate correlations for all i and j within the MPS unit cell. To evaluate the correlation function for a single i, you can use sites1=[i]. Alternatively, you can use term_correlation_function_right() (or term_correlation_function_left()):

>>> psi.correlation_function("Sz", "Sx")  
array([[ 0.  , -0.25,  0.  , -0.25,  0.  , -0.25],
       [ 0.  ,  0.  ,  0.  ,  0.  ,  0.  ,  0.  ],
       [ 0.  , -0.25,  0.  , -0.25,  0.  , -0.25],
       [ 0.  ,  0.  ,  0.  ,  0.  ,  0.  ,  0.  ],
       [ 0.  , -0.25,  0.  , -0.25,  0.  , -0.25],
       [ 0.  ,  0.  ,  0.  ,  0.  ,  0.  ,  0.  ]])
>>> psi.correlation_function("Sz", "Sx", [0])
array([[ 0.  , -0.25,  0.  , -0.25,  0.  , -0.25]])
>>> corr1 = psi.correlation_function("Sz", "Sx", [0], range(1, 10))
>>> corr2 = psi.term_correlation_function_right([("Sz", 0)], [("Sx", 0)], 0, range(1, 10))
>>> assert np.all(np.abs(corr2 - corr1) < 1.e-12)

For fermions, it auto-determines that/whether a Jordan Wigner string is needed:

>>> fermion = tenpy.networks.site.FermionSite(conserve='N')
>>> p_state = ['empty', 'full'] * 3
>>> psi = tenpy.networks.mps.MPS.from_product_state([fermion]*6, p_state, "finite")
>>> CdC = psi.correlation_function("Cd", "C")  # optionally: use `hermitian=True`
>>> psi.correlation_function("C", "Cd")[1, 2] == -CdC[2, 1]
True
>>> np.all(np.diag(CdC) == psi.expectation_value("Cd C"))  # "Cd C" is equivalent to "N"
True

See also

expectation_value_term

for a single combination of i and j of A_i B_j`.

term_correlation_function_right

for correlations between multi-site terms, fix left term.

term_correlation_function_left

for correlations between multi-site terms, fix right term.

expectation_value(ops, sites=None, axes=None)[source]

Expectation value <bra|ops|ket> of (n-site) operator(s).

Calculates n-site expectation values of operators sandwiched between bra and ket. For examples the contraction for a two-site operator on site i would look like:

|          .--S--B[i]--B[i+1]--.
|          |     |     |       |
|          |     |-----|       |
|          LP[i] | op  |       RP[i+1]
|          |     |-----|       |
|          |     |     |       |
|          .--S--B*[i]-B*[i+1]-.

Here, the B are taken from ket, the B* from bra. For MPS expectation values these are the same and LP/ RP are trivial.

Parameters:
  • ops ((list of) { Array | str }) – The operators, for which the expectation value should be taken, All operators should all have the same number of legs (namely 2 n). If less than len(sites) operators are given, we repeat them periodically. Strings (like 'Id', 'Sz') are translated into single-site operators defined by sites.

  • sites (list) – List of site indices. Expectation values are evaluated there. If None (default), the entire chain is taken (clipping for finite b.c.)

  • axes (None | (list of str, list of str)) – Two lists of each n leg labels giving the physical legs of the operator used for contraction. The first n legs are contracted with conjugated B, the second n legs with the non-conjugated B. None defaults to (['p'], ['p*']) for single site (n=1), or (['p0', 'p1', ... 'p{n-1}'], ['p0*', 'p1*', .... 'p{n-1}*']) for n > 1.

Returns:

exp_vals – Expectation values, exp_vals[i] = <bra|ops[i]|ket>, where ops[i] acts on site(s) j, j+1, ..., j+{n-1} with j=sites[i].

Warning

The MPSEnvironment variant of this method takes the accumulated MPS norm into account, which is non-trivial e.g. when you used apply_local_op with non-unitary operators.

In contrast, the MPS variant of this method ignores the norm, i.e. returns the expectation value for the normalized state.

Return type:

1D ndarray

Examples

Let’s prepare a state in alternating |+z>, |+x> states:

>>> spin_half = tenpy.networks.site.SpinHalfSite(conserve=None)
>>> p_state = ['up', [np.sqrt(0.5), -np.sqrt(0.5)]]*3
>>> psi = tenpy.networks.mps.MPS.from_product_state([spin_half]*6, p_state)

One site examples (n=1):

>>> Sz = psi.expectation_value('Sz')
>>> print(Sz)
[0.5 0.  0.5 0.  0.5 0. ]
>>> Sx = psi.expectation_value('Sx')
>>> print(Sx)
[ 0.  -0.5  0.  -0.5  0.  -0.5]
>>> print(psi.expectation_value(['Sz', 'Sx']))
[ 0.5 -0.5  0.5 -0.5  0.5 -0.5]
>>> print(psi.expectation_value('Sz', sites=[0, 3, 4]))
[0.5 0.  0.5]

Two site example (n=2), assuming homogeneous sites:

>>> SzSx = npc.outer(psi.sites[0].Sz.replace_labels(['p', 'p*'], ['p0', 'p0*']),
...                  psi.sites[1].Sx.replace_labels(['p', 'p*'], ['p1', 'p1*']))
>>> print(psi.expectation_value(SzSx))  # note: len L-1 for finite bc, or L for infinite
[-0.25  0.   -0.25  0.   -0.25]

Example measuring <psi|SzSx|psi> on each second site, for inhomogeneous sites:

>>> SzSx_list = [npc.outer(psi.sites[i].Sz.replace_labels(['p', 'p*'], ['p0', 'p0*']),
...                        psi.sites[i+1].Sx.replace_labels(['p', 'p*'], ['p1', 'p1*']))
...              for i in range(0, psi.L-1, 2)]
>>> print(psi.expectation_value(SzSx_list, range(0, psi.L-1, 2)))
[-0.25 -0.25 -0.25]

Expectation value with different bra and ket in an MPSEnvironment:

>>> spin_half = tenpy.networks.site.SpinHalfSite(conserve=None)
>>> p2_state = [[np.sqrt(0.5), -np.sqrt(0.5)], 'up']*3
>>> phi = tenpy.networks.mps.MPS.from_product_state([spin_half]*6, p2_state)
>>> env = tenpy.networks.mps.MPSEnvironment(phi, psi)
>>> Sz = env.expectation_value('Sz')
>>> print(Sz)
[0.0625 0.0625 0.0625 0.0625 0.0625 0.0625]
expectation_value_multi_sites(operators, i0)[source]

Expectation value <bra|op0_{i0}op1_{i0+1}...opN_{i0+N}|ket>.

Calculates the expectation value of a tensor product of single-site operators acting on different sites next to each other. In other words, evaluate the expectation value of a term op0_i0 op1_{i0+1} op2_{i0+2} ..., looking like this (with op short for operators, for len(operators)=3):

|          .--S--B[i0]---B[i0+1]--B[i0+2]--B[i0+3]--.
|          |     |       |        |        |        |
|         LP[i0] op[0]   op[1]    op[2]    op[3]    RP[i0+3]
|          |     |       |        |        |        |
|          .--S--B*[i0]--B*[i0+1]-B*[i0+2]-B*[i0+3]-.

Warning

This function does not automatically add Jordan-Wigner strings! For correct handling of fermions, use expectation_value_term() instead.

Parameters:
  • operators (List of { Array | str }) – List of one-site operators. This method calculates the expectation value of the n-sites operator given by their tensor product.

  • i0 (int) – The left most index on which an operator acts, i.e., operators[i] acts on site i + i0.

Returns:

exp_val – The expectation value of the tensorproduct of the given onsite operators, <bra|operators[0]_{i0} operators[1]_{i0+1} ... |ket>.

Warning

The MPSEnvironment variant of this method takes the accumulated MPS norm into account, which is non-trivial e.g. when you used apply_local_op with non-unitary operators.

In contrast, the MPS variant of this method ignores the norm, i.e. returns the expectation value for the normalized state.

Return type:

float/complex

expectation_value_term(term, autoJW=True)[source]

Expectation value <bra|op_{i0}op_{i1}...op_{iN}|ket>.

Calculates the expectation value of a tensor product of single-site operators acting on different sites i0, i1, … (not necessarily next to each other). In other words, evaluate the expectation value of a term op0_i0 op1_i1 op2_i2 ....

For example the contraction of three one-site operators on sites i0, i1=i0+1, i2=i0+3 would look like:

|          .--S--B[i0]---B[i0+1]--B[i0+2]--B[i0+3]--.
|          |     |       |        |        |        |
|         LP[i0]op1     op2       |       op3       RP[i0+3]
|          |     |       |        |        |        |
|          .--S--B*[i0]--B*[i0+1]-B*[i0+2]-B*[i0+3]-.
Parameters:
  • term (list of (str, int)) – List of tuples op, i where i is the MPS index of the site the operator named op acts on. The order inside term determines the order in which they act (in the mathematical convention: the last operator in term is right-most, so it acts first on a ket).

  • autoJW (bool) – If True (default), automatically insert Jordan Wigner strings for Fermions as needed.

Returns:

exp_val – The expectation value of the tensorproduct of the given onsite operators, <bra|op_i0 op_i1 ... op_iN |ket>.

Warning

The MPSEnvironment variant of this method takes the accumulated MPS norm into account, which is non-trivial e.g. when you used apply_local_op with non-unitary operators.

In contrast, the MPS variant of this method ignores the norm, i.e. returns the expectation value for the normalized state.

Return type:

float/complex

See also

correlation_function

efficient way to evaluate many correlation functions.

Examples

>>> a = psi.expectation_value_term([('Sx', 2), ('Sz', 4)])
>>> b = psi.expectation_value_term([('Sz', 4), ('Sx', 2)])
>>> c = psi.expectation_value_multi_sites(['Sz', 'Id', 'Sx'], i0=2)
>>> assert a == b == c
get_op(op_list, i)[source]

Given a list of operators, select the one corresponding to site i.

Parameters:
  • op_list ((list of) {str | npc.array}) – List of operators from which we choose. We assume that op_list[j] acts on site j. If the length is shorter than L, we repeat it periodically. Strings are translated using get_op() of site i.

  • i (int) – Index of the site on which the operator acts.

Returns:

op – One of the entries in op_list, not copied.

Return type:

npc.array

term_correlation_function_left(term_L, term_R, i_L=None, j_R=0, autoJW=True, opstr=None)[source]

Correlation function between (multi-site) terms, moving the left term, fix right term.

Same as term_correlation_function_right(), but vary index i of the left term instead of the j of the right term.

term_correlation_function_right(term_L, term_R, i_L=0, j_R=None, autoJW=True, opstr=None)[source]

Correlation function between (multi-site) terms, moving the right term, fix left term.

For term_L = [('A', 0), ('B', 1)] and term_R = [('C', 0), ('D', 1)], calculate the correlation function \(A_{i+0} B_{i+1} C_{j+0} D_{j+1}\) for fixed i and varying j according to i_L/j_R. The terms may not overlap. For fermions, the order of the terms is following the usual mathematical convention, where term_R acts first on a physical ket.

Warning

This function assumes that bra and ket are normalized, i.e. for MPSEnvironment. Thus you may want to take into account MPS.norm of both bra and ket.

Parameters:
  • term_L (list of (str, int)) – Each a term representing a sum of operators on different sites, e.g., [('Sz', 0), ('Sz', 1)] or [('Cd', 0), ('C', 1)].

  • term_R (list of (str, int)) – Each a term representing a sum of operators on different sites, e.g., [('Sz', 0), ('Sz', 1)] or [('Cd', 0), ('C', 1)].

  • i_L (int) – Offset added to the indices of term_L.

  • j_R (list of int | None) – List of offsets to be added to the indices of term_R. Is sorted before use, i.e. the order is ignored. For finite MPS, None defaults to range(j0, L), where j0 is chosen such that term_R starts one site right of the term_L. For infinite MPS, None defaults to range(L, 11*L, L), i.e., one term per MPS unit cell for a distance of up to 10 unit cells.

  • autoJW (bool) – Whether to automatically take care of Jordan-Wigner strings.

  • opstr (str) – Force an intermediate operator string to used inbetween the terms. Can only be used in combination with autoJW=False.

Returns:

corrs – Values of the correlation function, one for each entry in the list j_R.

Warning

The MPSEnvironment variant of this method takes the accumulated MPS norm into account, which is non-trivial e.g. when you used apply_local_op with non-unitary operators.

In contrast, the MPS variant of this method ignores the norm, i.e. returns the expectation value for the normalized state.

Return type:

1D array

See also

correlation_function

varying both i and j at once.

term_list_correlation_function_right

generalization to sums of terms on the left/right.

term_list_correlation_function_right(term_list_L, term_list_R, i_L=0, j_R=None, autoJW=True, opstr=None)[source]

Correlation function between sums of multi-site terms, moving the right sum of term.

Generalization of term_correlation_function_right() to the case where term_list_L and term_R are sums of terms. This function calculates <bra|term_list_L[i_L] term_list_R[j]|ket> for j in j_R.

Assumes that overall terms with an odd number of operators requiring a Jordan-Wigner string don’t contribute. (In systems conserving the fermionic particle number (parity), this is true.)

Parameters:
  • term_list_L (TermList) – Each a TermList representing the sum of terms to be applied.

  • term_list_R (TermList) – Each a TermList representing the sum of terms to be applied.

  • i_L (int) – Offset added to all the indices of term_list_L.

  • j_R (list of int | None) – List of offsets to be added to the indices of term_list_R. Is sorted before use, i.e. the order is ignored. For finite MPS, None defaults to range(j0, L), where j0 is chosen such that term_R starts one site right of the term_L. For infinite MPS, None defaults to range(L, 11*L, L), i.e., one term per MPS unit cell for a distance of up to 10 unit cells.

  • autoJW (bool) – Whether to automatically take care of Jordan-Wigner strings.

  • opstr (str) – Force an intermediate operator string to be used inbetween the terms. (Even used within the term_list_L/R for terms with smaller-than maximal support.) Can only be used in combination with autoJW=False.

Returns:

corrs – Values of the correlation function, one for each entry in the list j_R.

Warning

The MPSEnvironment variant of this method takes the accumulated MPS norm into account, which is non-trivial e.g. when you used apply_local_op with non-unitary operators.

In contrast, the MPS variant of this method ignores the norm, i.e. returns the expectation value for the normalized state.

Return type:

1D array

See also

term_correlation_function_right

version for a single term in both term_list_L/R.