Library documentation

States

MPSKit.FiniteMPSType
FiniteMPS{A<:GenericMPSTensor,B<:MPSBondTensor} <: AbstractFiniteMPS

Type that represents a finite Matrix Product State.

Fields

  • ALs – left-gauged MPS tensors
  • ARs – right-gauged MPS tensors
  • ACs – center-gauged MPS tensors
  • CLs – gauge tensors

Where each is entry can be a tensor or missing.

Notes

By convention, we have that:

  • AL[i] * CL[i+1] = AC[i] = CL[i] * AR[i]
  • AL[i]' * AL[i] = 1
  • AR[i] * AR[i]' = 1

Constructors

FiniteMPS([f, eltype], physicalspaces::Vector{<:Union{S,CompositeSpace{S}}},
          maxvirtualspaces::Union{S,Vector{S}};
          normalize=true, left=oneunit(S), right=oneunit(S)) where {S<:ElementarySpace}
FiniteMPS([f, eltype], N::Int, physicalspace::Union{S,CompositeSpace{S}},
          maxvirtualspaces::Union{S,Vector{S}};
          normalize=true, left=oneunit(S), right=oneunit(S)) where {S<:ElementarySpace}
FiniteMPS(As::Vector{<:GenericMPSTensor}; normalize=false, overwrite=false)

Construct an MPS via a specification of physical and virtual spaces, or from a list of tensors As. All cases reduce to the latter.

Arguments

  • As::Vector{<:GenericMPSTensor}: vector of site tensors

  • f::Function=rand: initializer function for tensor data

  • eltype::Type{<:Number}=ComplexF64: scalar type of tensors

  • physicalspaces::Vector{<:Union{S, CompositeSpace{S}}: list of physical spaces

  • N::Int: number of sites

  • physicalspace::Union{S,CompositeSpace{S}}: local physical space

  • virtualspaces::Vector{<:Union{S, CompositeSpace{S}}: list of virtual spaces

  • maxvirtualspace::S: maximum virtual space

Keywords

  • normalize: normalize the constructed state
  • overwrite=false: overwrite the given input tensors
  • left=oneunit(S): left-most virtual space
  • right=oneunit(S): right-most virtual space
source
MPSKit.InfiniteMPSType
InfiniteMPS{A<:GenericMPSTensor,B<:MPSBondTensor} <: AbtractMPS

Type that represents an infinite Matrix Product State.

Fields

  • AL – left-gauged MPS tensors
  • AR – right-gauged MPS tensors
  • AC – center-gauged MPS tensors
  • CR – gauge tensors

Notes

By convention, we have that:

  • AL[i] * CR[i] = AC[i] = CR[i-1] * AR[i]
  • AL[i]' * AL[i] = 1
  • AR[i] * AR[i]' = 1

Constructors

InfiniteMPS([f, eltype], physicalspaces::Vector{<:Union{S, CompositeSpace{S}},
            virtualspaces::Vector{<:Union{S, CompositeSpace{S}};
            kwargs...) where {S<:ElementarySpace}
InfiniteMPS(As::AbstractVector{<:GenericMPSTensor}; kwargs...)
InfiniteMPS(ALs::AbstractVector{<:GenericMPSTensor}, C₀::MPSBondTensor;
            kwargs...)

Construct an MPS via a specification of physical and virtual spaces, or from a list of tensors As, or a list of left-gauged tensors ALs.

Arguments

  • As::AbstractVector{<:GenericMPSTensor}: vector of site tensors

  • ALs::AbstractVector{<:GenericMPSTensor}: vector of left-gauged site tensors

  • C₀::MPSBondTensor: initial gauge tensor

  • f::Function=rand: initializer function for tensor data

  • eltype::Type{<:Number}=ComplexF64: scalar type of tensors

  • physicalspaces::AbstractVector{<:Union{S, CompositeSpace{S}}: list of physical spaces

  • virtualspaces::AbstractVector{<:Union{S, CompositeSpace{S}}: list of virtual spaces

Keywords

  • tol: gauge fixing tolerance
  • maxiter: gauge fixing maximum iterations
source
MPSKit.WindowMPSType
WindowMPS{A<:GenericMPSTensor,B<:MPSBondTensor} <: AbstractFiniteMPS

Type that represents a finite Matrix Product State embedded in an infinte Matrix Product State.

Fields

  • left_gs::InfiniteMPS – left infinite environment
  • window::FiniteMPS – finite window Matrix Product State
  • right_gs::InfiniteMPS – right infinite environment

Constructors

WindowMPS(left_gs::InfiniteMPS, window_state::FiniteMPS, [right_gs::InfiniteMPS])
WindowMPS(left_gs::InfiniteMPS, window_tensors::AbstractVector, [right_gs::InfiniteMPS])
WindowMPS([f, eltype], physicalspaces::Vector{<:Union{S, CompositeSpace{S}},
          virtualspaces::Vector{<:Union{S, CompositeSpace{S}}, left_gs::InfiniteMPS,
          [right_gs::InfiniteMPS])
WindowMPS([f, eltype], physicalspaces::Vector{<:Union{S,CompositeSpace{S}}},
          maxvirtualspace::S, left_gs::InfiniteMPS, [right_gs::InfiniteMPS])

Construct a WindowMPS via a specification of left and right infinite environment, and either a window state or a vector of tensors to construct the window. Alternatively, it is possible to supply the same arguments as for the constructor of FiniteMPS, followed by a left (and right) environment to construct the WindowMPS in one step.

Note

By default, the right environment is chosen to be equal to the left, however no copy is made. In this case, changing the left state will also affect the right state.

WindowMPS(state::InfiniteMPS, L::Int)

Construct a WindowMPS from an InfiniteMPS, by promoting a region of length L to a FiniteMPS.

source
MPSKit.MPSMultilineType
const MPSMultiline = Multiline{<:InfiniteMPS}

Type that represents multiple lines of InfiniteMPS objects.

Constructors

MPSMultiline(mpss::AbstractVector{<:InfiniteMPS})
MPSMultiline([f, eltype], physicalspaces::Matrix{<:Union{S, CompositeSpace{S}},
             virtualspaces::Matrix{<:Union{S, CompositeSpace{S}}) where
             {S<:ElementarySpace}
MPSMultiline(As::AbstractMatrix{<:GenericMPSTensor}; kwargs...)
MPSMultiline(ALs::AbstractMatrix{<:GenericMPSTensor}, 
             C₀::AbstractVector{<:MPSBondTensor}; kwargs...)

See also: Multiline

source

Operators

Missing docstring.

Missing docstring for SparseMPO. Check Documenter's build log for details.

MPSKit.MPOHamiltonianType
MPOHamiltonian

represents a general periodic quantum hamiltonian

really just a sparsempo, with some garantuees on its structure
source

Environments

MPSKit.FinEnvType
FinEnv keeps track of the environments for FiniteMPS / WindowMPS
It automatically checks if the queried environment is still correctly cached and if not - recalculates

if above is set to nothing, above === below.

opp can be a vector of nothing, in which case it'll just be the overlap
source
Missing docstring.

Missing docstring for MPSKit.IDMRGEnvs. Check Documenter's build log for details.

Generic actions

MPSKit.∂CFunction
Zero-site derivative (the C matrix to the right of pos)
source
Missing docstring.

Missing docstring for ∂∂C. Check Documenter's build log for details.

Missing docstring.

Missing docstring for ∂AC. Check Documenter's build log for details.

Missing docstring.

Missing docstring for ∂∂AC. Check Documenter's build log for details.

Missing docstring.

Missing docstring for ∂∂AC2. Check Documenter's build log for details.

Missing docstring.

Missing docstring for c_proj. Check Documenter's build log for details.

Missing docstring.

Missing docstring for ac_proj. Check Documenter's build log for details.

Missing docstring.

Missing docstring for ac2_proj. Check Documenter's build log for details.

MPSKit.expectation_valueFunction
expectation_value(ψ, O, [location], [environments])

Compute the expectation value of an operator O on a state ψ. If location is given, the operator is applied at that location. If environments is given, the expectation value might be computed more efficiently by re-using previously calculated environments.

Note

For MPOHamiltonian, the expectation value is not uniquely defined, as it is unclear to what site a given term belongs. For this reason, the returned value is half the expectation value of all terms that start and end on the site.

Arguments

  • ψ::AbstractMPS : the state on which to compute the expectation value
  • O : the operator to compute the expectation value of. This can either be an AbstractMPO, a single AbstractTensorMap or an array of AbstractTensorMaps.
  • location::Union{Int,AbstractRange{Int}} : the location at which to apply the operator. Only applicable for operators that act on a subset of all sites.
  • environments::Cache : the environments to use for the calculation. If not given, they will be calculated.
source

Algorithms

MPSKit.find_groundstateFunction
find_groundstate(ψ, H, [environments]; kwargs...)
find_groundstate(ψ, H, algorithm, environments)

Compute the groundstate for Hamiltonian H with initial guess ψ. If not specified, an optimization algorithm will be attempted based on the supplied keywords.

Arguments

  • ψ::AbstractMPS: initial guess
  • H::AbstractMPO: operator for which to find the groundstate
  • [environments]: MPS environment manager
  • algorithm: optimization algorithm

Keywords

  • tol::Float64: tolerance for convergence criterium
  • maxiter::Int: maximum amount of iterations
  • verbosity::Int: display progress information
source
MPSKit.timestepFunction
timestep(ψ₀, H, t, dt, [alg], [envs]; kwargs...)
timestep!(ψ₀, H, t, dt, [alg], [envs]; kwargs...)

Time-step the state ψ₀ with Hamiltonian H over a given time step dt at time t, solving the Schroedinger equation: $i ∂ψ/∂t = H ψ$.

Arguments

  • ψ₀::AbstractMPS: initial state
  • H::AbstractMPO: operator that generates the time evolution (can be time-dependent).
  • t::Number: starting time of time-step
  • dt::Number: time-step magnitude
  • [alg]: algorithm to use for the time evolution. Defaults to TDVP.
  • [envs]: MPS environment manager
source
Missing docstring.

Missing docstring for dynamicaldmrg. Check Documenter's build log for details.

MPSKit.excitationsFunction
excitations(H, algorithm::QuasiparticleAnsatz, ψ::FiniteQP, [left_environments],
            [right_environments]; num=1)
excitations(H, algorithm::QuasiparticleAnsatz, ψ::InfiniteQP, [left_environments],
            [right_environments]; num=1, solver=Defaults.solver)
excitations(H, algorithm::FiniteExcited, ψs::NTuple{<:Any, <:FiniteMPS};
            num=1, init=copy(first(ψs)))

Compute the first excited states and their energy gap above a groundstate.

Arguments

  • H::AbstractMPO: operator for which to find the excitations
  • algorithm: optimization algorithm
  • ψ::QP: initial quasiparticle guess
  • ψs::NTuple{N, <:FiniteMPS}: N first excited states
  • [left_environments]: left groundstate environment
  • [right_environments]: right groundstate environment

Keywords

  • num::Int: number of excited states to compute
  • solver: algorithm for the linear solver of the quasiparticle environments
  • init: initial excited state guess
source
MPSKit.approximateFunction
approximate(ψ₀, (O, ψ), algorithm, [environments]; kwargs...)
approximate!(ψ₀, (O, ψ), algorithm, [environments]; kwargs...)

Compute an approximation to the application of an operator O to the state ψ in the form of an MPS ψ₀.

Arguments

  • ψ₀::AbstractMPS: initial guess of the approximated state
  • (O::AbstractMPO, ψ::AbstractMPS): operator O and state ψ to be approximated
  • algorithm: approximation algorithm. See below for a list of available algorithms.
  • [environments]: MPS environment manager

Keywords

  • tol::Float64: tolerance for convergence criterium
  • maxiter::Int: maximum amount of iterations
  • verbosity::Int: display progress information

Algorithms

  • DMRG: Alternating least square method for maximizing the fidelity with a single-site scheme.

  • DMRG2: Alternating least square method for maximizing the fidelity with a two-site scheme.

  • IDMRG1: Variant of DMRG for maximizing fidelity density in the thermodynamic limit.

  • IDMRG2: Variant of DMRG2 for maximizing fidelity density in the thermodynamic limit.

  • VOMPS: Tangent space method for truncating uniform MPS.

source

Groundstate algorithms

MPSKit.VUMPSType
VUMPS{F} <: Algorithm

Variational optimization algorithm for uniform matrix product states, as introduced in https://arxiv.org/abs/1701.07035.

Fields

  • tol::Float64: tolerance for convergence criterium

  • maxiter::Int: maximum amount of iterations

  • finalize::F: user-supplied function which is applied after each iteration, with signature finalize(iter, ψ, H, envs) -> ψ, envs

  • verbosity::Int: display progress information

  • alg_gauge=Defaults.alg_gauge(): algorithm for gauging

  • alg_eigsolve=Defaults.alg_eigsolve(): algorithm for eigensolvers

  • alg_environments=Defaults.alg_environments(): algorithm for updating environments

source
MPSKit.IDMRG1Type
IDMRG1{A} <: Algorithm

Single site infinite DMRG algorithm for finding groundstates.

Fields

  • tol::Float64: tolerance for convergence criterium
  • tol_gauge::Float64: tolerance for gauging algorithm
  • eigalg::A: eigensolver algorithm
  • maxiter::Int: maximum number of outer iterations
  • verbosity::Int: display progress information
source
MPSKit.IDMRG2Type
IDMRG2{A} <: Algorithm

2-site infinite DMRG algorithm for finding groundstates.

Fields

  • tol::Float64: tolerance for convergence criterium
  • tol_gauge::Float64: tolerance for gauging algorithm
  • eigalg::A: eigensolver algorithm
  • maxiter::Int: maximum number of outer iterations
  • verbosity::Int: display progress information
  • trscheme::TruncationScheme: truncation algorithm for [tsvd]TensorKit.tsvd
source
MPSKit.DMRGType
DMRG{A,F} <: Algorithm

Single site DMRG algorithm for finding groundstates.

Fields

  • tol::Float64: tolerance for convergence criterium
  • eigalg::A: eigensolver algorithm
  • maxiter::Int: maximum number of outer iterations
  • verbosity::Int: display progress information
  • finalize::F: user-supplied function which is applied after each iteration, with signature finalize(iter, ψ, H, envs) -> ψ, envs
source
MPSKit.DMRG2Type
DMRG2{A,F} <: Algorithm

2-site DMRG algorithm for finding groundstates.

Fields

  • tol::Float64: tolerance for convergence criterium
  • eigalg::A: eigensolver algorithm
  • maxiter::Int: maximum number of outer iterations
  • verbosity::Int: display progress information
  • finalize::F: user-supplied function which is applied after each iteration, with signature finalize(iter, ψ, H, envs) -> ψ, envs
  • trscheme: truncation algorithm for [tsvd]TensorKit.tsvd
source
MPSKit.GradientGrassmannType
GradientGrassmann <: Algorithm

Variational gradient-based optimization algorithm that keeps the MPS in left-canonical form, as points on a Grassmann manifold. The optimization is then a Riemannian gradient descent with a preconditioner to induce the metric from the Hilbert space inner product.

Fields

  • method::OptimKit.OptimizationAlgorithm: algorithm to perform the gradient search
  • finalize!::Function: user-supplied function which is applied after each iteration, with signature finalize!(x::GrassmannMPS.ManifoldPoint, f, g, numiter) -> x, f, g

Constructors

GradientGrassmann(; kwargs...)

Keywords

  • method=ConjugateGradient: instance of optimization algorithm, or type of optimization algorithm to construct
  • finalize!: finalizer algorithm
  • tol::Float64: tolerance for convergence criterium
  • maxiter::Int: maximum amount of iterations
  • verbosity::Int: level of information display
source

Time evolution algorithms

MPSKit.TDVPType
TDVP{A} <: Algorithm

Single site TDVP algorithm for time evolution.

Fields

  • integrator::A: integration algorithm (defaults to Lanczos exponentiation)
  • tolgauge::Float64: tolerance for gauging algorithm
  • gaugemaxiter::Int: maximum amount of gauging iterations
  • finalize::F: user-supplied function which is applied after each timestep, with signature finalize(t, Ψ, H, envs) -> Ψ, envs
source
MPSKit.TDVP2Type
TDVP2{A} <: Algorithm

2-site TDVP algorithm for time evolution.

Fields

  • integrator::A: integrator algorithm (defaults to Lanczos exponentiation)
  • tolgauge::Float64: tolerance for gauging algorithm
  • gaugemaxiter::Int: maximum amount of gauging iterations
  • trscheme: truncation algorithm for [tsvd]TensorKit.tsvd
  • finalize::F: user-supplied function which is applied after each timestep, with signature finalize(t, Ψ, H, envs) -> Ψ, envs
source
Missing docstring.

Missing docstring for TaylorCluster. Check Documenter's build log for details.

Missing docstring.

Missing docstring for WII. Check Documenter's build log for details.

Leading boundary algorithms

Missing docstring.

Missing docstring for VUMPS. Check Documenter's build log for details.

MPSKit.VOMPSType
VOMPS{F} <: Algorithm

Power method algorithm for infinite MPS. SciPost:4.1.004

Fields

  • tol::Float64: tolerance for convergence criterium

  • maxiter::Int: maximum amount of iterations

  • finalize::F: user-supplied function which is applied after each iteration, with signature finalize(iter, ψ, toapprox, envs) -> ψ, envs

  • verbosity::Int: display progress information

  • alg_gauge=Defaults.alg_gauge(): algorithm for gauging

  • alg_environments=Defaults.alg_environments(): algorithm for updating environments

source
Missing docstring.

Missing docstring for GradientGrassmann. Check Documenter's build log for details.

Bond change algorithms

MPSKit.OptimalExpandType
struct OptimalExpand <: Algorithm end

An algorithm that expands the given mps using the algorithm given in the VUMPS paper, by selecting the dominant contributions of a two-site updated MPS tensor, orthogonal to the original ψ.

Fields

  • trscheme::TruncationScheme = truncdim(1) : The truncation scheme to use.
source
MPSKit.RandExpandType
struct RandExpand <: Algorithm end

An algorithm that expands the bond dimension by adding random unitary vectors that are orthogonal to the existing state. This is achieved by performing a truncated SVD on a random two-site MPS tensor, which is made orthogonal to the existing state.

Fields

  • trscheme::TruncationScheme = truncdim(1) : The truncation scheme to use.
source
MPSKit.VUMPSSvdCutType
struct VUMPSSvdCut <: Algorithm end

An algorithm that uses an IDMRG2 step to change the bond dimension of a state.

Fields

  • tol_gauge::Real = Defaults.tolgauge : The tolerance for the gauge.
  • tol::Real = Defaults.tol : The tolerance for the Galerkin truncation.
  • tol_eigenval::Real = Defaults.tol : The tolerance for the eigenvalue solver.
  • trscheme::TruncationScheme = notrunc() : The truncation scheme to use.
source
MPSKit.SvdCutType
struct SvdCut <: Algorithm end

An algorithm that uses truncated SVD to change the bond dimension of a ψ.

Fields

  • trscheme::TruncationScheme = notrunc() : The truncation scheme to use.
source

[Excitations]

MPSKit.QuasiparticleAnsatzType
QuasiparticleAnsatz <: Algorithm

Optimization algorithm for quasiparticle excitations on top of MPS groundstates, as introduced in this paper.

Fields

  • alg::A = Defaults.eigsolver: algorithm to use for the eigenvalue problem

Constructors

QuasiparticleAnsatz()
QuasiparticleAnsatz(; kwargs...)
QuasiparticleAnsatz(alg)

Create a QuasiparticleAnsatz algorithm with the given algorithm, or by passing the keyword arguments to Arnoldi.

source
MPSKit.FiniteExcitedType
FiniteExcited{A} <: Algorithm

Variational optimization algorithm for excitations of finite Matrix Product States by minimizing the energy of $H - λᵢ |ψᵢ><ψᵢ|$.

Fields

  • gsalg::A: optimization algorithm.
  • weight::Float64: energy penalty for previous states.
source