MPSEnvironment
full name: tenpy.networks.mps.MPSEnvironment
parent module:
tenpy.networks.mps
type: class
Inheritance Diagram
Methods
|
|
Apply signs on a virtual MPS leg equivalent to a Jordan-Wigner string on the left. |
|
Update short_term_keys for the cache and possibly preload tensors. |
|
Delete all partial contractions except the left-most LP and right-most RP. |
|
|
Correlation function |
Delete stored part strictly to the left of site i. |
|
Delete stored part strictly to the right of site i. |
|
|
Expectation value |
Expectation value |
|
Expectation value |
|
Calculate expectation values for a bunch of terms and sum them up. |
|
Calculate the overlap by a full contraction of the network. |
|
|
Calculate LP at given site from nearest available one. |
Return number of physical sites in the contractions of get_LP(i). |
|
|
Calculate RP at given site from nearest available one. |
Return number of physical sites in the contractions of get_RP(i). |
|
Return data for (re-)initialization of the environment. |
|
|
Given a list of operators, select the one corresponding to site i. |
Return True if LP left of site i is stored. |
|
Return True if RP right of site i is stored. |
|
|
Build initial left part |
|
Build initial right part |
(Re)initialize first LP and last RP from the given data. |
|
|
Store part to the left of site i. |
|
Store part to the right of site i. |
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. |
|
Correlation function between sums of multi-site terms, moving the right sum of term. |
|
Sanity check, raises ValueErrors, if something is wrong. |
- class tenpy.networks.mps.MPSEnvironment(bra, ket, cache=None, **init_env_data)[source]
Bases:
BaseEnvironment
,BaseMPSExpectationValue
Class storing partial contractions between two different MPS and providing expectation values.
The network for a contraction \(<bra|Op|ket>\) of a local operator Op, say exemplary at sites i, i+1 looks like:
| .-----M[0]--- ... --M[1]---M[2]--- ... ->--. | | | | | | | | | |------| | | LP[0] | | Op | RP[-1] | | | |------| | | | | | | | | .-----N[0]*-- ... --N[1]*--N[2]*-- ... -<--.
This class stores the partial contractions LP and RP coming from the left and right up to each bond, allowing efficient calculations of expectation values, correlation functions etc.
We use the following label convention (where arrows indicate qconj):
| .-->- vR vL ->-. | | | | LP RP | | | | .--<- vR* vL* -<-.
- full_contraction(i0)[source]
Calculate the overlap by a full contraction of the network.
The full contraction of the environments gives the overlap
<bra|ket>
, taking into account theMPS.norm
of both bra and ket. For this purpose, this function contractsget_LP(i0+1, store=False)
andget_RP(i0, store=False)
with appropriate singular values in between.- Parameters:
i0 (int) – Site index.
- 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 bycharge_to_JW_signs()
.Warning
We may loose an overall, global minus sign in the case that some B tensors have non-trivial qtotal!
- cache_optimize(short_term_LP=[], short_term_RP=[], preload_LP=None, preload_RP=None)[source]
Update short_term_keys for the cache and possibly preload tensors.
- Parameters:
short_term_LP (list of int) – i indices for
get_LP()
andget_RP()
, respectively, for which a repeated look-up could happen, i.e., for which tensors should be kept in RAM until the next call to this function.short_term_RP (list of int) – i indices for
get_LP()
andget_RP()
, respectively, for which a repeated look-up could happen, i.e., for which tensors should be kept in RAM until the next call to this function.preload_LP (int | None) – If not None, preload the tensors for the corresponding
get_LP()
andget_RP()
call, respectively, from disk.preload_RP (int | None) – If not None, preload the tensors for the corresponding
get_LP()
andget_RP()
call, respectively, from disk.
- 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 andstr_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. 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 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 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] = <bra|ops1[i] ops2[j]|ket>
, whereops1[i]
acts on sitei=sites1[x]
andops2[j]
on sitej=sites2[y]
. If opstr is given, it gives (forstr_on_first=True
): - Fori < j
:C[x, y] = <bra|ops1[i] prod_{i <= r < j} opstr[r] ops2[j]|ket>
. - Fori > j
:C[x, y] = <bra|prod_{j <= r < i} opstr[r] ops1[i] ops2[j]|ket>
. - Fori = j
:C[x, y] = <bra|ops1[i] ops2[j]|ket>
. The condition<= r
is replaced by a strict< r
, ifstr_on_first=False
.Warning
The
MPSEnvironment
variant of this method takes the accumulated MPSnorm
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 useterm_correlation_function_right()
(orterm_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` >>> bool(psi.correlation_function("C", "Cd")[1, 2] == -CdC[2, 1]) True >>> bool(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(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 thanlen(sites)
operators are given, we repeat them periodically. Strings (like'Id', 'Sz'
) are translated into single-site operators defined bysites
.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>
, whereops[i]
acts on site(s)j, j+1, ..., j+{n-1}
withj=sites[i]
.Warning
The
MPSEnvironment
variant of this method takes the accumulated MPSnorm
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]
- 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, forlen(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:
- 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 MPSnorm
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
- 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 MPSnorm
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 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.
- get_LP(i, store=True)[source]
Calculate LP at given site from nearest available one.
The returned
LP_i
corresponds to the following contraction, where the M’s and the N’s are in the ‘A’ form:| .-------M[0]--- ... --M[i-1]--->- 'vR' | | | | | LP[0] | | | | | | | .-------N[0]*-- ... --N[i-1]*--<- 'vR*'
- get_LP_age(i)[source]
Return number of physical sites in the contractions of get_LP(i).
Might be
None
.
- get_RP(i, store=True)[source]
Calculate RP at given site from nearest available one.
The returned
RP_i
corresponds to the following contraction, where the M’s and the N’s are in the ‘B’ form:| 'vL' ->---M[i+1]-- ... --M[L-1]----. | | | | | | | RP[L-1] | | | | | 'vL*' -<---N[i+1]*- ... --N[L-1]*---.
- get_RP_age(i)[source]
Return number of physical sites in the contractions of get_RP(i).
Might be
None
.
- get_initialization_data(first=0, last=None, include_bra=False, include_ket=False)[source]
Return data for (re-)initialization of the environment.
- Parameters:
first (int) – The first and last site, to the left and right of which we should return the environments. Defaults to 0 and
L
- 1.last (int) – The first and last site, to the left and right of which we should return the environments. Defaults to 0 and
L
- 1.include_bra (bool) – Whether to also return the
bra
andket
, respectively.include_ket (bool) – Whether to also return the
bra
andket
, respectively.
- Returns:
init_env_data – A dictionary with the following entries.
- init_LP, init_RP
Array
LP on the left of site first and RP on the right of site last, which can be used as init_LP and init_RP for the initialization of a new environment.
- age_LP, age_RPint
The number of physical sites involved into the contraction yielding init_LP and init_RP, respectively.
- bra, ket
MPS
References of
bra
andket
. Only included if include_bra and include_ket are True, respectively.
- init_LP, init_RP
- Return type:
- get_op(op_list, i)[source]
Given a list of operators, select the one corresponding to site i.
- Parameters:
- 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 ofop_list
is an array.
- init_LP(i, start_env_sites=0)[source]
Build initial left part
LP
.If bra and ket are the same and in left canonical form, this is the environment you get contracting the overlaps from the left infinity up to bond left of site i.
For segment MPS, the
segment_boundaries
are read out (if set).
- init_RP(i, start_env_sites=0)[source]
Build initial right part
RP
for an MPS/MPOEnvironment.If bra and ket are the same and in right canonical form, this is the environment you get contracting from the right infinity up to bond right of site i.
For segment MPS, the
segment_boundaries
are read out (if set).
- init_first_LP_last_RP(init_LP=None, init_RP=None, age_LP=0, age_RP=0, start_env_sites=0)[source]
(Re)initialize first LP and last RP from the given data.
- Parameters:
init_LP (
None
|Array
) – Initial very left partLP
. IfNone
, build one with :meth`init_LP`.init_RP (
None
|Array
) – Initial very right partRP
. IfNone
, build one withinit_RP()
.age_LP (int) – The number of physical sites involved into the contraction of init_LP.
age_RP (int) – The number of physical sites involved into the contraction of init_RP.
start_env_sites (int) – If init_LP and init_RP are not specified, contract each start_env_sites for them.
- 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_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.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 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.
Warning
The
MPSEnvironment
variant of this method takes the accumulated MPSnorm
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_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 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 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 MPSnorm
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.