MPO

class tenpy.networks.mpo.MPO(sites, Ws, bc='finite', IdL=None, IdR=None, max_range=None)[source]

Bases: object

Matrix product operator, finite (MPO) or infinite (iMPO).

Parameters
siteslist of Site

Defines the local Hilbert space for each site.

Wslist of Array

The matrices of the MPO. Should have labels wL, wR, p, p*.

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

Boundary conditions as described in mps. 'finite' requires Ws[0].get_leg('wL').ind_len = 1.

IdL(iterable of) {int | None}

Indices on the bonds, which correpond to ‘only identities to the left’. A single entry holds for all bonds.

IdR(iterable of) {int | None}

Indices on the bonds, which correpond to ‘only identities to the right’.

max_rangeint | np.inf | None

Maximum range of hopping/interactions (in unit of sites) of the MPO. None for unknown.

Attributes
Lint

Number of physical sites; for an iMPO the len of the MPO unit cell.

chinfoChargeInfo

The nature of the charge.

siteslist of Site

Defines the local Hilbert space for each site.

dtypetype

The data type of the _W.

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

Boundary conditions as described in mps. 'finite' requires Ws[0].get_leg('wL').ind_len = 1.

IdLlist of {int | None}

Indices on the bonds (length L`+1), which correpond to ‘only identities to the left’. ``None` for bonds where it is not set. In standard form, this is 0 (except for unset bonds in finite case)

IdRlist of {int | None}

Indices on the bonds (length L`+1), which correpond to ‘only identities to the right’. ``None` for bonds where it is not set. In standard form, this is the last index on the bond (except for unset bonds in finite case).

max_rangeint | np.inf | None

Maximum range of hopping/interactions (in unit of sites) of the MPO. None for unknown.

groupedint

Number of sites grouped together, see group_sites().

_Wlist of Array

The matrices of the MPO. Labels are 'wL', 'wR', 'p', 'p*'.

_valid_bctuple of str

Valid boundary conditions. The same as for an MPS.

Methods

dagger(self)

Return hermition conjugate copy of self.

expectation_value(self, psi[, tol, max_range])

Calculate <psi|self|psi>/<psi|psi>.

from_grids(sites, grids[, bc, IdL, IdR, …])

Initialize an MPO from grids.

get_IdL(self, i)

Return index of IdL at bond to the left of site i.

get_IdR(self, i)

Return index of IdR at bond to the right of site i.

get_W(self, i[, copy])

Return W at site i.

get_full_hamiltonian(self[, maxsize])

extract the full Hamiltonian as a d**L x d**L matrix.

get_grouped_mpo(self, blocklen)

Contract blocklen subsequent tensors into a single one and return result as a new MPO object.

group_sites(self[, n, grouped_sites])

Modify self inplace to group sites.

is_equal(self, other[, eps, max_range])

Check if self and other represent the same MPO to precision eps.

is_hermitian(self[, eps, max_range])

Check if self is a hermitian MPO.

set_W(self, i, W)

Set W at site i.

sort_legcharges(self)

Sort virtual legs by charges.

test_sanity(self)

Sanity check, raises ValueErrors, if something is wrong.

classmethod from_grids(sites, grids, bc='finite', IdL=None, IdR=None, Ws_qtotal=None, leg0=None, max_range=None)[source]

Initialize an MPO from grids.

Parameters
siteslist of Site

Defines the local Hilbert space for each site.

gridslist of list of list of entries

For each site (outer-most list) a matrix-grid (corresponding to wL, wR) with entries being or representing (see grid_insert_ops()) onsite-operators.

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

Boundary conditions as described in mps.

IdL(iterable of) {int | None}

Indices on the bonds, which correpond to ‘only identities to the left’. A single entry holds for all bonds.

IdR(iterable of) {int | None}

Indices on the bonds, which correpond to ‘only identities to the right’.

Ws_qtotal(list of) total charge

The qtotal to be used for each grid. Defaults to zero charges.

leg0LegCharge

LegCharge for ‘wL’ of the left-most W. By default, construct it.

max_rangeint | np.inf | None

Maximum range of hopping/interactions (in unit of sites) of the MPO. None for unknown.

See also

grid_insert_ops

used to plug in entries of the grid.

tenpy.linalg.np_conserved.grid_outer

used for final conversion.

test_sanity(self)[source]

Sanity check, raises ValueErrors, if something is wrong.

property L

Number of physical sites; for an iMPO the len of the MPO unit cell.

property dim

List of local physical dimensions.

property finite

Distinguish MPO vs iMPO.

True for an MPO (bc='finite', 'segment'), False for an iMPO (bc='infinite').

property chi

Dimensions of the virtual bonds.

get_W(self, i, copy=False)[source]

Return W at site i.

set_W(self, i, W)[source]

Set W at site i.

get_IdL(self, i)[source]

Return index of IdL at bond to the left of site i.

May be None.

get_IdR(self, i)[source]

Return index of IdR at bond to the right of site i.

May be None.

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.

sort_legcharges(self)[source]

Sort virtual legs by charges. In place.

The MPO seen as matrix of the wL, wR legs is usually very sparse. This sparsity is captured by the LegCharges for these bonds not being sorted and bunched. This requires a tensordot to do more block-multiplications with smaller blocks. This is in general faster for large blocks, but might lead to a larger overhead for small blocks. Therefore, this function allows to sort the virtual legs by charges.

expectation_value(self, psi, tol=1e-10, max_range=100)[source]

Calculate <psi|self|psi>/<psi|psi>.

For a finite MPS, simply contract the network <psi|self|psi>. For an infinite MPS, it assumes that self is the a of terms, with IdL and IdR defined on each site. Under this assumption, it calculates the expectation value of terms with the left-most non-trivial operator inside the MPO unit cell and returns the average value per site.

Parameters
psiMPS

State for which the expectation value should be taken.

tolfloat

Ignored for finite psi. For infinite MPO containing exponentially decaying long-range terms, stop evaluating further terms if the terms in LP have norm < tol.

max_rangeint

Ignored for finite psi. Contract at most self.L * max_range sites, even if tol is not reached. In that case, issue a warning.

Returns
exp_valfloat/complex

The expectation value of self with respect to the state psi. For an infinite MPS: the density per site.

dagger(self)[source]

Return hermition conjugate copy of self.

is_hermitian(self, eps=1e-10, max_range=None)[source]

Check if self is a hermitian MPO.

Shorthand for self.is_equal(self.dagger(), eps, max_range).

is_equal(self, other, eps=1e-10, max_range=None)[source]

Check if self and other represent the same MPO to precision eps.

To compare them efficiently we view self and other as MPS and compare the overlaps abs(<self|self> + <other|other> - 2 Re(<self|other>)) < eps*(<self|self>+<other|other>)

Parameters
otherMPO

The MPO to compare to.

epsfloat

Precision threshold what counts as zero.

max_rangeNone | int

Ignored for finite MPS; for finite MPS we consider only the terms contained in the sites with indices range(self.L + max_range). None defaults to max_range (or L in case this is infinite or None).

Returns
equalbool

Whether self equals other to the desired precision.

get_grouped_mpo(self, blocklen)[source]

Contract blocklen subsequent tensors into a single one and return result as a new MPO object.

get_full_hamiltonian(self, maxsize=1000000.0)[source]

extract the full Hamiltonian as a d**L x d**L matrix.