MPS and MPO

Types

ITensors.MPSType
MPS

A finite size matrix product state type. Keeps track of the orthogonality center.

source
ITensors.MPOType
MPO

A finite size matrix product operator type. Keeps track of the orthogonality center.

source

MPS Constructors

ITensors.MPSMethod
MPS(N::Int)

Construct an MPS with N sites with default constructed ITensors.

source
ITensors.MPSMethod
MPS([::Type{ElT} = Float64, ]sites)

Construct an MPS filled with Empty ITensors of type ElT from a collection of indices.

source
ITensors.randomMPSMethod
randomMPS(sites::Vector{<:Index}; linkdim=1)

Construct a random MPS with link dimension linkdim of type Float64.

source
ITensors.randomMPSMethod
randomMPS(::Type{ElT<:Number}, sites::Vector{<:Index}; linkdim=1)

Construct a random MPS with link dimension linkdim of type ElT.

source
ITensors.randomMPSMethod
randomMPS(sites::Vector{<:Index}, state; linkdim=1)

Construct a real, random MPS with link dimension linkdim, made by randomizing an initial product state specified by state. This version of randomMPS is necessary when creating QN-conserving random MPS (consisting of QNITensors). The initial state array provided determines the total QN of the resulting random MPS.

source
ITensors.productMPSMethod
productMPS(sites::Vector{<:Index},states)

Construct a product state MPS having site indices sites, and which corresponds to the initial state given by the array states. The states array may consist of either an array of integers or strings, as recognized by the state function defined for the relevant Index tag type.

Examples

N = 10
sites = siteinds("S=1/2",N)
states = [isodd(n) ? "Up" : "Dn" for n=1:N]
psi = productMPS(sites,states)
source
ITensors.productMPSMethod
productMPS(::Type{T},
           sites::Vector{<:Index},
           states::Union{Vector{String},
                         Vector{Int},
                         String,
                         Int})

Construct a product state MPS of element type T, having site indices sites, and which corresponds to the initial state given by the array states. The input states may be an array of strings or an array of ints recognized by the state function defined for the relevant Index tag type. In addition, a single string or int can be input to create a uniform state.

Examples

N = 10
sites = siteinds("S=1/2", N)
states = [isodd(n) ? "Up" : "Dn" for n=1:N]
psi = productMPS(ComplexF64, sites, states)
phi = productMPS(sites, "Up")
source
ITensors.productMPSMethod
productMPS(ivals::Vector{<:IndexVal})

Construct a product state MPS with element type Float64 and nonzero values determined from the input IndexVals.

source
ITensors.productMPSMethod
productMPS(::Type{T<:Number}, ivals::Vector{<:IndexVal})

Construct a product state MPS with element type T and nonzero values determined from the input IndexVals.

source

MPO Constructors

ITensors.MPOMethod
MPO(N::Int)

Make an MPO of length N filled with default ITensors.

source
ITensors.MPOMethod
MPO([::Type{ElT} = Float64}, ]sites, ops::Vector{String})

Make an MPO with pairs of sites s[i] and s[i]' and operators ops on each site.

source
ITensors.MPOMethod
MPO([::Type{ElT} = Float64, ]sites, op::String)

Make an MPO with pairs of sites s[i] and s[i]' and operator op on every site.

source
ITensors.MPOMethod
MPO(A::MPS; kwargs...)

For an MPS |A>, make the MPO |A><A|. Keyword arguments like cutoff can be used to truncate the resulting MPO.

source

Properties

ITensors.fluxMethod
flux(M::MPS)

flux(M::MPO)

totalqn(M::MPS)

totalqn(M::MPO)

For an MPS or MPO which conserves quantum numbers, compute the total QN flux. For a tensor network such as an MPS or MPO, the flux is the sum of fluxes of each of the tensors in the network. The name totalqn is an alias for flux.

source
ITensors.hasqnsMethod
hasqns(M::MPS)

hasqns(M::MPO)

Return true if the MPS or MPO has tensors which carry quantum numbers.

source

Obtaining and finding indices

ITensors.siteindsMethod
siteinds(commoninds, A::MPO, B::MPS, j::Integer; kwargs...)
siteinds(commonind, A::MPO, B::MPO, j::Integer; kwargs...)

Get the site index (or indices) of the jth MPO tensor of A that is shared with MPS/MPO B.

source
ITensors.siteindsMethod
siteinds(uniqueinds, A::MPO, B::MPS, j::Integer; kwargs...)
siteinds(uniqueind, A::MPO, B::MPS, j::Integer; kwargs...)

Get the site index (or indices) of MPO A that is unique to A (not shared with MPS/MPO B).

source
ITensors.findsiteFunction
findsite(M::Union{MPS, MPO}, is)

Return the first site of the MPS or MPO that has at least one Index in common with the Index or collection of indices is.

To find all sites with common indices with is, use the findsites function.

Examples

s = siteinds("S=1/2", 5)
ψ = randomMPS(s)
findsite(ψ, s[3]) == 3
findsite(ψ, (s[3], s[4])) == 3

M = MPO(s)
findsite(M, s[4]) == 4
findsite(M, s[4]') == 4
findsite(M, (s[4]', s[4])) == 4
findsite(M, (s[4]', s[3])) == 3
source
ITensors.findsitesFunction
findsites(M::Union{MPS, MPO}, is)

Return the sites of the MPS or MPO that have indices in common with the collection of site indices is.

Examples

s = siteinds("S=1/2", 5)
ψ = randomMPS(s)
findsites(ψ, s[3]) == [3]
findsites(ψ, (s[4], s[1])) == [1, 4]

M = MPO(s)
findsites(M, s[4]) == [4]
findsites(M, s[4]') == [4]
findsites(M, (s[4]', s[4])) == [4]
findsites(M, (s[4]', s[3])) == [3, 4]
source
ITensors.firstsiteindsFunction
firstsiteinds(M::MPO; kwargs...)

Get a Vector of the first site Index found on each site of M.

By default, it finds the first site Index with prime level 0.

source
ITensors.linkindMethod
linkind(M::MPS, j::Integer)
linkind(M::MPO, j::Integer)

Get the link or bond Index connecting the MPS or MPO tensor on site j to site j+1.

If there is no link Index, return nothing.

source
ITensors.siteindMethod
siteind(M::MPS, j::Int; kwargs...)

Get the first site Index of the MPS. Return nothing if none is found.

source
ITensors.siteindMethod
siteind(::typeof(first), M::Union{MPS,MPO}, j::Integer; kwargs...)

Return the first site Index found on the MPS or MPO (the first Index unique to the jth MPS/MPO tensor).

You can choose different filters, like prime level and tags, with the kwargs.

source
ITensors.siteindsMethod
siteinds(M::MPS)
siteinds(::typeof(first), M::MPS)

Get a vector of the first site Index found on each tensor of the MPS.

siteinds(::typeof(only), M::MPS)

Get a vector of the only site Index found on each tensor of the MPS. Errors if more than one is found.

siteinds(::typeof(all), M::MPS)

Get a vector of the all site Indices found on each tensor of the MPS. Returns a Vector of IndexSets.

source
ITensors.siteindMethod
siteind(M::MPO, j::Int; plev = 0, kwargs...)

Get the first site Index of the MPO found, by default with prime level 0.

source
ITensors.siteindsMethod
siteinds(M::MPO; kwargs...)

Get a Vector of IndexSets of all the site indices of M.

source
ITensors.siteindsMethod
siteinds(M::Union{MPS, MPO}}, j::Integer; kwargs...)

Return the site Indices found of the MPO or MPO at the site j as an IndexSet.

Optionally filter prime tags and prime levels with keyword arguments like plev and tags.

source

Priming and tagging

ITensors.primeMethod
prime[!](M::MPS, args...; kwargs...)
prime[!](M::MPO, args...; kwargs...)

Apply prime to all ITensors of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors. Alternatively apply the function in-place.

source
ITensors.primeMethod
prime[!](siteinds, M::MPS, args...; kwargs...)
prime[!](siteinds, M::MPO, args...; kwargs...)

Apply prime to all site indices of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors.

source
ITensors.primeMethod
prime[!](linkinds, M::MPS, args...; kwargs...)
prime[!](linkinds, M::MPO, args...; kwargs...)

Apply prime to all link indices of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors.

source
ITensors.primeMethod
prime[!](siteinds, commoninds, M1::MPO, M2::MPS, args...; kwargs...)
prime[!](siteinds, commoninds, M1::MPO, M2::MPO, args...; kwargs...)

Apply prime to the site indices that are shared by M1 and M2.

Returns new MPSs/MPOs. The ITensors of the MPSs/MPOs will be a view of the storage of the original ITensors.

source
ITensors.primeMethod
prime[!](siteinds, uniqueinds, M1::MPO, M2::MPS, args...; kwargs...)

Apply prime to the site indices of M1 that are not shared with M2. Returns new MPSs/MPOs.

The ITensors of the MPSs/MPOs will be a view of the storage of the original ITensors.

source
ITensors.setprimeMethod
setprime[!](M::MPS, args...; kwargs...)
setprime[!](M::MPO, args...; kwargs...)

Apply setprime to all ITensors of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors. Alternatively apply the function in-place.

source
ITensors.setprimeMethod
setprime[!](siteinds, M::MPS, args...; kwargs...)
setprime[!](siteinds, M::MPO, args...; kwargs...)

Apply setprime to all site indices of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors.

source
ITensors.setprimeMethod
setprime[!](linkinds, M::MPS, args...; kwargs...)
setprime[!](linkinds, M::MPO, args...; kwargs...)

Apply setprime to all link indices of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors.

source
ITensors.setprimeMethod
setprime[!](siteinds, commoninds, M1::MPO, M2::MPS, args...; kwargs...)
setprime[!](siteinds, commoninds, M1::MPO, M2::MPO, args...; kwargs...)

Apply setprime to the site indices that are shared by M1 and M2.

Returns new MPSs/MPOs. The ITensors of the MPSs/MPOs will be a view of the storage of the original ITensors.

source
ITensors.setprimeMethod
setprime[!](siteinds, uniqueinds, M1::MPO, M2::MPS, args...; kwargs...)

Apply setprime to the site indices of M1 that are not shared with M2. Returns new MPSs/MPOs.

The ITensors of the MPSs/MPOs will be a view of the storage of the original ITensors.

source
ITensors.noprimeMethod
noprime[!](M::MPS, args...; kwargs...)
noprime[!](M::MPO, args...; kwargs...)

Apply noprime to all ITensors of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors. Alternatively apply the function in-place.

source
ITensors.noprimeMethod
noprime[!](siteinds, M::MPS, args...; kwargs...)
noprime[!](siteinds, M::MPO, args...; kwargs...)

Apply noprime to all site indices of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors.

source
ITensors.noprimeMethod
noprime[!](linkinds, M::MPS, args...; kwargs...)
noprime[!](linkinds, M::MPO, args...; kwargs...)

Apply noprime to all link indices of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors.

source
ITensors.noprimeMethod
noprime[!](siteinds, commoninds, M1::MPO, M2::MPS, args...; kwargs...)
noprime[!](siteinds, commoninds, M1::MPO, M2::MPO, args...; kwargs...)

Apply noprime to the site indices that are shared by M1 and M2.

Returns new MPSs/MPOs. The ITensors of the MPSs/MPOs will be a view of the storage of the original ITensors.

source
ITensors.noprimeMethod
noprime[!](siteinds, uniqueinds, M1::MPO, M2::MPS, args...; kwargs...)

Apply noprime to the site indices of M1 that are not shared with M2. Returns new MPSs/MPOs.

The ITensors of the MPSs/MPOs will be a view of the storage of the original ITensors.

source
ITensors.addtagsMethod
addtags[!](M::MPS, args...; kwargs...)
addtags[!](M::MPO, args...; kwargs...)

Apply addtags to all ITensors of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors. Alternatively apply the function in-place.

source
ITensors.addtagsMethod
addtags[!](siteinds, M::MPS, args...; kwargs...)
addtags[!](siteinds, M::MPO, args...; kwargs...)

Apply addtags to all site indices of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors.

source
ITensors.addtagsMethod
addtags[!](linkinds, M::MPS, args...; kwargs...)
addtags[!](linkinds, M::MPO, args...; kwargs...)

Apply addtags to all link indices of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors.

source
ITensors.addtagsMethod
addtags[!](siteinds, commoninds, M1::MPO, M2::MPS, args...; kwargs...)
addtags[!](siteinds, commoninds, M1::MPO, M2::MPO, args...; kwargs...)

Apply addtags to the site indices that are shared by M1 and M2.

Returns new MPSs/MPOs. The ITensors of the MPSs/MPOs will be a view of the storage of the original ITensors.

source
ITensors.addtagsMethod
addtags[!](siteinds, uniqueinds, M1::MPO, M2::MPS, args...; kwargs...)

Apply addtags to the site indices of M1 that are not shared with M2. Returns new MPSs/MPOs.

The ITensors of the MPSs/MPOs will be a view of the storage of the original ITensors.

source
ITensors.removetagsMethod
removetags[!](M::MPS, args...; kwargs...)
removetags[!](M::MPO, args...; kwargs...)

Apply removetags to all ITensors of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors. Alternatively apply the function in-place.

source
ITensors.removetagsMethod
removetags[!](siteinds, M::MPS, args...; kwargs...)
removetags[!](siteinds, M::MPO, args...; kwargs...)

Apply removetags to all site indices of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors.

source
ITensors.removetagsMethod
removetags[!](linkinds, M::MPS, args...; kwargs...)
removetags[!](linkinds, M::MPO, args...; kwargs...)

Apply removetags to all link indices of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors.

source
ITensors.removetagsMethod
removetags[!](siteinds, commoninds, M1::MPO, M2::MPS, args...; kwargs...)
removetags[!](siteinds, commoninds, M1::MPO, M2::MPO, args...; kwargs...)

Apply removetags to the site indices that are shared by M1 and M2.

Returns new MPSs/MPOs. The ITensors of the MPSs/MPOs will be a view of the storage of the original ITensors.

source
ITensors.removetagsMethod
removetags[!](siteinds, uniqueinds, M1::MPO, M2::MPS, args...; kwargs...)

Apply removetags to the site indices of M1 that are not shared with M2. Returns new MPSs/MPOs.

The ITensors of the MPSs/MPOs will be a view of the storage of the original ITensors.

source
ITensors.replacetagsMethod
replacetags[!](M::MPS, args...; kwargs...)
replacetags[!](M::MPO, args...; kwargs...)

Apply replacetags to all ITensors of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors. Alternatively apply the function in-place.

source
ITensors.replacetagsMethod
replacetags[!](siteinds, M::MPS, args...; kwargs...)
replacetags[!](siteinds, M::MPO, args...; kwargs...)

Apply replacetags to all site indices of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors.

source
ITensors.replacetagsMethod
replacetags[!](linkinds, M::MPS, args...; kwargs...)
replacetags[!](linkinds, M::MPO, args...; kwargs...)

Apply replacetags to all link indices of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors.

source
ITensors.replacetagsMethod
replacetags[!](siteinds, commoninds, M1::MPO, M2::MPS, args...; kwargs...)
replacetags[!](siteinds, commoninds, M1::MPO, M2::MPO, args...; kwargs...)

Apply replacetags to the site indices that are shared by M1 and M2.

Returns new MPSs/MPOs. The ITensors of the MPSs/MPOs will be a view of the storage of the original ITensors.

source
ITensors.replacetagsMethod
replacetags[!](siteinds, uniqueinds, M1::MPO, M2::MPS, args...; kwargs...)

Apply replacetags to the site indices of M1 that are not shared with M2. Returns new MPSs/MPOs.

The ITensors of the MPSs/MPOs will be a view of the storage of the original ITensors.

source
ITensors.settagsMethod
settags[!](M::MPS, args...; kwargs...)
settags[!](M::MPO, args...; kwargs...)

Apply settags to all ITensors of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors. Alternatively apply the function in-place.

source
ITensors.settagsMethod
settags[!](siteinds, M::MPS, args...; kwargs...)
settags[!](siteinds, M::MPO, args...; kwargs...)

Apply settags to all site indices of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors.

source
ITensors.settagsMethod
settags[!](linkinds, M::MPS, args...; kwargs...)
settags[!](linkinds, M::MPO, args...; kwargs...)

Apply settags to all link indices of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors.

source
ITensors.settagsMethod
settags[!](siteinds, commoninds, M1::MPO, M2::MPS, args...; kwargs...)
settags[!](siteinds, commoninds, M1::MPO, M2::MPO, args...; kwargs...)

Apply settags to the site indices that are shared by M1 and M2.

Returns new MPSs/MPOs. The ITensors of the MPSs/MPOs will be a view of the storage of the original ITensors.

source
ITensors.settagsMethod
settags[!](siteinds, uniqueinds, M1::MPO, M2::MPS, args...; kwargs...)

Apply settags to the site indices of M1 that are not shared with M2. Returns new MPSs/MPOs.

The ITensors of the MPSs/MPOs will be a view of the storage of the original ITensors.

source

Operations

ITensors.dagMethod
dag[!](M::MPS, args...; kwargs...)
dag[!](M::MPO, args...; kwargs...)

Apply dag to all ITensors of an MPS/MPO, returning a new MPS/MPO.

The ITensors of the MPS/MPO will be a view of the storage of the original ITensors. Alternatively apply the function in-place.

source
NDTensors.denseMethod
dense(::MPS/MPO)

Given an MPS (or MPO), return a new MPS (or MPO) having called dense on each ITensor to convert each tensor to use dense storage and remove any QN or other sparse structure information, if it is not dense already.

source
ITensors.movesiteMethod
movesite(::Union{MPS, MPO}, n1n2::Pair{Int, Int})

Create a new MPS/MPO where the site at n1 is moved to n2, for a pair n1n2 = n1 => n2.

This is done with a series a pairwise swaps, and can introduce a lot of entanglement into your state, so use with caution.

source
ITensors.orthogonalize!Function
orthogonalize!(M::MPS, j::Int; kwargs...)
orthogonalize(M::MPS, j::Int; kwargs...)

orthogonalize!(M::MPO, j::Int; kwargs...)
orthogonalize(M::MPO, j::Int; kwargs...)

Move the orthogonality center of the MPS to site j. No observable property of the MPS will be changed, and no truncation of the bond indices is performed. Afterward, tensors 1,2,...,j-1 will be left-orthogonal and tensors j+1,j+2,...,N will be right-orthogonal.

Either modify in-place with orthogonalize! or out-of-place with orthogonalize.

source
ITensors.replacebond!Method
replacebond!(M::MPS, b::Int, phi::ITensor; kwargs...)

Factorize the ITensor phi and replace the ITensors b and b+1 of MPS M with the factors. Choose the orthogonality with ortho="left"/"right".

source
ITensors.sampleMethod
sample(m::MPS)

Given a normalized MPS m with orthocenter(m)==1, returns a Vector{Int} of length(m) corresponding to one sample of the probability distribution defined by squaring the components of the tensor that the MPS represents

source
ITensors.sample!Method
sample!(m::MPS)

Given a normalized MPS m, returns a Vector{Int} of length(m) corresponding to one sample of the probability distribution defined by squaring the components of the tensor that the MPS represents. If the MPS does not have an orthogonality center, orthogonalize!(m,1) will be called before computing the sample.

source
ITensors.sampleMethod
sample(M::MPO)

Given a normalized MPO M, returns a Vector{Int} of length(M) corresponding to one sample of the probability distribution defined by the MPO, treating the MPO as a density matrix.

The MPO M should have an (approximately) positive spectrum.

source
NDTensors.truncate!Function
truncate!(M::MPS; kwargs...)
truncate!(M::MPO; kwargs...)

Perform a truncation of all bonds of an MPS/MPO, using the truncation parameters (cutoff,maxdim, etc.) provided as keyword arguments.

source

Gate evolution

ITensors.productMethod
product(As::Vector{<:ITensor}, M::Union{MPS, MPO}; <keyword arguments>)
apply([...])

Apply the ITensors As to the MPS or MPO M, treating them as gates or matrices from pairs of prime or unprimed indices.

Examples

Apply one-site gates to an MPS:

N = 3

ITensors.op(::OpName"σx", ::SiteType"S=1/2", s::Index) =
  2*op("Sx", s)

ITensors.op(::OpName"σz", ::SiteType"S=1/2", s::Index) =
  2*op("Sz", s)

# Make the operator list.
os = [("σx", n) for n in 1:N]
append!(os, [("σz", n) for n in 1:N])

@show os

s = siteinds("S=1/2", N)
gates = ops(os, s)

# Starting state |↑↑↑⟩
ψ0 = productMPS(s, "↑")

# Apply the gates.
ψ = apply(gates, ψ0; cutoff = 1e-15)

# Test against exact (full) wavefunction
prodψ = apply(gates, prod(ψ0))
@show prod(ψ) ≈ prodψ

# The result is:
# σz₃ σz₂ σz₁ σx₃ σx₂ σx₁ |↑↑↑⟩ = -|↓↓↓⟩
@show inner(ψ, productMPS(s, "↓")) == -1

Apply nonlocal two-site gates and one-site gates to an MPS:

# 2-site gate
function ITensors.op(::OpName"CX", ::SiteType"S=1/2", s1::Index, s2::Index)
  mat = [1 0 0 0
         0 1 0 0
         0 0 0 1
         0 0 1 0]
  return itensor(mat, s2', s1', s2, s1)
end

os = [("CX", 1, 3), ("σz", 3)]

@show os

# Start with the state |↓↑↑⟩
ψ0 = productMPS(s, n -> n == 1 ? "↓" : "↑")

# The result is:
# σz₃ CX₁₃ |↓↑↑⟩ = -|↓↑↓⟩
ψ = apply(ops(os, s), ψ0; cutoff = 1e-15)
@show inner(ψ, productMPS(s, n -> n == 1 || n == 3 ? "↓" : "↑")) == -1

Perform TEBD-like time evolution:

# Define the nearest neighbor term `S⋅S` for the Heisenberg model
function ITensors.op(::OpName"expS⋅S", ::SiteType"S=1/2",
                     s1::Index, s2::Index; τ::Number)
  O = 0.5 * op("S+", s1) * op("S-", s2) +
      0.5 * op("S-", s1) * op("S+", s2) +
            op("Sz", s1) * op("Sz", s2)
  return exp(τ * O)
end

τ = -0.1im
os = [("expS⋅S", (1, 2), (τ = τ,)),
      ("expS⋅S", (2, 3), (τ = τ,))]
ψ0 = productMPS(s, n -> n == 1 ? "↓" : "↑")
expτH = ops(os, s)
ψτ = apply(expτH, ψ0)
source

Algebra Operations

LinearAlgebra.dotMethod
dot(A::MPS, B::MPS; make_inds_match = true)
inner(A::MPS, B::MPS; make_inds_match = true)

dot(A::MPO, B::MPO)
inner(A::MPO, B::MPO)

Compute the inner product <A|B>. If A and B are MPOs, computes the Frobenius inner product.

If make_inds_match = true, the function attempts to make the site indices match before contracting (so for example, the inputs can have different site indices, as long as they have the same dimensions or QN blocks).

For now, make_inds_match is only supported for MPSs.

See also logdot/loginner.

source
ITensors.logdotMethod
logdot(A::MPS, B::MPS; make_inds_match = true)
loginner(A::MPS, B::MPS; make_inds_match = true)

logdot(A::MPO, B::MPO)
loginner(A::MPO, B::MPO)

Compute the logarithm of the inner product <A|B>. If A and B are MPOs, computes the logarithm of the Frobenius inner product.

This is useful for larger MPS/MPO, where in the limit of large numbers of sites the inner product can diverge or approach zero.

If make_inds_match = true, the function attempts to make the site indices match before contracting (so for example, the inputs can have different site indices, as long as they have the same dimensions or QN blocks).

For now, make_inds_match is only supported for MPSs.

source
ITensors.lognormMethod
lognorm(A::MPS)
lognorm(A::MPO)

Compute the logarithm of the norm of the MPS or MPO.

This is useful for larger MPS/MPO that are not gauged, where in the limit of large numbers of sites the norm can diverge or approach zero.

See also norm and loginner/logdot.

source
Base.:+Method
+(A::MPS/MPO...; kwargs...)
add(A::MPS/MPO...; kwargs...)

Add arbitrary numbers of MPS/MPO with each other, with some optional truncation.

A cutoff of 1e-15 is used by default, and in general users should set their own cutoff for their particular application.

In the future we will give an interface for returning the truncation error.

Examples

N = 10

s = siteinds("S=1/2", N; conserve_qns = true)

state = n -> isodd(n) ? "↑" : "↓"
ψ₁ = randomMPS(s, state, 2)
ψ₂ = randomMPS(s, state, 2)
ψ₃ = randomMPS(s, state, 2)

ψ = +(ψ₁, ψ₂; cutoff = 1e-8)

# Can use:
#
# ψ = ψ₁ + ψ₂
#
# but generally you want to set a custom `cutoff` and `maxdim`.

println()
@show inner(ψ, ψ)
@show inner(ψ₁, ψ₂) + inner(ψ₁, ψ₂) + inner(ψ₂, ψ₁) + inner(ψ₂, ψ₂)

# Computes ψ₁ + 2ψ₂
ψ = ψ₁ + 2ψ₂

println()
@show inner(ψ, ψ)
@show inner(ψ₁, ψ₁) + 2 * inner(ψ₁, ψ₂) + 2 * inner(ψ₂, ψ₁) + 4 * inner(ψ₂, ψ₂)

# Computes ψ₁ + 2ψ₂ + ψ₃
ψ = ψ₁ + 2ψ₂ + ψ₃

println()
@show inner(ψ, ψ)
@show inner(ψ₁, ψ₁) + 2 * inner(ψ₁, ψ₂) + inner(ψ₁, ψ₃) +
      2 * inner(ψ₂, ψ₁) + 4 * inner(ψ₂, ψ₂) + 2 * inner(ψ₂, ψ₃) +
      inner(ψ₃, ψ₁) + 2 * inner(ψ₃, ψ₂) + inner(ψ₃, ψ₃)
source
Base.:*Method
contract(::MPS, ::MPO; kwargs...)
*(::MPS, ::MPO; kwargs...)

contract(::MPO, ::MPS; kwargs...)
*(::MPO, ::MPS; kwargs...)

Contract the MPO with the MPS, returning an MPS with the unique site indices of the MPO.

Choose the method with the method keyword, for example "densitymatrix" and "naive".

source