MPS

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

Bases: object

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

Parameters
siteslist of Site

Defines the local Hilbert space for each site.

Bslist of Array

The ‘matrices’ of the MPS. Labels are vL, vR, p (in any order).

SVslist 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 tabel 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.

Attributes
L

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

chi

Dimensions of the (nontrivial) virtual bonds.

finite

Distinguish MPS vs iMPS.

nontrivial_bonds

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

siteslist of Site

Defines the local Hilbert space for each site.

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

Boundary conditions as described in above table.

formlist of {None | tuple(float, float)}

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).

chinfoChargeInfo

The nature of the charge.

dtypetype

The data type of the _B.

normfloat

The norm of the state, i.e. sqrt(<psi|psi>). Ignored for (normalized) expectation_value(), but important for overlap().

groupedint

Number of sites grouped together, see group_sites().

_Blist of npc.Array

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.

_Slist of (None | 1D array)

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.

_valid_formsdict

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

_valid_bctuple of str

Valid boundary conditions.

_transfermatrix_keepint

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

Methods

add(self, other, alpha, beta[, cutoff])

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

apply_local_op(self, i, op[, unitary, …])

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

average_charge(self[, bond])

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

canonical_form(self[, renormalize])

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

canonical_form_finite(self[, renormalize, …])

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

canonical_form_infinite(self[, renormalize, …])

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

charge_variance(self[, bond])

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

compute_K(self, perm[, swap_op, trunc_par, …])

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

convert_form(self[, new_form])

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

copy(self)

Returns a copy of self.

correlation_function(self, ops1, ops2[, …])

Correlation function <psi|op1_i op2_j|psi>/<psi|psi> of single site operators.

correlation_length(self[, target, tol_ev0, …])

Calculate the correlation length by diagonalizing the transfer matrix.

entanglement_entropy(self[, n, bonds, …])

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

entanglement_entropy_segment(self[, …])

Calculate entanglement entropy for general geometry of the bipartition.

entanglement_spectrum(self[, by_charge])

return entanglement energy spectrum.

expectation_value(self, ops[, sites, axes])

Expectation value <psi|ops|psi>/<psi|psi> of (n-site) operator(s).

expectation_value_multi_sites(self, …)

Expectation value <psi|op0_{i0}op1_{i0+1}...opN_{i0+N}|psi>/<psi|psi>.

expectation_value_term(self, term[, autoJW])

Expectation value <psi|op_{i0}op_{i1}...op_{iN}|psi>/<psi|psi>.

expectation_value_terms_sum(self, term_list)

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

from_Bflat(sites, Bflat[, SVs, bc, dtype, …])

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

from_full(sites, psi[, form, cutoff, …])

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

from_product_state(sites, p_state[, bc, …])

Construct a matrix product state from a given product state.

from_singlets(site, L, pairs[, up, down, …])

Create an MPS of entangled singlets.

gauge_total_charge(self[, qtotal, vL_leg, …])

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

get_B(self, i[, form, copy, cutoff, label_p])

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

get_SL(self, i)

Return singular values on the left of site i

get_SR(self, i)

Return singular values on the right of site i

get_grouped_mps(self, blocklen)

contract blocklen subsequent tensors into a single one and return result as a new MPS.

get_op(self, op_list, i)

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

get_rho_segment(self, segment)

Return reduced density matrix for a segment.

get_theta(self, i[, n, cutoff, formL, formR])

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

get_total_charge(self[, only_physical_legs])

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

group_sites(self[, n, grouped_sites])

Modify self inplace to group sites.

group_split(self[, trunc_par])

Modify self inplace to split previously grouped sites.

increase_L(self[, new_L])

Modify self inplace to enlarge the unit cell.

mutinf_two_site(self[, max_range, n])

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

norm_test(self)

Check that self is in canonical form.

overlap(self, other[, charge_sector, …])

Compute overlap <self|other>.

permute_sites(self, perm[, swap_op, …])

Applies the permutation perm to the state (inplace).

probability_per_charge(self[, bond])

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

set_B(self, i, B[, form])

Set B at site i.

set_SL(self, i, S)

Set singular values on the left of site i

set_SR(self, i, S)

Set singular values on the right of site i

swap_sites(self, i[, swap_op, trunc_par])

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

test_sanity(self)

Sanity check, raises ValueErrors, if something is wrong.

test_sanity(self)[source]

Sanity check, raises ValueErrors, if something is wrong.

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
siteslist of Site

The sites defining the local Hilbert space.

p_stateiterable of {int | str | 1D array}

Defines the product state to be represented. If p_state[i] is str, then site i is in state self.sites[i].state_labels(p_state[i]). If p_state[i] is int, then site i is in state p_state[i]. If p_state[i] is an array, then site i wavefunction is p_state[i].

bc{‘infinite’, ‘finite’, ‘segmemt’}

MPS boundary conditions. See docstring of MPS.

dtypetype or string

The data type of the array entries.

permutebool

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 argument 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.

chargeLcharges

Leg charge at bond 0, which are purely conventional.

Returns
product_mpsMPS

An MPS representing the specified product state.

Examples

Example to get a Neel state for a TIChain:

>>> M = TFIChain({'L': 10})
>>> 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:

>>> M = TFIChain({'L': 8, 'conserve': None})
>>> p_state = ["up", "down"] * (L//2)  # repeats entries L/2 times
>>> 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(M.lat.mps_sites(), p_state, bc=M.lat.bc_MPS, dtype=np.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 orderd ascending from 'down' to 'up'. The SpinHalfSite on the other hand uses the order 'up', 'down' where that the Pauli matrices look as usual.

Moreover, 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.

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
siteslist of Site

The sites defining the local Hilbert space.

Bflatiterable of numpy ndarrays

The matrix defining the MPS on each site, with legs 'p', 'vL', 'vR' (physical, virtual left/right).

SVslist 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’, ‘segmemt’}

MPS boundary conditions. See docstring of MPS.

dtypetype or string

The data type of the array entries. Defaults to the common dtype of Bflat.

permutebool

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 p_state argument 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_LLegCharge | None

Leg charges at bond 0, which are purely conventional. If None, use trivial charges.

Returns
mpsMPS

An MPS with the matrices Bflat converted to npc arrays.

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
siteslist of Site

The sites defining the local Hilbert space.

psiArray

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.

cutofffloat

Cutoff of singular values used in the SVDs.

normalizebool

Whether the resulting MPS should have ‘norm’ 1.

bc‘finite’ | ‘segment’

Boundary conditions.

outer_SNone | (array, array)

For ‘semgent’ bc the singular values on the left and right of the considered segment, None for ‘finite’ boundary conditions.

Returns
psi_mpsMPS

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

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

Create an MPS of entangled singlets.

Parameters
siteSite

The site defining the local Hilbert space, taken uniformly for all sites.

Lint

The number of sites.

pairslist of (int, int)

Pairs of sites to be entangled; the returned MPS will have a singlet for each pair in pairs.

up, downint | 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.

lonelylist of int

Sites which are not included into a singlet pair.

lonely_stateint | str

The state for the lonely sites.

bc{‘infinite’, ‘finite’, ‘segmemt’}

MPS boundary conditions. See docstring of MPS.

Returns
singlet_mpsMPS

An MPS representing singlets on the specified pairs of sites.

copy(self)[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.

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(self, i, form='B', copy=False, cutoff=1e-16, label_p=None)[source]

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

Parameters
iint

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.

copybool

Whether to return a copy even if form matches the current form.

cutofffloat

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_pNone | 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
BArray

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).

Raises
ValueErrorif self is not in canoncial form and form is not None.
set_B(self, i, B, form='B')[source]

Set B at site i.

Parameters
iint

Index choosing the site.

BArray

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.

get_SL(self, i)[source]

Return singular values on the left of site i

get_SR(self, i)[source]

Return singular values on the right of site i

set_SL(self, i, S)[source]

Set singular values on the left of site i

set_SR(self, i, S)[source]

Set singular values on the right of site i

get_op(self, 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.

iint

Index of the site on which the operator acts.

Returns
opnpc.array

One of the entries in op_list, not copied.

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

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

Parameters
iint

Site index.

nint

Number of sites. The result lives on sites[i:i+n].

cutofffloat

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.

formLfloat

Exponent for the singular values to the left.

formRfloat

Exponent for the singular values to the right.

Returns
thetaArray

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.

convert_form(self, new_form='B')[source]

Tranform 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
ValueErrorif trying to convert from a None form. Use canonical_form() instead!
increase_L(self, new_L=None)[source]

Modify self inplace to enlarge the unit cell.

For an infinite MPS, we have unit cells.

Parameters
new_Lint

New number of sites. Defaults to twice the number of current sites.

group_sites(self, 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
nint

Number of sites to be grouped together.

grouped_sitesNone | list of GroupedSite

The sites grouped together.

See also

group_split

Reverts the grouping.

group_split(self, trunc_par=None)[source]

Modify self inplace to split previously grouped sites.

Parameters
trunc_pardict

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

Returns
trunc_errTruncationError

The error introduced by the truncation for the splitting.

See also

group_sites

Should have been used before to combine sites.

get_grouped_mps(self, blocklen)[source]

contract blocklen subsequent tensors into a single one and return result as a new MPS.

blocklen = number of subsequent sites to be combined.

Returns
new MPS object with bunched sites.
get_total_charge(self, only_physical_legs=False)[source]

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

Parameters
only_physical_legsbool

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
qtotalcharges

The sum of the qtotal of the individual B tensors.

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

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

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_legNone | LegCharge

Desired new virtual leg on the very left. Needs to have the same block strucuture as current leg, but can have shifted charge entries.

vR_legNone | LegCharge

Desired new virtual leg on the very right. Needs to have the same block strucuture as current leg, but can have shifted charge entries. Should be vL_leg.conj() for infinite MPS, if qtotal is not given.

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

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

Consider a bipartition of the sytem 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
nint/float

Selects which entropy to calculate; n=1 (default) is the ususal von-Neumann entanglement entropy.

bondsNone | (iterable of) int

Selects the bonds at which the entropy should be calculated. None defaults to range(0, L+1)[self.nontrivial_bonds].

for_matrix_Sbool

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

Returns
entropies1D ndarray

Entanglement entropies for half-cuts. entropies[j] contains the entropy for a cut at bond bonds[j] (i.e. left to site bonds[j]).

entanglement_entropy_segment(self, 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 acchieved by explicitly calculating the reduced density matrix of A and thus works only for small segments.

Parameters
segmentlist of int

Given a first site i, the region A_i is defined to be [i+j for j in segment].

first_siteNone | (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.

nint | float

Selects which entropy to calculate; n=1 (default) is the ususal von-Neumann entanglement entropy, otherwise the n-th Renyi entropy.

Returns
entropies1D ndarray

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

entanglement_spectrum(self, by_charge=False)[source]

return entanglement energy spectrum.

Parameters
by_chargebool

Wheter we should sort the spectrum on each bond by the possible charges.

Returns
ent_spectrumlist

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.

get_rho_segment(self, 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
segmentiterable of int

Sites for which the reduced density matrix is to be calculated. Assumed to be sorted.

Returns
rhoArray

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

probability_per_charge(self, bond=0)[source]

Return probabilites 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 probabilty that this combination occurs in the given state.

Parameters
bondint

The bond to be considered. The returned charges are summed on the left of this bond.

Returns
charge_values2D array

Columns correspond to the different charges in self.chinfo. Rows are the different charge fluctuations at this bond

probabilities1D array

For each row of charge_values the probablity for these values of charge fluctuations.

average_charge(self, 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
bondint

The bond to be considered. The returned charges are summed over the sites left of bond.

Returns
average_charge1D array

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

charge_variance(self, 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
bondint

The bond to be considered. The returned charges are summed over the sites left of bond.

Returns
average_charge1D array

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

mutinf_two_site(self, 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_rangeint

Maximal distance |i-j| for which the mutual information should be calculated. None defaults to L-1.

nfloat

Selects the entropy to use, see entropy().

Returns
coords2D array

Coordinates for the mutinf array.

mutinf1D array

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

overlap(self, other, charge_sector=0, ignore_form=False, **kwargs)[source]

Compute overlap <self|other>.

Parameters
otherMPS

An MPS with the same physical sites.

charge_sectorNone | charges | 0

Selects the charge sector in which the dominant eigenvector of the TransferMatrix is. 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.

ignore_formbool

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!)

**kwargs :

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

Returns
overlapdtype.type

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.

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

Expectation value <psi|ops|psi>/<psi|psi> of (n-site) operator(s).

Given the MPS in canonical form, it calculates n-site expectation values. For example the contraction for a two-site (n = 2) operator on site i would look like:

|          .--S--B[i]--B[i+1]--.
|          |     |     |       |
|          |     |-----|       |
|          |     | op  |       |
|          |     |-----|       |
|          |     |     |       |
|          .--S--B*[i]-B*[i+1]-.
Parameters
ops(list of) { Array | str }

The operators, for wich the expectation value should be taken, All operators should all have the same number of legs (namely 2 n). If less than self.L operators are given, we repeat them periodically. Strings (like 'Id', 'Sz') are translated into single-site operators defined by sites.

sitesNone | list of int

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

axesNone | (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 operators (n = 1), or (['p0', 'p1', ... 'p{n-1}'], ['p0*', 'p1*', .... 'p{n-1}*']) for n > 1.

Returns
exp_vals1D ndarray

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

Examples

One site examples (n=1):

>>> psi.expectation_value('Sz')
[Sz0, Sz1, ..., Sz{L-1}]
>>> psi.expectation_value(['Sz', 'Sx'])
[Sz0, Sx1, Sz2, Sx3, ... ]
>>> psi.expectation_value('Sz', sites=[0, 3, 4])
[Sz0, Sz3, Sz4]

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*']))
>>> psi.expectation_value(SzSx)
[Sz0Sx1, Sz1Sx2, Sz2Sx3, ... ]   # with len L-1 for finite bc, or L for infinite

Example measuring <psi|SzSx|psi2> 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)]
>>> psi.expectation_value(SzSx_list, range(0, psi.L-1, 2))
[Sz0Sx1, Sz2Sx3, Sz4Sx5, ...]
expectation_value_term(self, term, autoJW=True)[source]

Expectation value <psi|op_{i0}op_{i1}...op_{iN}|psi>/<psi|psi>.

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]--.
|          |     |       |        |        |        |
|          |    op1     op2       |       op3       |
|          |     |       |        |        |        |
|          .--S--B*[i0]--B*[i0+1]-B*[i0+2]-B*[i0+3]-.
Parameters
termlist 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).

autoJWbool

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

Returns
exp_valfloat/complex

The expectation value of the tensorproduct of the given onsite operators, <psi|op_i0 op_i1 ... op_iN |psi>/<psi|psi>, where |psi> is the represented MPS.

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', 'Sz'], i0=2)
>>> assert a == b == c
expectation_value_multi_sites(self, operators, i0)[source]

Expectation value <psi|op0_{i0}op1_{i0+1}...opN_{i0+N}|psi>/<psi|psi>.

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} ....

Parameters
operatorsList of { Array | str }

List of one-site operators. This method calculates the expectation value of the n-sites operator given by their tensor product.

i0int

The left most index on which an operator acts, i.e., operators[i] acts on site i + i0.

Returns
exp_valfloat/complex

The expectation value of the tensorproduct of the given onsite operators, <psi|operators[0]_{i0} operators[1]_{i0+1} ... |psi>/<psi|psi>, where |psi> is the represented MPS.

expectation_value_terms_sum(self, 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 effiency, the term_list is converted to an MPO and the expectation value of the MPO is evaluated.

Note

Due to the way MPO expectation values are evaluated for infinite systems, it works only if all terms in the term_list start within the MPS unit cell.

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_listTermList

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_sumlist 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.

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

Correlation function <psi|op1_i op2_j|psi>/<psi|psi> of single site operators.

Given the MPS in canonical form, it calculates 2-site correlation functions. For examples the contraction for a two-site operator on site i would look like:

|          .--S--B[i]--B[i+1]--...--B[j]---.
|          |     |     |            |      |
|          |     |     |            op2    |
|          |     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    |            |      |
|          |     |      |            |      |     |     |      |            |      |
|          |     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*'.

Parameters
ops1(list of) { Array | str }

First operator of the correlation function (acting after ops2). ops1[x] acts on site sites1[x]. If less than len(sites1) operators are given, we repeat them periodically.

ops2(list of) { Array | str }

Second operator of the correlation function (acting before ops1). ops2[y] acts on site sites2[y]. If less than len(sites2) operators are given, we repeat them periodically.

sites1None | int | list of int

List of site indices; 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.

sites2None | 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.

opstrNone | (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_firstbool

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!)

hermitianbool

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.

Returns
C2D ndarray

The correlation function C[x, y] = <psi|ops1[i] ops2[j]|psi>, 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] = <psi|ops1[i] prod_{i <= r < j} opstr[r] ops2[j]|psi>.

  • For i > j: C[x, y] = <psi|prod_{j <= r < i} opstr[r] ops1[i] ops2[j]|psi>.

  • For i = j: C[x, y] = <psi|ops1[i] ops2[j]|psi>.

The condition <= r is replaced by a strict < r, if str_on_first=False.

norm_test(self)[source]

Check that self is in canonical form.

Returns
norm_error: array, shape (L, 2)

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_errror[i, 1] is the norm-difference of:

|   .--theta[i]---         .--s[i+1]--
|   |    |          vs     |
|   .--theta*[i]--         .--s[i+1]--
canonical_form(self, renormalize=True)[source]

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

Simply calls canonical_form_finite() or canonical_form_infinite().

canonical_form_finite(self, renormalize=True, cutoff=0.0)[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.

cutofffloat | None

Cutoff of singular values used in the SVDs.

Returns
U_L, V_RArray

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'.

canonical_form_infinite(self, 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.

tol_xifloat

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

correlation_length(self, target=1, tol_ev0=1e-08, charge_sector=0)[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
targetint

We look for the target + 1 largest eigenvalues.

tol_ev0float

Print warning if largest eigenvalue deviates from 1 by more than tol_ev0.

charge_sectorNone | charges | 0

Selects the charge sector in which the dominant eigenvector of the TransferMatrix is. 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.

Returns
xifloat | 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 lattice spacing in the MPS language, see the warning above.

add(self, 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
otherMPS

Another MPS of the same length to be added with self.

alpha, betacomplex float

Prefactors for self and other. We calculate alpha * |self> + beta * |other>

cutofffloat | None

Cutoff of singular values used in the SVDs.

Returns
sumMPS

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

U_L, V_RArray

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'.

apply_local_op(self, i, op, unitary=None, renormalize=False, cutoff=1e-13)[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.

Parameters
iint

(Left-most) index of the site(s) on which the operator should act.

opstr | 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.

unitaryNone | 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.

renormalizebool

Whether the final state should keep track of the norm (False, default) or be renormalized to have norm 1 (True).

cutofffloat

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

swap_sites(self, i, swap_op='auto', trunc_par={})[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 sign for fermions.

Parameters
iint

Swap the two sites at positions i and i+1.

swap_opNone | 'auto' | Array

The operator used to swap the phyiscal legs of the two-site wave function theta. For None, just transpose/relabel the legs, for 'auto' also take care of fermionic signs. Alternative give an npc Array which represents the full operator used for the swap. Should have legs ['p0', 'p1', 'p0*', 'p1*'] whith 'p0', 'p1*' contractible.

trunc_pardict

Parameters for truncation, see truncate(). chi_max defaults to max(self.chi).

Returns
trunc_errTruncationError

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

permute_sites(self, perm, swap_op='auto', trunc_par={}, verbose=0)[source]

Applies the permutation perm to the state (inplace).

Parameters
permndarray[ndim=1, int]

The applied permutation, such that psi.permute_sites(perm)[i] = psi[perm[i]] (where [i] indicates the i-th site).

swap_opNone | 'auto' | Array

The operator used to swap the phyiscal legs of a two-site wave function theta, see swap_sites().

trunc_pardict

Parameters for truncation, see truncate(). chi_max defaults to max(self.chi).

verbosefloat

Level of verbosity, print status messages if verbose > 0.

Returns
trunc_errTruncationError

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

compute_K(self, perm, swap_op='auto', trunc_par=None, canonicalize=1e-06, verbose=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 [PollmannTurner2012], see also (the appendix and Fig. 11 in the arXiv version of) [CincioVidal2013].

Parameters
perm1D ndarray | Lattice

Permuation 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_opNone | 'auto' | Array

The operator used to swap the phyiscal legs of a two-site wave function theta, see swap_sites().

trunc_pardict

Parameters for truncation, see truncate(). If not set, chi_max defaults to max(self.chi).

canonicalizefloat

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

verbosefloat

Level of verbosity, print status messages if verbose > 0.

Returns
UArray

Unitary representation of the applied permutation on left Schmidt states.

Wndarray

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 Schmidt values S and momenta K from W.

qLegCharge

LegCharge corresponding to W.

ovcomplex

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_errTruncationError

The error of the represented state introduced by the truncation after swaps when performing the truncation.