seemps.state.CanonicalMPS#

class seemps.state.CanonicalMPS#

Bases: MPS

Canonical MPS class.

This implements a Matrix Product State object with open boundary conditions, that is always on canonical form with respect to a given site. The tensors have three indices, A[α,i,β], where α,β are the internal labels and i is the physical state of the given site.

Parameters:
dataIterable[Tensor3]

A set of tensors that will be orthogonalized. It can be an MPS state.

centerint, optional

The center for the canonical form. Defaults to the first site center = 0.

normalizebool, optional

Whether to normalize the state to compensate for truncation errors. Defaults to the value set by strategy.

strategyStrategy, optional

The truncation strategy for the orthogonalization and later algorithms. Defaults to DEFAULT_STRATEGY.

Renyi_entropy(site: int | None = None, alpha: float = 2.0) float[source]#

Compute the Renyi entropy of the MPS for a bipartition around site.

Parameters:
siteint, optional

Site in the range [0, self.size), defaulting to self.center. The system is diveded into [0, self.site) and [self.site, self.size).

alphafloat, default = 2

Power of the Renyi entropy.

Returns:
float

Von Neumann entropy of bipartition.

Schmidt_weights(site: int | None = None) Vector[source]#

Return the Schmidt weights for a bipartition around site.

Parameters:
siteint, optional

Site in the range [0, self.size), defaulting to self.center. The system is diveded into [0, self.site) and [self.site, self.size).

Returns:
numbers: np.ndarray

Vector of non-negative Schmidt weights.

copy()[source]#

Return a shallow copy of the CanonicalMPS, preserving the tensors.

entanglement_entropy(site: int | None = None) float[source]#

Compute the entanglement entropy of the MPS for a bipartition around site.

Parameters:
siteint, optional

Site in the range [0, self.size), defaulting to self.center. The system is diveded into [0, self.site) and [self.site, self.size).

Returns:
float

Von Neumann entropy of bipartition.

classmethod from_vector(ψ: VectorLike, dimensions: Sequence[int], strategy: Strategy = DEFAULT_STRATEGY, normalize: bool = True, center: int = 0) CanonicalMPS[source]#

Create an MPS in canonical form starting from a state vector.

Parameters:
ψVectorLike

Real or complex vector of a wavefunction.

dimensionsSequence[int]

Sequence of integers representing the dimensions of the quantum systems that form this state.

strategyStrategy, default = DEFAULT_STRATEGY

Default truncation strategy for algorithms working on this state.

normalizebool, default = True

Whether the state is normalized to compensate truncation errors.

centerint, default = 0

Center for the canonical form of this decomposition.

Returns:
CanonicalMPS

A valid matrix-product state approximating this state vector.

See also

from_vector()
left_environment(site: int) ndarray[tuple[Any, ...], dtype[_ScalarT]][source]#

Optimized version of left_environment()

norm_squared() float[source]#

Norm-2 squared \(\Vert{\psi}\Vert^2\) of this MPS.

normalize_inplace() CanonicalMPS[source]#

Normalize the state by updating the central tensor.

recenter(center: int, strategy: TypeAliasForwardRef('~seemps.state.Strategy') | None = None) CanonicalMPS[source]#

Update destructively the state to be in canonical form with respect to a different site.

Parameters:
centerint

The new site for orthogonalization in [0, self.size)

strategyStrategy, optional

Truncation strategy. Defaults to self.strategy

Returns:
CanonicalMPS

This same object.

right_environment(site: int) ndarray[tuple[Any, ...], dtype[_ScalarT]][source]#

Optimized version of right_environment()

update_2site_left(AA: ndarray[tuple[Any, ...], dtype[_ScalarT]], site: int, strategy: Strategy) None[source]#

Split a two-site tensor into two one-site tensors by left orthonormalization and insert the tensor in canonical form into the MPS Ψ at the given site and the site on its right. Update the neighboring sites in the process.

Parameters:
AATensor4

Two-site tensor A[a,i,j,b]

siteint

The index of the site whose quantum number is i. The new center will be self.site.

strategyStrategy

Truncation strategy, including relative tolerances and maximum bond dimensions

update_2site_right(AA: ndarray[tuple[Any, ...], dtype[_ScalarT]], site: int, strategy: Strategy) None[source]#

Split a two-site tensor into two one-site tensors by right orthonormalization and insert the tensor in canonical form into the MPS at the given site and the site on its right. Update the neighboring sites in the process.

Parameters:
AATensor4

Two-site tensor A[a,i,j,b]

siteint

The index of the site whose quantum number is i. The new center will be self.site+1.

strategyStrategy

Truncation strategy, including relative tolerances and maximum bond dimensions

update_canonical(A: ndarray[tuple[Any, ...], dtype[_ScalarT]], direction: int, truncation: Strategy) None[source]#

Update the state, replacing the tensor at self.center and moving the center to self.center + direction.

Parameters:
ATensor3

The new tensor.

direction{ +1, -1 }

Direction in which the update is performed.

truncationStrategy

Truncation parameters such as tolerance or maximum bond dimension.

Returns:
float

The truncation error of this update.

zero_state() CanonicalMPS[source]#

Return a zero wavefunction with the same physical dimensions.