FiniteMPS.MPSDefaultModule

module MPSDefault D::Int64 = 512 tol::Float64 = 1e-8 end

Default maximum bond dimension D and truncation tolerance tol for MPS, used as default truncation scheme in some 2-site update algorithms.

FiniteMPS.SU₂SpinModule
 module SU₂Spin

Prepare the local space of SU₂ spin-1/2.

Fields

 pspace::VectorSpace

Local d = 2 Hilbert space.

 SS::NTuple{2, TensorMap}

Two rank-3 operators of Heisenberg S⋅S interaction.

FiniteMPS.U₁SU₂FermionModule
 module U₁SU₂Fermion

Prepare some commonly used objects for U₁×SU₂ fermions.

Nothing is exported, please call U₁SU₂Fermion.xxx to use xxx.

Fields

 pspace::VectorSpace

Local d = 4 Hilbert space.

 Z::TensorMap

Rank-2 fermion parity operator Z = (-1)^n.

 n::TensorMap

Rank-2 particle number operator n = n↑ + n↓.

 nd::TensorMap

Rank-2 double occupancy operator nd = n↑n↓.

 SS::NTuple{2, TensorMap}

Two rank-3 operators of Heisenberg S⋅S interaction.

 FdagF::NTuple{2, TensorMap}

Two rank-3 operators of hopping c↑^dag c↑ + c↓^dag c↓.

 FFdag::NTuple{2, TensorMap}

Two rank-3 operators of hopping c↑ c↑^dag + c↓ c↓^dag.

 ΔₛdagΔₛ::NTuple{4, TensorMap}

Four operators of singlet pairing correlation Δₛ^dagΔₛ, where Δₛ = (c↓c↑ - c↑c↓)/√2. Rank = (3, 4, 4, 3).

 ΔₜdagΔₜ::NTuple{4, TensorMap}

Four operators of triplet pairing correlation Δₜ^dag⋅Δₜ, where Δₜ is the triplet pairing operator that carries 2 charge and 1 spin quantum numbers. Rank = (3, 4, 4, 3).

 Δₛ::NTuple{2, TensorMap}
 Δₛdag::NTuple{2, TensorMap}

Singlet pairing operators Δₛ and Δₛ^dag. Rank = (4, 3). Note the first operator has nontrivial left bond index.

 CpCm::NTuple{2, TensorMap}

Two rank-3 operators of C+C- correlation where C+ = c↑^dag c↓^dag and C- = C+^dag = c↓c↑. Note both operators are bosonic.

FiniteMPS.U₁SU₂tJFermionModule
 module U₁SU₂tJFermion

Prepare some commonly used objects for U₁×SU₂ tJ fermions, i.e. local d = 3 Hilbert space without double occupancy.

Behaviors of all operators are the same as U₁SU₂Fermion up to the projection, details please see U₁SU₂Fermion.

FiniteMPS.U₁SpinlessFermionModule
 module U₁SpinlessFermion

Prepare some commonly used objects for U₁ spinless fermions.

Fields

 pspace::VectorSpace

Local d = 2 Hilbert space.

 Z::TensorMap

Rank-2 fermion parity operator Z = (-1)^n.

 n::TensorMap

Rank-2 particle number operator.

 FdagF::NTuple{2, TensorMap}

Two rank-3 operators of hopping c^dag c.

 FFdag::NTuple{2, TensorMap}

Two rank-3 operators of hopping c c^dag.

FiniteMPS.U₁U₁FermionModule
 module U₁U₁Fermion

Prepare the local space of d = 4 spin-1/2 fermions with U₁ charge and U₁ spin symmetry.

Fields

 pspace::VectorSpace

Local d = 4 Hilbert space.

 Z::TensorMap

Rank-2 fermion parity operator Z = (-1)^n.

 n₊::TensorMap
 n₋::TensorMap
 n::TensorMap

Rank-2 particle number operators. and denote spin up and down as and cannot be used in variable names.

 nd::TensorMap

Rank-2 double occupancy operator nd = n↑n↓.

 Sz::TensorMap

Rank-2 spin-z operator Sz = (n↑ - n↓)/2.

 S₊₋::NTuple{2, TensorMap}
 S₋₊::NTuple{2, TensorMap}

Two rank-3 operators of S₊₋ and S₋₊ interaction. Note Heisenberg S⋅S = SzSz + (S₊₋ + S₋₊)/2.

 FdagF₊::NTuple{2, TensorMap}
 FdagF₋::NTuple{2, TensorMap}

Two rank-3 operators of hopping c↑^dag c↑ and c↓^dag c↓.

 FFdag₊::NTuple{2, TensorMap}
 FFdag₋::NTuple{2, TensorMap}

Two rank-3 operators of hopping c↑ c↑^dag and c↓ c↓^dag.

 ΔdagΔ₊₊::NTuple{4, TensorMap}
 ΔdagΔ₋₋::NTuple{4, TensorMap}
 ΔdagΔ₊₋::NTuple{4, TensorMap}

Rank-4 operators of pairing correlation. Note ΔdagΔ₊₋ means c↑^dag c↓^dag c↓ c↑.

FiniteMPS.ℤ₂SU₂FermionModule
 module ℤ₂SU₂Fermion

Prepare some commonly used objects for ℤ₂×SU₂ fermions. Basis convention in (0, 0) sector is {|0⟩, |↑↓⟩}.

Each operator has the same name and behavior as U₁SU₂Fermion, details please see U₁SU₂Fermion.

FiniteMPS.SimpleLeftTensorType
 const SimpleLeftTensor = Union{Nothing, LocalLeftTensor}

Type of left environment tensor of simple MPO, i.e. a channel of sparse MPO.

FiniteMPS.SimpleRightTensorType
 const SimpleRightTensor = Union{Nothing, LocalRightTensor}

Type of right environment tensor of simple MPO, i.e. a channel of sparse MPO.

FiniteMPS.AbstractDistributionType
 abstract type AbstractDistribution{F<:Union{Float64, ComplexF64}}

Abstract type of probability distributions on Stiefel manifold Vₙₖ.

FiniteMPS.AbstractLocalOperatorType
 abstract type AbstractLocalOperator <: AbstractTensorWrapper

Wrapper type for classifying differnet local operators in order to accelerate contractions.

FiniteMPS.AbstractMPSTensorType
abstract type AbstractMPSTensor <: AbstractTensorWrapper

Elements of MPS, note a MPO is nothing but a MPS with rank-4 tensors, hence we using this type to deal with both MPS and MPO

FiniteMPS.AbstractStoreTypeType
 abstract type AbstractStoreType

Only has 2 concrete types StoreMemory and StoreDisk, to determine how the collections such as MPS, MPO and Environment store the local tensors.

FiniteMPS.AbstractTensorWrapperType
 abstract type AbstractTensorWrapper

Wrapper type for classifying different Tensors.

Note each concrete subtype must have a field A::AbstractTensorMap to save the Tensor.

FiniteMPS.AdjointMPSType
 struct AdjointMPS{L} <: AbstractMPS{L}
      parent::DenseMPS{L}
 end

Lazy wrapper type for adjoint of MPS.

 adjoint(::DenseMPS) -> ::AdjointMPS
 adjoint(::AdjointMPS) -> ::DenseMPS

Functions to be directly propagated to the parent:

 lastindex, length, keys, norm, normalize!, Center, iterate, canonicalize!

Functions to be propagated to the parent with some adaptations:

 getindex, setindex!, coef
FiniteMPS.AdjointMPSTensorType
 struct AdjointMPSTensor{R} <: AbstractMPSTensor
      A::AbstractTensorMap
 end

Lazy wrapper type for tensors of adjoint MPS.

 adjoint(::MPSTensor) -> ::AdjointMPSTensor
 adjoint(::AdjointMPSTensor) -> ::MPSTensor

Convention (' marks codomain):

                   3              4                   R
                   |              |                   |
 2-- A --1'     2--A--1'       3--A--2'       (R-1)-- A  --(R-2)'        
                                  |                 / | \
                                  1'              1' ... (R-3)'
FiniteMPS.BondInfoType
 struct BondInfo
      D::Int64
      DD::Int64
      TrunErr::Float64
      SE::Float64
 end

Type for storing the information of a bond.

Constructors

 BondInfo(s::AbstractTensorMap, ϵ::Float64 = 0.0)

Outer constructor via giving the s tensor and ϵ form tsvd.

 BondInfo(A::AbstractTensorMap, direction::Symbol)
 BondInfo(A::MPSTensor, direction::Symbol)

Outer constructor via giving a tensor A and direction = :L or :R. We cannot get truncation error and singular values hence TrunErr and SE are set to 0.0 and NaN, respectively.

FiniteMPS.CBEAlgorithmType
 abstract type CBEAlgorithm{T <: SweepDirection}

Abstract type of all (controlled bond expansion) CBE algorithms.

FiniteMPS.CBEInfoType
 struct CBEInfo{N}
      Alg::CBEAlgorithm
      info::NTuple{N, BondInfo}
      D₀::NTuple{2, Int64}
      D::NTuple{2, Int64}
      ϵ::Float64
 end

Information of CBE. Alg is the algorithm used. info contain the truncation info of N times svd. D₀ and D are the initial and final bond dimension, respectively. ϵ = |Al*Ar - Al_ex*Ar_ex| if calculated (otherwise NaN).

FiniteMPS.CheapCBEType
 struct CheapCBE{T<:SweepDirection} <: CBEAlgorithm{T}
      D::Int64
      tol::Float64
      check::Bool
 end

Similar to StandardCBE, but not weigthed by the opposite environment to save computational cost. Note this trick will affect the accuracy of the expanded space. It should only be used to redistribute bond dimensions to different sectors, instead of increasing bond dimension adaptively.

FiniteMPS.CompositeMPSTensorType
 struct CompositeMPSTensor{N, T <: NTuple{N, MPSTensor}} <: AbstractMPSTensor
      A::AbstractTensorMap
 end

Wrapper type for multi-site local tensors.

The 2nd parameter indicates the types of the N original on-site tensors.

Constructors

 CompositeMPSTensor{N, T}(::AbstractTensorMap) where T <: NTuple{N, MPSTensor}

Directly construct.

 CompositeMPSTensor(::NTuple{N, MPSTensor})
 CompositeMPSTensor(::MPSTensor, ::MPSTensor, ...)

Contract N on-site tensors to get the N-site tensor.

FiniteMPS.DMRGInfoType
 struct DMRGInfo
      Eg::Float64
      Lanczos::LanczosInfo
      Bond::BondInfo
 end

Information of each DMRG update.

FiniteMPS.DenseMPSType
 abstract type DenseMPS{L, T <:Union{Float64, ComplexF64}} <: AbstractMPS{L}

Abstract type of dense MPS/MPO with length L.

FiniteMPS.FullCBEType
 struct FullCBE{T <: SweepDirection} <: CBEAlgorithm{T} 
      check::Bool
 end

Special case of CBE algorithm, directly keep the full bond space, usually used near boundary.

Constructor

 FullCBE(direction::SweepDirection = AnyDirection(); check::Bool = false)
FiniteMPS.GaussianDistributionType
 struct GaussianDistribution{F} <: AbstractDistribution{F}
      σ::Float64
 end

Type of the gaussian distribution around [Iₖₖ, 0] on Stiefel manifold Vₙₖ.

Note this is not strictly gaussian distribution on Vₙₖ, in practical we will first generate a gaussian distribution N(0,σ) on so(n) and then exp to SO(n).

FiniteMPS.IdentityOperatorType
 mutable struct IdentityOperator <: AbstractLocalOperator
      pspace::Union{Nothing, VectorSpace}
      si::Int64
      strength::Number 
 end

Lazy type of identity operator, used for skipping some tensor contractions.

Constructors

 IdentityOperator([pspace::VectorSpace,] si::Int64, strength::Number = NaN)
FiniteMPS.IdentityProjectiveHamiltonianType
 struct IdentityProjectiveHamiltonian{N} <: AbstractProjectiveHamiltonian
      El::SimpleLeftTensor
      Er::SimpleRightTensor
      si::Vector{Int64}
 end

Special type to deal with the cases which satisfy ⟨Ψ₁|Id|Ψ₂⟩ == ⟨Ψ₁|Ψ₂⟩, thus the environment is a 2-layer simple one.

FiniteMPS.InteractionTreeType
 struct InteractionTree{N, T}
      Root::InteractionTreeNode{Nothing}
 end

Wrapper type of the root node in interaction tree.

N is the children number of root node. When generating sparse MPO, N is exactly the length of boundary left environment tensor. Usually N == 1, except for some cases when we want to deal with several sparse MPOs independently, e.g. H and N in tanTRG when using fixed particle number technique.

T is the type of value of all non-root nodes. Note value of Root is always of type Nothing.

FiniteMPS.InteractionTreeNodeType
 mutable struct InteractionTreeNode{T} 
      Op::Union{Nothing, AbstractLocalOperator}
      value::Union{Nothing, T}
      parent::Union{Nothing,InteractionTreeNode}
      children::Vector{InteractionTreeNode}
 end

Type of node of interaction tree, with field Op to store the local operator, value to store anything others.

Constructors

 InteractionTreeNode(Op::Union{Nothing,AbstractLocalOperator},
      value::T,
      parent::Union{Nothing,InteractionTreeNode} = nothing,
      children::Vector{InteractionTreeNode} = []) -> ::InteractionTreeNode{T}
FiniteMPS.LanczosInfoType
 struct LanczosInfo
      converged::Bool
      normres::Vector{Float64}
      numiter::Int64
      numops::Int64
 end

Similar to KrylovKit.ConvergenceInfo but delete residuals to save memory.

FiniteMPS.LeftPreFuseTensorType
 struct LeftPreFuseTensor{R} <: AbstractEnvironmentTensor
      A::AbstractTensorMap
 end

Wrapper type for the prefused left environment tensors, i.e., the combination of a left environment tensor and a local operator.

Convention (' marks codomain):

  --R-1  R
 |       |
 El------Hl--(R-2)' 
 |       |
  --1'   2'
FiniteMPS.LocalLeftTensorType
 struct LocalLeftTensor{R} <: AbstractEnvironmentTensor
      A::AbstractTensorMap
      tag::NTuple{R, String}
 end

Wrapper type for rank-R left environment tenosr, with an additional field tag to distinguish legs of different channels.

Convention (' marks codomain):

  --R
 |
 | --R-1(') 
 El  ⋮
 | --2(')
 |
  --1'

Constructors

 LocalLeftTensor(A::AbstractTensorMap [, tag::NTuple{R, String}])

Default tag = ("", "", ..., "").

 LocalLeftTensor{R}(A::AbstractTensorMap)

Used for automatic converting, only support default tag.

FiniteMPS.LocalOperatorType
 mutable struct LocalOperator{R₁,R₂} <: AbstractLocalOperator
      A::AbstractTensorMap
      name::String
      si::Int64
      strength::Number
      tag::tag2Tuple{R₁,R₂}
 end

Warpper type for local operators, the building blocks of sparse MPO.

R₁ and R₂ indicate the rank corresponding to codomain and domain, respectively.

Warning: this warpper type does not support automatic converting.

Convention (' marks codomain):

 2      2          3         3
 |      |          |         |
 A      A--3   1'--A     1'--A--4               
 |      |          |         |
 1'     1'         2'        2'

Constructors

 LocalOperator(O::AbstractTensorMap,
      name::Union{String,Symbol},
      si::Int64,
      [,strength::Number = NaN]
      [, tag::tag2Tuple{R₁,R₂}];
      swap::Bool=false)

Default tag: "phys" for physical indices and name for virtual indices.

If swap == ture, it will swap the left and right virtual indices.

FiniteMPS.LocalRightTensorType
 struct LocalRightTensor{R} <: AbstractEnvironmentTensor
      A::AbstractTensorMap
      tag::NTuple{R, String}
 end

Wrapper type for rank-R right environment tenosr, with an additional field tag to distinguish legs of different channels.

Convention (' marks codomain):

     1'--
         |
   2(')--| 
   ⋮     Er  
 R-1(')--| 
         |
      R--

Constructors

 LocalRigthTensor(A::AbstractTensorMap [, tag::NTuple{R, String}])

Default tag = ("", "", ..., "").

 LocalRightTensor{R}(A::AbstractTensorMap)

Used for automatic converting, only support default tag.

FiniteMPS.MPOType
 mutable struct MPO{L, T <:Union{Float64, ComplexF64}, C} <: DenseMPS{L, T}
      const A::AbstractVector{AbstractMPSTensor}
      const Center::Vector{Int64} 
      c::T 
 end

All the fields and constructors are exactly the same to those of MPS, we redefine the type MPO just for using multiple dispacth when implementing the algebra between MPS and MPO.

Details of constructors please see MPS.

FiniteMPS.MPSType
 mutable struct MPS{L, T <:Union{Float64, ComplexF64}, C} <: DenseMPS{L, T}
      const A::AbstractVector{AbstractMPSTensor}
      const Center::Vector{Int64} 
      c::T 
 end

Concrete type of MPS, where L is the length, T == Float64 or ComplexF64 is the number type of local tensors.

C <: AbstractStoreType is the type to determine the storage of local tensors, usually C == StoreMemory. We can use package SerializedElementArrays to store local tensors in disk, in which case C == StoreDisk.

Fields

 const A::AbstractVector{AbstractMPSTensor}

Length L vector to store the local tensors. Note the vector A is immutable while the local tensors in it are mutable.

 const Center::Vector{Int64}

Length 2 vector to label the canonical form. [a, b] means left-canonical from 1 to a-1 and right-canonical from b+1 to L.

 c::T

The global coefficient, i.e. we represented a MPS with L local tensors and an additional scalar c, in order to avoid too large/small local tensors.

Constructors

 MPS{L, T}(A::AbstractVector{<:AbstractMPSTensor}, Center::Vector{Int64} = [1, L], c::T = one(T))

Standard constructor. Note we will promote all local tensors if T == ComplexF64 while the input local tensors in A are of Float64.

 MPS(A::AbstractVector{<:AbstractMPSTensor}, args...)

Deduce type parameters via L = length(A) and T = typeof(c). For default c cases, assume T = Float64 if all local tensors are real, otherwise T = ComplexF64.

 MPS(A::AbstractVector{<:AbstractTensorMap}, args...)

Automatically convert local tensors to warpper type MPSTensor.

 MPS{L, T}() 
 MPS(L, T = Float64)

Initialize an MPS{L, T} with the local tensors A to be filled. Note we initialize Center = [1, L] and c = one(T).

Newest update, kwargs disk::Bool = false can be added to each constructor to control how to store the local tensors.

FiniteMPS.MPSTensorType
 struct MPSTensor{R} <: AbstractMPSTensor
      A::AbstractTensorMap
 end

Wrapper type for rank-R MPS local tensors.

Convention (' marks codomain):

      3 ... (R-1)
      \ | /  
 1'--   A  ---R         1'-- A -- 2
        | 
        2'

In particular, R == 2 for bond tensor.

Constructors

 MPSTensor(::AbstractTensorMap) 
 MPSTensor{R}(::AbstractTensorMap)
FiniteMPS.NoCBEType
 struct NoCBE{T <: SweepDirection} <: CBEAlgorithm{T}
FiniteMPS.ObservableTreeType
 struct ObservableTree{N}
      Root::InteractionTreeNode{Tuple{String, Vararg{Int64}}}
 end

Similar to InteractionTree but specially used for calculation of observables. The value field is used to tell which observable corresponds to the current node.

Constructors

 ObservableTree{N}()

Initialize an empty object.

FiniteMPS.PreFuseProjectiveHamiltonianType
 struct PreFuseProjectiveHamiltonian{N, Tl, Tr} <: AbstractProjectiveHamiltonian 
      El::Tl
      Er::Tr
      si::Vector{Int64}
      E₀::Float64
 end

Prefused N-site projective Hamiltonian. Note El and Er can be a original environment tensor or a prefused one, depending on N. If N == 1, only one of them will be prefused.

FiniteMPS.SimpleEnvironmentType
 struct SimpleEnvironment{L, N, T<:NTuple{N, AbstractMPS{L}}, C<:AbstractStoreType} <: AbstractEnvironment{L}
      layer::T
      El::AbstractVector{SimpleLeftTensor}
      Er::AbstractVector{SimpleRightTensor}
      Center::Vector{Int64} 
 end

Type of all environments whose local left/right tensors are simple, i.e. a single tensor instead of a vector storing several tensors.

Center = [a, b] means El[1:a] and Er[b:end] are valid.

FiniteMPS.SparseEnvironmentType
 struct SparseEnvironment{L, N, T<:NTuple{N, AbstractMPS{L}}, C<:AbstractStoreType} <: AbstractEnvironment{L}
      layer::T
      El::AbstractVector{SparseLeftTensor}
      Er::AbstractVector{SparseRightTensor}
      Center::Vector{Int64} 
 end

Type of all environments whose local left/right tensors are vectors storing several tensors, usually due to a SparseMPO in a layer.

Center = [a, b] means El[1:a] and Er[b:end] are valid.

FiniteMPS.SparseLeftTensorType
 const SparseLeftTensor = Vector{SimpleLeftTensor}

Type of left environment tensor of sparse MPO.

FiniteMPS.SparseMPOType
 struct SparseMPO{L} <: AbstractMPS{L}
      A::Vector{SparseMPOTensor}
 end

Concrete type of sparse MPO.

Note an instance of this type is usually a Hamiltonian, which will not cost too much memory, therefore we always store the local tensors in memory.

Constructors

 SparseMPO(A::AbstractVector{SparseMPOTensor})
FiniteMPS.SparseProjectiveHamiltonianType
 struct SparseProjectiveHamiltonian{N} <: AbstractProjectiveHamiltonian  
      El::SparseLeftTensor
      Er::SparseRightTensor
      H::NTuple{N, SparseMPOTensor}
      si::Vector{Int64}
      validIdx::Vector{Tuple}
      E₀::Float64
 end

N-site projective Hamiltonian, sparse version. Note we shift H to H - E₀ to avoid numerical overflow.

Convention: – – – – | | | | | | | El– i – H1 – j –Er El– i – H1 – j – H2 – k –Er ... | | | | | | | – – – –

validIdx stores all tuples (i, j, ...) which are valid, i.e. all El[i], H1[i, j] and Er[j] are not nothing (N == 1).

FiniteMPS.SparseRightTensorType
 const SparseRightTensor = Vector{SimpleRightTensor}

Type of right environment tensor of sparse MPO.

FiniteMPS.StandardCBEType
 struct StandardCBE{T <: SweepDirection} <: CBEAlgorithm{T}
      D::Int64
      tol::Float64
      check::Bool
 end

Standard CBE algorithm, details see PhysRevLett.130.246402. Note some detailed implementations are modified so that it can be compatible with the sparse MPO.

Constructor

 StandardCBE(D::Int64, tol::Float64, direction::SweepDirection = AnyDirection(); check::Bool = false)
FiniteMPS.StoreDiskType
 struct StoreDisk <: AbstractStoreType

Tell the collection to store local tensors in disk.

FiniteMPS.StoreMemoryType
 struct StoreMemory <: AbstractStoreType

Tell the collection to store local tensors in memory.

FiniteMPS.SweepDirectionType
 abstract type SweepDirection end

Direction of sweeps, wrapped as a type to be passed as a type parameter.

FiniteMPS.TDVPInfoType
 struct TDVPInfo{N,T}
      dt::T
      Lanczos::LanczosInfo
      Bond::BondInfo
 end

Information of each N-site TDVP update.

Constructors

 TDVPInfo{N}(dt::Number, Lanczos::LanczosInfo, Bond::BondInfo)
FiniteMPS.TDVPIntegratorType
 struct TDVPIntegrator{N} 
      dt::NTuple{N, Float64}
      direction::NTuple{N, SweepDirection}
 end

Type of TDVP integrators using composition methods. Note dt is relative, thus sum(dt) == 1.

Constructors

 TDVPIntegrator(dt::NTuple{N, Rational}, direction::NTuple{N, SweepDirection})

Standard constructor with checking dt and direction.

 TDVPIntegrator(dt::Number...)

Assume direction is repeated as L-R-L-R-..., thus length(dt) must be even.

FiniteMPS.UniformDistributionType
 struct UniformDistribution{F} <: AbstractDistribution{F}

Type of the uniform distribution on Stiefel manifold Vₙₖ, induced by the Haar measure on O(n) == Vₙₙ.

FiniteMPS.tag2TupleType
 const tag2Tuple{R₁,R₂} = Tuple{NTuple{R₁,String}, NTuple{R₂,String}}

Type of field tag of LocalOperator.

Base.:*Method
 *(A::LocalOperator{R₁, R₂}, B::LocalOperator{R₃,R₄}) -> ::LocalOperator{R₁ + R₃ - 1, R₂ + R₄ - 1}

The multiplication of two local operators.

Since we only use this function when generating InteractionTree, field strength of A and B must be NaN.

Field name of output obj is "A.name" * "B.name".

Warning: we write this function case by case via multiple dispatch, hence it may throw a "no method matching" error for some interactions.

Base.:*Method
 *(A::LocalLeftTensor{R}, B::LocalRightTensor{R}) -> ::Number

Contract a left environment tenosr and a right environment tensor with the same virtual spaces to get a scalar.

Base.:*Method
 *(El::LocalLeftTensor, A::MPSTensor)

Conctract a MPS tensor to a left environment tensor.

Base.:*Method
 *(A::MPSTensor, Er::LocalRightTensor)

Conctract a MPS tensor to a right environment tensor.

Base.:*Method
 *(A::MPSTensor, B::MPSTensor) -> ::AbstractTensorMap

Contract the virtual bond between 2 neighbor local tensors.

Base.:*Method
 *(O::LocalOperator{R₁,R₂}, A::MPSTensor{R₃}) -> ::MPSTensor{R₁ + R₂ +R₃ - 2}

Apply a local operator on a MPS tensor.

Convention: R₁+2 R₁+3 ... | / 1–-A– end | 2–-O –- ... / | ... R₁+1 i.e., merge the domain and codomain, legs of A first.

Base.:+Method
 +(A::LocalOperator{R₁,R₂}, B::LocalOperator{R₁,R₂}) -> ::LocalOperator{R₁,R₂}

Plus of two local operators on the same site. Note the strength of each one must not be NaN.

Field name of output obj is "A.name(A.strength) + B.name(B.strength)".

Base.:==Method
 ==(A::AbstractLocalOperator, B::AbstractLocalOperator) -> ::Bool

Test if two LocalOperator objects are equal. Note we do not consider the field strength.

Base.complexMethod
 complex(obj::DenseMPS{L}) -> ::DenseMPS{L, ComplexF64}

Return a copy of given MPS but with ComplexF64 as basic field.

FiniteMPS.AutomataMPOFunction
 AutomataMPO(Tree::InteractionTree,
      L::Int64 = treeheight(Tree.Root) - 1)
      -> ::SparseMPO

Convert an interaction tree to a sparse MPO.

FiniteMPS.CBEMethod
 CBE(PH::SparseProjectiveHamiltonian{2},
      Al::MPSTensor, Ar::MPSTensor,
      Alg::CBEAlgorithm;
      kwargs...) -> Al_ex::MPSTensor, Ar_ex::MPSTensor, info::CBEInfo

Return the two local tensors Al_ex and Ar_ex after CBE.

FiniteMPS.CenterMethod
 Center(obj::AbstractEnvironment) -> Vector (length 2)

Interface of environment, return the info of canonical center. Center = [a, b] means El[1:a] and Er[b:end] are valid.

FiniteMPS.CenterMethod
 Center(obj::DenseMPS) -> Vector (length 2)

Interface of DenseMPS, return the info of canonical center. [a, b] means left-canonical from 1 to a-1 and right-canonical from b+1 to L.

FiniteMPS.DMRGSweep1!Method
 DMRGSweep1!(Env::SparseEnvironment{L,3,T}, ::SweepDirection; kwargs...)
      -> info::Vector{DMRGInfo}, Timer::TimerOutput

1-site DMRG sweep from left to right or sweep back from right to left.

Kwargs

 krylovalg::KrylovKit.KrylovAlgorithm = DMRGDefaultLanczos

Krylov algorithm used in DMRG update.

 GCstep::Bool = false

GC.gc() manually after each step if true.

 GCsweep::Bool = false

GC.gc() manually after each (left to right or right to left) sweep if true.

 verbose::Int64 = 0

Print the TimerOutput after each sweep or each local update if verbose = 1 or 2, respectively.

 CBEAlg::CBEAlgorithm = NoCBE()

CBE algorithm for 1-DMRG.

 trunc::TruncationScheme = notrunc()

Control the truncation after each update, only used together with CBE. Details see tsvd.

FiniteMPS.DMRGSweep2!Method
 DMRGSweep2!(Env::SparseEnvironment{L,3,T}, ::SweepDirection; kwargs...) 
      -> info::Vector{DMRGInfo}, Timer::TimerOutput

2-site DMRG sweep from left to right or sweep back from right to left.

Kwargs

 krylovalg::KrylovKit.KrylovAlgorithm = DMRGDefaultLanczos

Krylov algorithm used in DMRG update.

 trunc::TruncationScheme = truncbelow(MPSDefault.tol) & truncdim(MPSDefault.D)

Control the truncation in svd after each 2-site update. Details see tsvd.

 GCstep::Bool = false

GC.gc() manually after each step if true.

 GCsweep::Bool = false

GC.gc() manually after each (left to right or right to left) sweep if true.

 verbose::Int64 = 0

Print the TimerOutput after each sweep or each local update if verbose = 1 or 2, respectively.

 noise::Real = 0

Add noise to the 2-site local tensor after each update.

FiniteMPS.EnvironmentMethod
 Environment(M::AbstractMPS{L}...; kwargs...)

Generic constructor of environments.

Generate an environment object according to the input MPS/MPO objects. The boundary environment tensor, i.e. El[1] and Er[L] will be initialized.

Kwargs

 disk::Bool = false

Store the local environment tensor in disk(true) or in memory(false).

 El::Union{SimpleLeftTensor, SparseLeftTensor}
 Er::Union{SimpleLeftTensor, SparseLeftTensor}

Initialize boundary El or Er manually. Default value is generated by function _defaultEl or _defaultEr, respectively.

FiniteMPS.ProjHamMethod
 ProjHam(Env::SparseEnvironment, siL::Int64 [, siR::Int64 = siL]; E₀::Number = 0.0)

Generic constructor for N-site projective Hamiltonian, where N = siR - siL + 1.

 ProjHam(Env::SimpleEnvironment, siL::Int64 [, siR::Int64 = siL])

Construct the special IdentityProjectiveHamiltonian from a simple environment.

FiniteMPS.SETTNMethod
 SETTN(H::SparseMPO, β::Number; kwargs...) -> ρ::MPO, lsF::Vector{Float64}

Use series-expansion thermal tensor network (SETTN)[https://doi.org/10.1103/PhysRevB.95.161104] method to initialize a high-temperature MPO ρ = e^(-βH/2). Note ρ is unnormalized. The list of free energy F = -lnTr[ρρ^†]/β with different expansion orders is also returned.

Kwargs

 trunc::TruncationScheme = truncbelow(D) (this keyword argument is necessary!) 
 disk::Bool = false
 maxorder::Int64 = 4
 tol::Float64 = 1e-8
 bspace::VectorSpace (details please see identityMPO)
 compress::Float64 = 1e-16 (finally compress ρ with `tol = compress`)
 ρ₀::MPO (initial MPO, default = Id)

Note we use mul! and axpby! to implement H^n -> H^(n+1) and ρ -> ρ + (-βH/2)^n / n!, respectively. All kwargs of these two functions are valid and will be propagated to them appropriately. We find that it is unstable if truncating the bond dimension with tol, thus the input keyword argument trunc must be a TruncationDimension object (NoTruncation is also unallowed to avoid infinitely growing bond dimension).

FiniteMPS.SymmetricIntegratorMethod
 SymmetricIntegrator(p::Int64) -> TDVPIntegrator

Construct predefined symmetric integrators with p-th order. Only p = 2, 3, 4 are supported.

FiniteMPS.TDVPSweep1!Method
 TDVPSweep1!(Env::SparseEnvironment{L,3,T},
      dt::Number,
      direction::SweepDirection;
      kwargs...) -> info, TimerSweep

Apply 1-site TDVP[https://doi.org/10.1103/PhysRevB.94.165116] sweep to perform time evolution for DenseMPS (both MPS and MPO) with step length dt. Env is the 3-layer environment ⟨Ψ|H|Ψ⟩.

 TDVPSweep1!(Env::SparseEnvironment{L,3,T}, dt::Number; kwargs...)

Wrap TDVPSweep1! with a symmetric integrator, i.e., sweeping from left to right and then from right to left with the same step length dt / 2.

Kwargs

 krylovalg::KrylovKit.KrylovAlgorithm = TDVPDefaultLanczos
 trunc::TruncationType = notrunc()
 GCstep::Bool = false
 GCsweep::Bool = false
 verbose::Int64 = 0
 CBEAlg::CBEAlgorithm = NoCBE()
FiniteMPS.TDVPSweep2!Method
 TDVPSweep2!(Env::SparseEnvironment{L,3,T},
      dt::Number,
      direction::SweepDirection;
      kwargs...) -> info, TimerSweep

Apply left-to-right or right-to-left 2-site TDVP[https://doi.org/10.1103/PhysRevB.94.165116] sweep to perform time evolution for DenseMPS (both MPS and MPO) with step length dt. Env is the 3-layer environment ⟨Ψ|H|Ψ⟩.

 TDVPSweep2!(Env::SparseEnvironment{L,3,T}, dt::Number; kwargs...)

Wrap TDVPSweep2! with a symmetric integrator, i.e., sweeping from left to right and then from right to left with the same step length dt / 2.

Kwargs

 krylovalg::KrylovKit.KrylovAlgorithm = TDVPDefaultLanczos
 trunc::TruncationType = truncbelow(MPSDefault.tol) & truncdim(MPSDefault.D)
 GCstep::Bool = false
 GCsweep::Bool = false
 verbose::Int64 = 0
FiniteMPS.action0Method
 action0(obj::SparseProjectiveHamiltonian{0}, x::MPSTensor{2}; kwargs...) -> ::MPSTensor

Action of 0-site projective Hamiltonian on the rank-2 bond local tensors.

FiniteMPS.action1Method
 action1(obj::SparseProjectiveHamiltonian{1}, x::MPSTensor; kwargs...) -> ::MPSTensor

Action of 1-site projective Hamiltonian on the 1-site local tensors.

FiniteMPS.action2Method
 action2(obj::SparseProjectiveHamiltonian{2}, x::CompositeMPSTensor{2, T}; kwargs...) -> ::CompositeMPSTensor{2, T}

Action of 2-site projective Hamiltonian on the 2-site local tensors, wrapped by CompositeMPSTensor{2, T} where T<:NTuple{2,MPSTensor}.

 action2(obj::IdentityProjectiveHamiltonian{2}, x::CompositeMPSTensor{2, T}; kwargs...) -> ::CompositeMPSTensor{2, T}

Special case for IdentityProjectiveHamiltonian.

FiniteMPS.addIntr!Method
 addIntr!(Root::InteractionTreeNode,
      Op::NTuple{N,AbstractTensorMap},
      si::NTuple{N,Int64},
      strength::Number;
      kwargs...)

Generic function to add an N-site interaction via addIntr1!, addIntr2! and addIntr4!.

Kwargs

 Z::Union{Nothing,AbstractTensorMap}=nothing
 name::NTuple{N,Union{Symbol,String}}

Detailed usage of kwargs see addIntr1!, addIntr2! and addIntr4!.

 Obs::Bool = false

Obs == true means this interaction is used for calculating observables, and thus the name and si information will be stored in the last node additionally.

FiniteMPS.addIntr1!Method
 addIntr1!(Root::InteractionTreeNode,
      Op::AbstractTensorMap,
      si::Int64,
      strength::Number;
      Obs::Bool = false,
      name::Union{Symbol,String} = :O) -> nothing

 addIntr1!(Tree::InteractionTree, args...) = addIntr1!(Tree.Root.children[1], args...)

Add an on-site term Op at site si to a given interaction tree.

 addIntr1!(Root::InteractionTreeNode,
      O::LocalOperator,
      strength::Number;
      value = nothing) -> nothing

Expert version, each method finally reduces to this one. The value will be stored in the last node.

FiniteMPS.addIntr2!Method
 addIntr2!(Root::InteractionTreeNode,
      Op::NTuple{2,AbstractTensorMap},
      si::NTuple{2,Int64},
      strength::Number;
      Obs::Bool = false,
      Z::Union{Nothing,AbstractTensorMap} = nothing,
      name::NTuple{2,Union{Symbol,String}} = (:A, :B)) -> nothing

 addIntr2!(Tree::InteractionTree, args...) = addIntr2!(Tree.Root.children[1], args...)

Add a two-site interaction Op at site si (2tuple) to a given interaction tree. If Z is given, assume Op is ferminic operator and add Z automatically.

 addIntr2!(Root::InteractionTreeNode,
      OL::LocalOperator,
      OR::LocalOperator,
      strength::Number,
      Z::Union{Nothing,AbstractTensorMap};
      value = nothing) -> nothing

Expert version, each method finally reduces to this one. The value will be stored in the last node.

Note if OL.si == OR.si, it will recurse to addIntr1! automatically.

FiniteMPS.addIntr4!Method
 addIntr4!(Root::InteractionTreeNode,
      Op::NTuple{4,AbstractTensorMap},
      si::NTuple{4,Int64},
      strength::Number;
      Obs::Bool = false,
      Z::Union{Nothing,AbstractTensorMap} = nothing,
      name::NTuple{4,Union{Symbol,String}} = (:A, :B, :C, :D)) -> nothing

 addIntr4!(Tree::InteractionTree, args...) = addIntr4!(Tree.Root.children[1], args...)

Add a 4-site interaction Op at site si (4tuple) to a given interaction tree. If Z is given, assume each operator in tuple Op is ferminic operator and add Z automatically.

 addIntr4!(Root::InteractionTreeNode,
      A::LocalOperator,
      B::LocalOperator,
      C::LocalOperator,
      D::LocalOperator,
      strength::Number,
      Z::Union{Nothing,AbstractTensorMap};
      value = nothing) -> nothing

Expert version, each method finally reduces to this one.

Note if there exist repeated si, it will recurse to addIntr2! or addIntr3!(TODO) automatically.

FiniteMPS.addObs!Method
 addObs!(Tree::ObservableTree{M},
      Op::NTuple{N,AbstractTensorMap},
      si::NTuple{N,Int64}
      n::Int64 = 1; 
      Z::Union{Nothing,AbstractTensorMap} = nothing,
      name::NTuple{N,Union{Symbol,String}} = (:A, :B, ...))

Add a term to the n-th root of ObservableTree{M}. Detailed usage see addIntr!.

Warning: a same name can be given to two local operators iff they are exactly the same, otherwise, it will confuse the convert function when trying to extract the values stored in the tree to a dictionary. For example, you can simply name SzSz correlation as (:S, :S). However this name is inappropriate for S+S- correlation.

Warning: there is a known issue that only one permutation will be calculated if there exist multiple permutations with the same operator and name. For example, if you add (:Sz,:Sz) with sites (1, 2) and (2, 1), only one of them will appear in the final result. Changing the name to distinguish the 1st and 2nd Sz by using e.g. (:Sz1, :Sz2) can solve this issue. However, a more elegant solution is to avoid adding both of them as we know the expected value must be the same. Thus, we will not prioritize fixing this issue in the near future.

FiniteMPS.addchild!Method
 addchild!(node::InteractionTreeNode, child::InteractionTreeNode) -> nothing

Add a child to a given node.

 addchild!(node::InteractionTreeNode, Op::AbstractLocalOperator [, value]) -> nothing

Initialize a child node with Op, and add it to the given node.

Note we will use similar(for array-like types) or zero(for other types) to initialize the field value in default case, which may lead to a "no method matching" error.

FiniteMPS.calObs!Method
 calObs!(Tree::ObservableTree, Ψ::AbstractMPS{L}; kwargs...) -> Timer::TimerOutput

Calculate observables respect to state Ψ, the info to tell which observables to calculate is stored in Tree. The results are stored in each leaf node of Tree. Note the value in each node will be in-place updated, so do not call this function twice with the same Tree object.

Kwargs

 serial::Bool = false

Force to compute in serial mode, usually used for debugging.

FiniteMPS.canonicalize!Method
 canonicalize!(obj::AbstractEnvironment, 
      siL::Int64
      [, siR::Int64 = siL]; kwargs...) -> obj::AbstractEnvironment

Canonicalize the environment s.t. at least El[i ≤ siL] and Er[i ≥ siR] are valid.

Kwargs

 free::Bool = true

If true, call free!(obj) to free the local environment tensors which are no longer required. Details see free!.

FiniteMPS.canonicalize!Method
 canonicalize!(obj::AbstractMPS, 
      siL::Int64
      [, siR::Int64 = siL]; kwargs...) -> obj::AbstractMPS

Canonicalize the MPS s.t. all sites ≤ siL are left-canonical, all sites ≥ siR are right-canonical.

kwargs will be propagated to leftorth and rightorth to determine how to truncate the SVD spectra.

FiniteMPS.cleanup!Method
 cleanup!(d::SerializedElementArray, n::Int64) -> nothing

Use rm to cleanup the file corresponding to the n-th element of d::SerializedElementArray in disk manually.

 cleanup!(d::SerializedElementArray, lsn::AbstractArray{Int64}) -> nothing

Vector version of the above usage.

 cleanup!(d::SerializedElementArray) -> nothing

Cleanup all files of d::SerializedElementArray manually. Note setindex! cannot be applied to d after this operation.

FiniteMPS.coefMethod
 coef(obj::DenseMPS) -> ::F

Interface of DenseMPS, return the global coefficient, where F is the number type of given MPS.

FiniteMPS.connection!Method
 connection!(obj::SparseEnvironment; 
      kwargs...) -> ::Matrix

Return the connection ⟨∇⟨Hᵢ⟩, ∇⟨Hⱼ⟩⟩ where H₀, H₁, ⋯, Hₙ are the components of the total Hamiltonian, decomposed according to the boundary left environment.

Note the state obj[3] must be right canonicalized. After this funcation, the canonical center will move to the right boundary.

kwargs

 moveback::Bool = false

Move the canonical center back to the left boundary if true.

FiniteMPS.dataMethod
 data(A::AbstractTensorMap) -> collection of data
 data(A::AdjointTensorMap) = data(A.parent)

Interface of AbstractTensorMap, return the data of a given tensor.

FiniteMPS.free!Method
 free!(obj::AbstractEnvironment; 
      siL::AbstractVector{Int64} = obj.Center[1] + 1 : L,
      siR::AbstractVector{Int64} = 1 : obj.Center[2] - 1     
 ) -> nothing

Free the local tensors in El[siL] and Er[siR].

FiniteMPS.getOpNameMethod
 getOpName(::AbstractLocalOperator) -> ::String

Interface of AbstractLocalOperator, return "I" for IdentityOperator and O.name for O::LocalOperator.

FiniteMPS.getPhysSpaceMethod
 getPhysSpace(O::AbstractLocalOperator) -> ::VectorSpace

Interface of AbstractLocalOperator, return the local physical space.

FiniteMPS.hastagMethod
 hastag(::AbstractLocalOperator) -> ::Bool

Check if this LocalOperator has field tag or not.

FiniteMPS.identityMPOMethod
 identityMPO(::Type{T} = Float64, L::Int64, pspace::AbstractVector; kwargs...)

Construct an identity MPO where the physical spaces are informed by a length L vertor of VectorSpace.

 identityMPO(::Type{T} = Float64, L::Int64, pspace::VectorSpace; kwargs...)

Assume the all the physical spaces are the same.

 identityMPO(obj::DenseMPS{L, T}; kwargs...)

Deduce the scalar type T and physical spaces from a MPS/MPO.

Kwargs

 disk::Bool = false
 bspace::VectorSpace

Space of left boundary bond, default = trivial space.

FiniteMPS.initialize!Method
 initialize!(obj::AbstractEnvironment; kwargs...)

Initialize the boundary environment tensors, i.e. El[1] and Er[L].

Kwargs

 El::Union{SimpleLeftTensor, SparseLeftTensor}
 Er::Union{SimpleRightTensor, SparseRightTensor}

Directly give El or Er, otherwise, use _defaultEl or _defaultEr to generate one.

 free::Bool = false

If true, call free!(obj) to free the local environment tensors which are no longer required. Details see free!.

FiniteMPS.issparseMethod
 issparse(::AbstractMPS) -> ::Bool

Check if a MPS/MPO object is sparse, e.g. ::MPS -> false, ::SparseMPO -> true.

FiniteMPS.istrivialMethod
 istrivial(V::VectorSpace) -> ::Bool

Check if a given VectorSpace is trivial.

Examples

julia> istrivial(Rep[U₁ × SU₂]((1, 1/2) => 2))
false

julia> istrivial(Rep[U₁ × SU₂]((0, 0) => 2))
false

julia> istrivial(Rep[U₁ × SU₂]((0, 0) => 1))
true
FiniteMPS.manualGCMethod
 manualGC([T::TimerOutput]) -> nothing

Manually call GC.gc() on all workers and use @timeit to collect the time cost if T::TimerOutput is provided.

FiniteMPS.noise!Method
 noise!(A::CompositeMPSTensor{2}, σ::Real)

Apply noise to a given 2-site local tensor by contracting a d×d random isometry to it.

FiniteMPS.oplusEmbedMethod
 oplusEmbed(lsV::Vector{<:GradedSpace};
      reverse::Bool=false) -> lsEmbed::Vector{<:AbstractTensorMap}

Return the embedding maps from vectors in lsV to their direct sum space, with the same order as lsV. If reverse == true, return the submersions from the direct sum space to the vectors instead.

 oplusEmbed(A::AbstractTensorMap,
      B::AbstractTensorMap,
      idx::Int64) -> EmbA::TensorMap, EmbB::TensorMap

Return the 2 embedding maps (from A and B) to the direct sum space (or their adjoint) corresponding to idx.

FiniteMPS.pushleft!Method
 pushleft!(::AbstractEnvironment)

Push left the given environment object, i.e. Center == [i, j] to [i, j - 1].

FiniteMPS.pushright!Method
 pushright!(::AbstractEnvironment)

Push right the given environment object, i.e. Center == [i, j] to [i + 1, j].

FiniteMPS.randMPSMethod
 randMPS([::Type{T},]  
      pspace::Vector{VectorSpace},
      aspace::Vector{VectorSpace};
      kwargs...) -> MPS{L}

Generate a length L random MPS with given length L vector pspace and aspace. T = Float64(default) or ComplexF64 is the number type. Note the canonical center is initialized to the first site.

 randMPS([::Type{T},] L::Int64, pspace::VectorSpace, apsace::VectorSpace; kwargs...) -> MPS{L}

Assume the same pspace and aspace, except for the boundary bond, which is assumed to be trivial.

FiniteMPS.randStiefelMethod
 randStiefel([Dist::UniformDistribution, ] n::Int64, k::Int64) -> ::Matrix (n × k)
 randStiefel(Dist::GaussianDistribution) -> ::Matrix (n × k)

Random sample on Stiefel manifold Vₙₖ with given distribution.

Default Dist = UniformDistribution{Float64}.

FiniteMPS.randisometry!Method
 randisometry!(A::AbstractTensorMap; kwargs...)

Randomize the data of tensor A based on randStiefel.

Keywords

 Dist => ::Symbol

:Uniform or :Gaussian ( == :Normal). If σ is given, default = :Gaussian, otherwise :Uniform.

 σ => ::Real

σ of gaussian distribution. No default value. It will throw an error if σ is not given and Dist = :Gaussian.

FiniteMPS.randisometryMethod
 randisometry([::Type{T},] codom::VectorSpace, dom::VectorSpace = codom; kwargs...)
 randisometry([::Type{T},] A::AbstractTensorMap; kwargs...)

Generate random tensors based on randStiefel. Valid kwargs please see the mutating version randisometry!.

FiniteMPS.scalar!Method
 scalar!(obj::AbstractEnvironment; kwargs...) -> ::Number

Fully contract the total tensor network to get a scalar.

Note this may change the Center of the environment obj.

Kwargs

 normalize::Bool = false

If true, calculate ⟨Ψ₁|H|Ψ₂⟩/⟨Ψ₁|Ψ₂⟩ instead of ⟨Ψ₁|H|Ψ₂⟩ for example.

 split::Bool = false

Split the value into each contribution of each left boundary environment if true. Thus, return a vector instead of a scalar in this case.

 tmp::Bool = false

tmp == true means the environment is temporary, and thus we will free the local environment tensors which are no longer required.

FiniteMPS.trivialMethod
 trivial(V::VectorSpace) -> ::VectorSpace

Return the trivial space of a given vector space.

Examples

julia> trivial(Rep[U₁ × SU₂]((1, 1/2) => 2))
Rep[U₁ × SU₂]((0, 0)=>1)

julia> trivial(ℂ^2)
ℂ^1
LinearAlgebra.axpby!Method
 axpby!(α::Number, x::DenseMPS, β::Number, y::DenseMPS; kwargs...)

Compute y = α*x + β*y variationally via 2-site update, where x and y are dense MPS/MPO. Note 'x' cannot reference to the same MPS/MPO with y.

Kwargs

 trunc::TruncationScheme = truncbelow(MPSDefault.tol) & truncdim(MPSDefault.D)
 GCstep::Bool = false
 GCsweep::Bool = false
 maxiter::Int64 = 8
 disk::Bool = false
 tol::Float64 = 1e-8
 verbose::Int64 = 0
LinearAlgebra.mul!Method
 mul!(C::DenseMPS, A::SparseMPO, B::DenseMPS, α::Number, β::Number; kwargs...)

Compute C = α A*B + β C variationally via 2-site update, where A is a sparse MPO, B and C are dense MPS/MPO. Note 'B' cannot reference to the same MPS/MPO with C.

 mul!(C::DenseMPS, A::SparseMPO, B::DenseMPS; kwargs...)

Compute C = A*B by letting α = 1 and β = 0.

Kwargs

 trunc::TruncationScheme = truncbelow(MPSDefault.tol) & truncdim(MPSDefault.D)
 GCstep::Bool = false
 GCsweep::Bool = false
 maxiter::Int64 = 8
 disk::Bool = false
 tol::Float64 = 1e-8
 verbose::Int64 = 0
 lsnoise::AbstractVector{Float64} = Float64[]
LinearAlgebra.normMethod
 norm(obj::DenseMPS) -> ::Float64

Return the inner-induced norm. Note we assume the MPS satisfies a canonical form and the center tensor is normalized, hence the norm is just abs(c).

LinearAlgebra.normalize!Method
 normalize!(obj::DenseMPS) -> obj

Normalize a given MPS according to inner-induced norm.

Note we assume the MPS satisfies a canonical form and the center tensor is normalized, hence we only normalize c.

LinearAlgebra.rankMethod
 rank(A::AbstractTensorMap) -> ::Int64

Return the rank of a given tensor.

 rank(A::AbstractTensorMap, idx::Int64) -> ::Int64

Return the rank corresponding to codomain (idx = 1) or domain (idx = 2).

TensorKit.dimMethod
 dim(A::AbstractTensorMap, idx::Int64) -> (D, DD)::NTuple{2, Int64}

Return the dimension of a given index of tensor A.

D is the number of multiplets, DD is the number of equivalent no symmetry states. Note D == DD for abelian groups.

TensorKit.fuseMethod
 fuse(El::SimpleLeftTensor) -> iso::AbstractTensorMap
 fuse(El::SimpleRightTensor) -> iso::AbstractTensorMap

Return the isometry to fuse the top 2 legs.

 fuse(lsEl::SparseLeftTensor) -> Vector{AbstractTensorMap} 
 fuse(lsEl::SparseRightTensor) -> Vector{AbstractTensorMap}

Additionally embed the isometry to the direct sum space of all channels.

TensorKit.leftorthMethod
 leftorth(::CompositeMPSTensor{2, ...};
      trunc = notrunc(),
      kwargs...) -> Q::AbstractTensorMap, R::AbstractTensorMap, info::BondInfo

Split a 2-site local tensor s.t. the left one is canonical.

TensorKit.leftorthMethod
 leftorth(A::MPSTensor; 
      trunc = notrunc(),
      kwargs...) -> Q::AbstractTensorMap, R::AbstractTensorMap, info::BondInfo

Left canonicalize a on-site MPS tensor.

If trunc = notrunc(), use TensorKit.leftorth, otherwise, use TensorKit.tsvd. Propagate kwargs to the TensorKit functions.

TensorKit.rightorthMethod
 rightorth(::CompositeMPSTensor{2, ...};
      trunc = notrunc(),
      kwargs...) -> L::AbstractTensorMap, Q::AbstractTensorMap, info::BondInfo

Split a 2-site local tensor s.t. the right one is canonical.

TensorKit.rightorthMethod
 rightorth(A::MPSTensor;
      trunc = notrunc(),
      kwargs...) -> L::AbstractTensorMap, Q::AbstractTensorMap, info::BondInfo

Right canonicalize a on-site MPS tensor.

If trunc = notrunc(), use TensorKit.rightorth, otherwise, use TensorKit.tsvd. Propagate kwargs to the TensorKit functions.

TensorKit.tsvdMethod
 tsvd(::CompositeMPSTensor{2, ...}; kwargs...) 
      -> u::AbstractTensorMap, s::AbstractTensorMap, vd::AbstractTensorMap, info::BondInfo

Use SVD to split a 2-site local tensor, details see TensorKit.tsvd.

TensorKit.tsvdMethod
 tsvd(A::AbstractTensorWrapper,
      p₁::NTuple{N₁,Int64},
      p₂::NTuple{N₂,Int64};
      kwargs...) 
      -> u::AbstractTensorMap, s::AbstractTensorMap, vd::AbstractTensorMap, info::BondInfo

Wrap TensorKit.tsvd, return BondInfo struct instead of truncation error ϵ.

VectorInterface.innerMethod
 inner(A::DenseMPS, B::DenseMPS)

Return the inner product ⟨A, B⟩ between MPS/MPO A and B.

VectorInterface.scalartypeMethod
 scalartype(obj::SparseMPO) -> Float64/ComplexF64

Return the scalar type of given SparseMPO. Note return Float64 iff all local tensors are real.

FiniteMPS.U₁U₁tJFermionModule
 module U₁U₁tJFermion

Prepare some commonly used objects for U₁×U₁ tJ fermions, i.e. local d = 3 Hilbert space without double occupancy.

Behaviors of all operators are the same as U₁U₁Fermion up to the projection, details please see U₁U₁Fermion.

FiniteMPS.NoSymSpinOneHalfModule
 module NoSymSpinOneHalf

Prepare the local space of U₁ spin-1/2. Basis convention is {|↑⟩, |↓⟩}.

Fields

 pspace::VectorSpace

Local d = 2 Hilbert space.

Sz::TensorMap
Sx::TensorMap
Sy::TensorMap

Rank-2 spin-1/2 operators.

S₊::TensorMap

Rank-2 spin-plus operator S₊ = Sx + iSy. S₋::TensorMap Rank-2 spin-minus operator S₋ = Sx - iSy.

FiniteMPS.U₁SpinModule
 module U₁Spin

Prepare the local space of U₁ spin-1/2.

Fields

 pspace::VectorSpace

Local d = 2 Hilbert space.

Sz::TensorMap

Rank-2 spin-z operator Sz = (n↑ - n↓)/2.

S₊₋::NTuple{2, TensorMap}
S₋₊::NTuple{2, TensorMap}

Two rank-3 operators of S₊₋ and S₋₊ interaction. Note Heisenberg S⋅S = SzSz + (S₊₋ + S₋₊)/2.