BaseMPSExpectationValue

  • full name: tenpy.networks.mps.BaseMPSExpectationValue

  • parent module: tenpy.networks.mps

  • type: class

Inheritance Diagram

Inheritance diagram of tenpy.networks.mps.BaseMPSExpectationValue

Methods

BaseMPSExpectationValue.__init__()

BaseMPSExpectationValue.apply_JW_string_left_of_virt_leg(...)

Apply signs on a virtual MPS leg equivalent to a Jordan-Wigner string on the left.

BaseMPSExpectationValue.correlation_function(...)

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

BaseMPSExpectationValue.expectation_value(ops)

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

BaseMPSExpectationValue.expectation_value_multi_sites(...)

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

BaseMPSExpectationValue.expectation_value_term(term)

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

BaseMPSExpectationValue.expectation_value_terms_sum(...)

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

BaseMPSExpectationValue.get_op(op_list, i)

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

BaseMPSExpectationValue.term_correlation_function_left(...)

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

BaseMPSExpectationValue.term_correlation_function_right(...)

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

BaseMPSExpectationValue.term_list_correlation_function_right(...)

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

class tenpy.networks.mps.BaseMPSExpectationValue[source]

Bases: object

Base class providing unified expectation value framework for MPS and MPSEnvironment.

For general expectation values of operators ‘ops’ between different states <bra|ops|ket> we need to include the left/ right environments LP and RP. These are calculated in MPSEnvironment. For “standard” expectation values <psi|ops|psi>, the environments are trivial identities due to the canonical from.

Subclasses need to have the attributes sites, L, bc, finite. See MPS for details.

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

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

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

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

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

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

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

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

Returns:

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

Warning

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

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

Return type:

1D ndarray

Examples

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

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

One site examples (n=1):

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

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

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

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

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

Expectation value with different bra and ket in an MPSEnvironment:

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

Apply signs on a virtual MPS leg equivalent to a Jordan-Wigner string on the left.

If we conserve the (parity of the) total fermion particle number, each Schmidt state |alpha> on a given bond (here left of site i) has a well-defined fermion parity number, so we can simply transform |alpha> --> (-1)**parity[alpha] |alpha>. The corresponding signs (-1)**parity[alpha] are extracted by charge_to_JW_signs().

Warning

We may loose an overall, global minus sign in the case that some B tensors have non-trivial qtotal!

Parameters:
  • theta (Array) – Tensor with virtual leg

  • virtual_leg_index (int | str) – Index of the virtual leg on the left of which we want to apply the JW string.

  • i (int) – Site index of theta.

expectation_value_multi_sites(operators, i0)[source]

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

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

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

Warning

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

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

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

Returns:

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

Warning

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

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

Return type:

float/complex

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

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

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

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

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

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

Warning

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

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

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

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

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

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

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

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

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

Returns:

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

Warning

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

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

Return type:

2D ndarray

Examples

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

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

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

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

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

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

See also

expectation_value_term

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

term_correlation_function_right

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

term_correlation_function_left

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

expectation_value_term(term, autoJW=True)[source]

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

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

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

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

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

Returns:

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

Warning

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

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

Return type:

float/complex

See also

correlation_function

efficient way to evaluate many correlation functions.

Examples

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

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

This is equivalent to the following expression:

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

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

Warning

This function works only for finite bra and ket and does not include normalization factors.

Parameters:

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

Returns:

  • terms_sum ((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.

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

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

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

Warning

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

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

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

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

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

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

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

Returns:

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

Warning

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

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

Return type:

1D array

See also

correlation_function

varying both i and j at once.

term_list_correlation_function_right

generalization to sums of terms on the left/right.

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

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

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

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

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

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

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

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

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

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

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

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

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

Returns:

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

Warning

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

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

Return type:

1D array

See also

term_correlation_function_right

version for a single term in both term_list_L/R.

get_op(op_list, i)[source]

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

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

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

Returns:

  • op (npc.array) – One of the entries in op_list, not copied.

  • needs_JW (bool) – If the operator needs a JW string. Always False if the entry of op_list is an array.