A TreeTensorNetwork (alias TTN) is an ITensorNetwork whose underlying graph is a
tree (no cycles). This additional structure enables exact, efficient canonical gauges
via QR decomposition — a key ingredient in variational algorithms such as DMRG and TDVP.
A TreeTensorNetwork carries an extra piece of metadata: the ortho_region, which
records which vertices currently form the orthogonality center of the network. Algorithms
update this field as the gauge changes.
MPS (matrix product states) are the special case of a TreeTensorNetwork on a
1D path graph. The [mps](@ref ITensorNetworks.mps) constructor enforces this topology
and provides a convenient interface for 1D calculations.
using Graphs: vertices
using ITensorNetworks: ITensorNetwork, TreeTensorNetwork, mps, ortho_region, orthogonalize,
siteinds, ttn
using ITensors: ITensors
using LinearAlgebra: norm
using NamedGraphs.NamedGraphGenerators: named_comb_tree
# Comb-tree TTN (a popular tree topology for 2D-like systems)
g = named_comb_tree((3, 2))
sites = siteinds("S=1/2", g)
psi = ttn(sites) # zero-initialised
psi = ttn(v -> "Up", sites) # product state
# 1D MPS
s1d = siteinds("S=1/2", 6)
mps_state = mps(v -> "Up", s1d) # product MPS
ITensorNetworks.ttn
ITensorNetworks.mps
The TreeTensorNetwork struct wraps an ITensorNetwork and records the current
orthogonality region. Use the TreeTensorNetwork constructor to convert a plain
ITensorNetwork with tree topology into a TTN, and ITensorNetwork to strip the
gauge metadata when you need a plain network again.
itn = ITensorNetwork(psi) # TTN → ITensorNetwork
psi = TreeTensorNetwork(itn) # ITensorNetwork → TTN
ITensorNetworks.TreeTensorNetwork
ITensorNetworks.ITensorNetwork(::ITensorNetworks.TreeTensorNetwork)
A dense tensor can be decomposed into a TTN by successive QR/SVD factorisations along the
tree edges. Truncation parameters (e.g. cutoff, maxdim) are forwarded to the
factorisation step.
g = named_comb_tree((3, 1))
sites = siteinds("S=1/2", g)
A = ITensors.random_itensor(only(sites[(1, 1)]), only(sites[(2, 1)]), only(sites[(3, 1)]))
ttn_A = ttn(A, sites)
ITensorNetworks.ttn(::ITensors.ITensor, ::ITensorNetworks.IndsNetwork)
One of the most powerful features of tree tensor networks is the ability to bring the
network into an orthogonal gauge in linear time. When the network is in a gauge
centered on vertex v, all tensors away from v are isometric with respect to the bond
pointing toward v. This makes computing local observables, inner products, and
eigenvalue problems numerically efficient and stable.
The current orthogonality center is tracked by the ortho_region field.
v = collect(vertices(psi))[1]
v1 = collect(vertices(psi))[1]
v2 = collect(vertices(psi))[2]
vs = [v]
psi = orthogonalize(psi, v) # QR-sweep to put ortho center at vertex v
psi = orthogonalize(psi, [v1, v2]) # two-site center (for nsites=2 sweeps)
ortho_region(psi) # query current ortho region (returns an index set)
ITensorNetworks.orthogonalize
ITensorNetworks.ortho_region
After algorithms that grow the bond dimension (e.g. addition, subspace expansion), use
truncate to recompress the network. For TreeTensorNetwork there are two forms:
- Whole-network recompression (TTN-specific):
truncate(ttn; kwargs...)sweeps from the leaves to the root, orthogonalising and truncating every bond in sequence. This is the preferred form after addition or DMRG expansion. - Single-bond truncation (available for any
ITensorNetwork):truncate(tn, edge; kwargs...)truncates one bond by SVD — see the ITensor Networks page.
psi = truncate(psi; cutoff = 1e-10, maxdim = 50)
The sweep-based form orthogonalises each bond before truncating it, so truncation errors
are controlled. All keyword arguments accepted by ITensors.svd (e.g. cutoff, maxdim,
mindim) are forwarded.
Base.truncate(::ITensorNetworks.AbstractTreeTensorNetwork)
Two TTNs with the same graph and site indices can be summed. The result has bond
dimension equal to the sum of the two inputs, and can be recompressed with truncate.
psi1, psi2 = psi, psi
psi3 = psi1 + psi2 # or add(psi1, psi2)
psi3 = truncate(psi3; cutoff = 1e-10, maxdim = 50)
2 * psi # scalar multiplication
psi / norm(psi) # manual normalisation
ITensorNetworks.add(::ITensorNetworks.AbstractTreeTensorNetwork, ::ITensorNetworks.AbstractTreeTensorNetwork)
Base.:+(::ITensorNetworks.AbstractTreeTensorNetwork...)