mps

  • full name: tenpy.networks.mps

  • parent module: tenpy.networks

  • type: module

Classes

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

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

MPSEnvironment(bra, ket[, init_LP, init_RP, …])

Stores partial contractions of \(<bra|Op|ket>\) for local operators Op.

TransferMatrix(bra, ket[, shift_bra, …])

Transfer matrix of two MPS (bra & ket).

Functions

build_initial_state(size, states, filling[, …])

Build an “initial state” list.

Module description

This module contains a base class for a Matrix Product State (MPS).

An MPS looks roughly like this:

|   -- B[0] -- B[1] -- B[2] -- ...
|       |       |      |

We use the following label convention for the B (where arrows indicate qconj):

|  vL ->- B ->- vR
|         |
|         ^
|         p

We store one 3-leg tensor _B[i] with labels 'vL', 'vR', 'p' for each of the L sites 0 <= i < L. Additionally, we store L+1 singular value arrays _S[ib] on each bond 0 <= ib <= L, independent of the boundary conditions. _S[ib] gives the singlur values on the bond i-1, i. However, be aware that e.g. chi returns only the dimensions of the nontrivial_bonds depending on the boundary conditions.

The matrices and singular values always represent a normalized state (i.e. np.linalg.norm(psi._S[ib]) == 1 up to roundoff errors), but (for finite MPS) we keep track of the norm in norm (which is respected by overlap(), …).

Valid MPS boundary conditions (not to confuse with bc_coupling of tenpy.models.model.CouplingModel) are the following:

bc

description

‘finite’

Finite MPS, G0 s1 G1 ... s{L-1} G{l-1}. This is acchieved by using a trivial left and right bond s[0] = s[-1] = np.array([1.]).

‘segment’

Generalization of ‘finite’, describes an MPS embedded in left and right environments. The left environment is described by chi[0] orthonormal states which are weighted by the singular values s[0]. Similar, s[L] weight some right orthonormal states. You can think of the left and right states to be generated by additional MPS, such that the overall structure is something like ... s L s L [s0 G0 s1 G1 ... s{L-1} G{L-1} s{L}] R s R s R ... (where we save the part in the brackets [ ... ] ).

‘infinite’

infinite MPS (iMPS): we save a ‘MPS unit cell’ [s0 G0 s1 G1 ... s{L-1} G{L-1}] which is repeated periodically, identifying all indices modulo self.L. In particular, the last bond L is identified with 0. (The MPS unit cell can differ from a lattice unit cell). bond is identified with the first one.

An MPS can be in different ‘canonical forms’ (see [Vidal2004], [Schollwoeck2011]). To take care of the different canonical forms, algorithms should use functions like get_theta(), get_B() and set_B() instead of accessing them directly, as they return the B in the desired form (which can be chosen as an argument). The values of the tuples for the form correspond to the exponent of the singular values on the left and right. To keep track of a “mixed” canonical form A A A s B B, we save the tuples for each site of the MPS in MPS.form.

form

tuple

description

'B'

(0, 1)

right canonical: _B[i] = -- Gamma[i] -- s[i+1]-- The default form, which algorithms asssume.

'C'

(0.5, 0.5)

symmetric form: _B[i] = -- s[i]**0.5 -- Gamma[i] -- s[i+1]**0.5--

'A'

(1, 0)

left canonical: _B[i] = -- s[i] -- Gamma[i] --.

'G'

(0, 0)

Save only _B[i] = -- Gamma[i] --.

'Th'

(1, 1)

Form of a local wave function theta with singular value on both sides. psi.get_B(i, 'Th') is equivalent to ``psi.get_theta(i, n=1).

None

None

General non-canoncial form. Valid form for initialization, but you need to call canonical_form() (or similar) before using algorithms.