MPS¶
full name: tenpy.networks.mps.MPS
parent module:
tenpy.networks.mps
type: class
Inheritance Diagram
Methods
|
Initialize self. |
|
Return an MPS which represents |
|
Apply a local (one or multi-site) operator to self. |
|
Return the average charge for the block on the left of a given bond. |
|
Bring self into canonical ‘B’ form, (re-)calculate singular values. |
|
Bring a finite (or segment) MPS into canonical form (in place). |
|
Bring an infinite MPS into canonical form (in place). |
|
Return the charge variance on the left of a given bond. |
|
Compresss an MPS. |
|
Compress self with a single sweep of SVDs; in place. |
|
Compute the momentum quantum numbers of the entanglement spectrum for 2D states. |
|
Tranform self into different canonical form (by scaling the legs with singular values). |
|
Returns a copy of self. |
|
Correlation function |
|
Calculate the correlation length by diagonalizing the transfer matrix. |
|
Repeat the unit cell for infinite MPS boundary conditions; in place. |
|
Calculate the (half-chain) entanglement entropy for all nontrivial bonds. |
|
Calculate entanglement entropy for general geometry of the bipartition. |
|
return entanglement energy spectrum. |
|
Expectation value |
|
Expectation value |
|
Expectation value |
|
Calculate expectation values for a bunch of terms and sum them up. |
|
Construct a matrix product state from a set of numpy arrays Bflat and singular vals. |
|
Construct an MPS from a single tensor psi with one leg per physical site. |
|
Load instance from a HDF5 file. |
|
Construct an MPS from a product state given in lattice coordinates. |
|
Construct a matrix product state from a given product state. |
|
Create an MPS of entangled singlets. |
|
Gauge the legcharges of the virtual bonds such that the MPS has a total qtotal. |
|
Return (view of) B at site i in canonical form. |
|
Return singular values on the left of site i |
|
Return singular values on the right of site i |
|
contract blocklen subsequent tensors into a single one and return result as a new MPS. |
|
Given a list of operators, select the one corresponding to site i. |
|
Return reduced density matrix for a segment. |
|
Calculates the n-site wavefunction on |
|
Calculate and return the qtotal of the whole MPS (when contracted). |
|
Modify self inplace to group sites. |
|
Modify self inplace to split previously grouped sites. |
|
Modify self inplace to enlarge the MPS unit cell; in place. |
|
Calculate the two-site mutual information \(I(i:j)\). |
Check that self is in canonical form. |
|
|
Compute overlap |
|
Applies the permutation perm to the state (inplace). |
|
Return probabilites of charge value on the left of a given bond. |
|
Shift the section we define as unit cellof an infinite MPS; in place. |
|
Export self into a HDF5 file. |
|
Set B at site i. |
|
Set singular values on the left of site i |
|
Set singular values on the right of site i |
|
SVD a two-site wave function theta and save it in self. |
Perform a spatial inversion along the MPS. |
|
|
Swap the two neighboring sites i and i+1 (inplace). |
|
Correlation function between (multi-site) terms, moving the left term, fix right term. |
|
Correlation function between (multi-site) terms, moving the right term, fix left term. |
Sanity check, raises ValueErrors, if something is wrong. |
Class Attributes and Properties
Number of physical sites; for an iMPS the len of the MPS unit cell. |
|
Dimensions of the (nontrivial) virtual bonds. |
|
List of local physical dimensions. |
|
Distinguish MPS vs iMPS. |
|
Slice of the non-trivial bond indices, depending on |
-
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
sites (list of
Site
) – Defines the local Hilbert space for each site.Bs (list of
Array
) – The ‘matrices’ of the MPS. Labels arevL, 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 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.
-
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. Forform = (nuL, nuR)
, the stored_B[i]
ares**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
-
norm
¶ The norm of the state, i.e.
sqrt(<psi|psi>)
. Ignored for (normalized)expectation_value()
, but important foroverlap()
.- Type
-
grouped
¶ Number of sites grouped together, see
group_sites()
.- Type
-
_B
¶ The ‘matrices’ of the MPS. Labels are
vL, vR, p
(in any order). We recommend usingget_B()
andset_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 beNone
if the MPS is not in canonical form. Otherwise,_S[i]
is to the left of_B[i]
. We recommend usingget_SL()
,get_SR()
,set_SL()
,set_SR()
, which takes proper care of the boundary conditions.- Type
list of (
None
| 1D array)
-
_valid_forms
¶ Class attribute. Mapping for canonical forms to a tuple
(nuL, nuR)
indicating thatself._Bs[i] = s[i]**nuL -- Gamma[i] -- s[i]**nuR
is saved.- Type
-
_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
-
_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 thePurification_MPS
if just the number of physical legs changed.- Type
list of str
-
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"
, andform
converted to a single array of shape (L, 2) as"canonical_form"
, Moreover, it savesnorm
,L
,grouped
and_transfermatrix_keep
(as “transfermatrix_keep”) as HDF5 attributes, as well as the maximum ofchi
under the name “max_bond_dimension”.
-
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, **kwargs)[source]¶ Construct an MPS from a product state given in lattice coordinates.
This is a wrapper around
from_product_state()
. The purpuse 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 shapelat.shape
, if it is smaller.**kwargs – Other keyword arguments as definied in
from_product_state()
. bc is set by default fromlat.bc_MPS
.
- Returns
product_mps – An MPS representing the specified product state.
- Return type
Examples
Let’s first consider a
Ladder
composed of aSpinHalfSite
and aFermionSite
.>>> spin_half = tenpy.networks.site.SpinHalfSite() >>> fermion = tenpy.networks.site.FermionSite() >>> ladder_i = tenpy.models.lattice.Ladder(2, [spin_half, fermion], bc_MPS="infinite")
To initialize a state of up-spins on the spin sites and half-filled ferions, you can use:
>>> p_state = [["up", "empty"], ["up", "full"]] >>> psi = tenpy.networks.MPS.from_lat_product_state(ladder_i, p_state)
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:>>> ladder_f = tenpy.models.lattice.Ladder(10, [spin_half, fermion], bc_MPS="finite") >>> psi = tenpy.networks.MPS.from_lat_product_state(ladder_f, p_state)
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).>>> honeycomb = tenpy.models.lattice.Honeycomb([4, 4], [fermion, fermion], bc_MPS="finite") >>> p_state_only_A = [[["empty", "full"]]] >>> psi_only_A = tenpy.networks.MPS.from_lat_product_state(honeycomb, p_state_only_A) >>> p_state_stripe_x = [[["empty", "empty"], ... ["full", "full"]]] >>> psi_stripe_x = tenpy.networks.MPS.from_lat_product_state(honeycomb, p_state_stripe_x) >>> p_state_stripe_y = [[["empty", "empty"]], ... [["full", "full"]]] >>> psi_stripe_y = tenpy.networks.MPS.from_lat_product_state(honeycomb, p_state_stripe_y)
-
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', 'segmemt'}) – 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’sperm
. 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
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 theSite
, in this example aSpinHalfSite
.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 thebloch_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 generalSpinSite
, where the states are orderd ascending from'down'
to'up'
. TheSpinHalfSite
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
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 ofnontrivial_bonds
are ignored.bc ({'infinite', 'finite', 'segmemt'}) – 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’sperm
. 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_L (LegCharge |
None
) – Leg charges at bond 0, which are purely conventional. IfNone
, use trivial charges.
- Returns
mps – An MPS with the matrices Bflat converted to npc arrays.
- Return type
-
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 ‘semgent’ 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
-
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.
up (int | str) – A singlet is defined as
(|up down> - |down up>)/2**0.5
,up
anddown
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
anddown
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', 'segmemt'}) – MPS boundary conditions. See docstring of
MPS
.
- Returns
singlet_mps – An MPS representing singlets on the specified pairs of sites.
- Return type
-
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. ForNone
, 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 (ifcopy=False
), or a copy (if the form changed orcopy=True
).- Return type
:raises ValueError : if self is not in canoncial 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
. IfNone
, no truncation is done.update_norm (bool) – If
True
, multiply the norm of theta intonorm
.
-
get_op
(op_list, i)[source]¶ Given a list of operators, select the one corresponding to site i.
- Parameters
- Returns
op – One of the entries in op_list, not copied.
- Return type
npc.array
-
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
-
convert_form
(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 ValueError : if trying to convert from a
None
form. Usecanonical_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 ofpsi.increase_L(new_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 cellof an infinite MPS; in place.
Suppose we have a unit cell with tensors
[A, B, C, D]
(repeated on both sites). Withshift = 1
, the new unit cell will be[D, A, B, C]
, whereasshift = -1
will give[B, C, D, A]
.- Parameters
shift (int) – By how many sites to move the tensors to the right.
-
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
See also
group_sites()
Should have been used before to combine sites.
-
get_grouped_mps
(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
- Return type
new MPS object with bunched sites.
-
get_total_charge
(only_physical_legs=False)[source]¶ Calculate and return the qtotal of the whole MPS (when contracted).
- 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.
- 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. Needs to have the same block strucuture as current leg, but can have shifted charge entries.
vR_leg (None | 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
(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
n (int/float) – Selects which entropy to calculate; n=1 (default) is the ususal von-Neumann entanglement entropy.
bonds (
None
| (iterable of) int) – Selects the bonds at which the entropy should be calculated.None
defaults torange(0, L+1)[self.nontrivial_bonds]
.for_matrix_S (bool) – 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
entropies – Entanglement entropies for half-cuts. entropies[j] contains the entropy for a cut at bond
bonds[j]
(i.e. left to sitebonds[j]
).- 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 acchieved by explicitly calculating the reduced density matrix of A and thus works only for small segments.
- 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 torange(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 ususal von-Neumann entanglement entropy, otherwise the n-th Renyi entropy.
- Returns
entropies –
entropies[i]
contains the entropy for the the regionA_i
defined above.- Return type
1D ndarray
-
entanglement_spectrum
(by_charge=False)[source]¶ return entanglement energy spectrum.
- Parameters
by_charge (bool) – Wheter 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
-
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*'
withk=len(segment)
.- Return type
-
probability_per_charge
(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
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 probablity 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>\).
-
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\).
-
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
- Returns
coords (2D array) – Coordinates for the mutinf array.
mutinf (1D array) –
mutinf[k]
is the mutual information \(I(i:j)\) between the sitesi, j = coords[k]
.
-
overlap
(other, charge_sector=None, ignore_form=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 formform
at each site. IfTrue
, 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
overlap – The contraction
<self|other> * self.norm * other.norm
(i.e., taking into account thenorm
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
(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 bysites
.sites (None | 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.)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 operators (n = 1), or(['p0', 'p1', ... 'p{n-1}'], ['p0*', 'p1*', .... 'p{n-1}*'])
for n > 1.
- Returns
exp_vals – Expectation values,
exp_vals[i] = <psi|ops[i]|psi>
, whereops[i]
acts on site(s)j, j+1, ..., j+{n-1}
withj=sites[i]
.- Return type
1D ndarray
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
(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
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,
<psi|op_i0 op_i1 ... op_iN |psi>/<psi|psi>
, where|psi>
is the represented MPS.- 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', 'Sz'], i0=2) >>> assert a == b == c
-
expectation_value_multi_sites
(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} ...
, looking like this (with op short for operators, forlen(operators)=3
):.–S–B[i0]—B[i0+1]–B[i0+2]–B[i0+3]–.| | | | | || op[0] op[1] op[2] op[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
- Returns
exp_val – 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.- Return type
float/complex
-
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 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
- 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 MPSterms_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
(ops1, ops2, sites1=None, sites2=None, opstr=None, str_on_first=True, hermitian=False, autoJW=True)[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. Forstr_on_first=False
, the opstr on sitemin(i, j)
is always left out.Strings (like
'Id', 'Sz'
) in the arguments are translated into single-site operators defined by theSite
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 and while fixing the right end. In that case, you might be better of (=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.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 sitesrange(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 sitesrange(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 betweenops1
andops2
. If less thanL
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: choosestr_on_first=True
for fermions!)hermitian (bool) – Optimization flag: if
sites1 == sites2
andOps1[i]^\dagger == Ops2[i]
(which is not checked explicitly!), the resultingC[x, y]
will be hermitian. We can use that to avoid calculations, sohermitian=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] = <psi|ops1[i] ops2[j]|psi>
, whereops1[i]
acts on sitei=sites1[x]
andops2[j]
on sitej=sites2[y]
. If opstr is given, it gives (forstr_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
, ifstr_on_first=False
.- Return type
2D ndarray
Examples
For a spin chain:
>>> psi.correlation_function("A", "B") [[A0B0, A0B1, ..., A0B{L-1}], [A1B0, A1B1, ..., A1B{L-1]], ..., [A{L-1}B0, ALB1, ..., A{L-1}B{L-1}], ]
To evaluate the correlation function for a single i, you can use
sites1=[i]
:>>> psi.correlation_function("A", "B", [3]) [[A3B0, A3B1, ..., A3B{L-1}]]
Alternatively, you can use
term_correlation_function_right()
(orterm_correlation_function_left()
):>>> corr1 = psi.correlation_function("A", "B", [0], range(1, 10)) >>> corr2 = psi.term_correlation_function_right([("A", 0), [("B", 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:
>>> CdC = psi.correlation_function("Cd", "C") # optionally: use `hermitian=True` >>> psi.correlation_function("C", "Cd")[1, 2] == -CdC[1, 2] 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.
-
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)]
andterm_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.- 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 torange(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.
- Return type
1D array
See also
correlation_function()
varying both i and j at once.
-
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.
For
term_L = [('A', 0), ('B', 1)]
andterm_R = [('C', 0), ('D', 1)]
, calculate the correlation function \(A_{i+0} B_{i+1} C_{j+0} D_{j+1}\) for varying i and fixed 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.- 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 (list of int | None) – List of offsets to be added to the indices of term_L. Is sorted descending before use, i.e., the order is ignored. 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.j_R (int) – Offset added to the indices of term_R.
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.
- Return type
1D array
See also
correlation_function()
varying both i and j at once.
-
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_errror[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
(renormalize=True)[source]¶ Bring self into canonical ‘B’ form, (re-)calculate singular values.
Simply calls
canonical_form_finite()
orcanonical_form_infinite()
.
-
canonical_form_finite
(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
- 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
-
canonical_form_infinite
(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.
-
correlation_length
(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
target (int) – We look for the target + 1 largest eigenvalues.
tol_ev0 (float) – Print warning if largest eigenvalue deviates from 1 by more than tol_ev0.
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 zero-charge sector. Defaults to0
, i.e., assumes the dominant eigenvector is in charge sector 0.
- Returns
xi – 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.
- Return type
float | 1D array
-
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
-
apply_local_op
(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
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 bysites
.unitary (None | bool) – Whether op is unitary, i.e., whether the canonical form is preserved (
True
) or whether we should callcanonical_form()
(False
).None
checks whethernorm(op dagger(op) - identity)
is smaller than cutoff.renormalize (bool) – Whether the final state should keep track of the norm (False, default) or be renormalized to have norm 1 (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.)
-
swap_sites
(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
i (int) – Swap the two sites at positions i and i+1.
swap_op (
None
|'auto'
|Array
) – The operator used to swap the phyiscal legs of the two-site wave function theta. ForNone
, just transpose/relabel the legs, for'auto'
also take care of fermionic signs. Alternative give an npcArray
which represents the full operator used for the swap. Should have legs['p0', 'p1', 'p0*', 'p1*']
whith'p0', 'p1*'
contractible.trunc_par (dict) – Parameters for truncation, see
truncation
. Defaults to{'chi_max': max(self.chi)}
.
- Returns
trunc_err – The error of the represented state introduced by the truncation after the swap.
- Return type
-
permute_sites
(perm, swap_op='auto', trunc_par={}, verbose=0)[source]¶ Applies the permutation perm to the state (inplace).
- 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'
|Array
) – The operator used to swap the phyiscal legs of a two-site wave function theta, seeswap_sites()
.trunc_par (dict) – Parameters for truncation, see
truncation
. Defaults to{'chi_max': max(self.chi)}
.verbose (float) – Level of verbosity, print status messages if verbose > 0.
- Returns
trunc_err – The error of the represented state introduced by the truncation after the swaps.
- Return type
-
compute_K
(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 [pollmann2012], see also (the appendix and Fig. 11 in the arXiv version of) [cincio2013].
- Parameters
perm (1D ndarray |
Lattice
) – Permuation to be applied to the physical indices, seepermute_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 aCouplingModel
, give its lat attribute for this argument)swap_op (
None
|'auto'
|Array
) – The operator used to swap the phyiscal legs of a two-site wave function theta, seeswap_sites()
.trunc_par (dict) – Parameters for truncation, see
truncation
. If not set, chi_max defaults tomax(self.chi)
.canonicalize (float) – Check that self is in canonical form; call
canonical_form()
ifnorm_test()
yieldsnp.linalg.norm(self.norm_test()) > canonicalize
.verbose (float) – Level of verbosity, print status messages if verbose > 0.
- 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 usenp.abs()
andnp.angle()
to extract the 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]¶ Compresss an MPS.
Options
-
config
MPS_compress
¶ option summary Whether to combine legs into pipes. This combines the virtual and [...]
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 [...]
N_sweeps (from VariationalCompression) in VariationalCompression
Number of sweeps to perform.
orthogonal_to (from Sweep) in DMRGEngine.init_env
List of other matrix product states to orthogonalize against. [...]
Number of sweeps to be performed without optimization to update [...]
Number of sweeps that have already been performed.
Truncation parameters as described in :cfg:config:`truncation`.
Level of verbosity (i.e. how much status information to print); higher=more [...]
-
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
.
-
option
-
config
-
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
.