Reference
ITensorNetworks.ITensorNetwork
— TypeITensorNetwork
ITensorNetworks.ProjTTN
— TypeProjTTN
ITensorNetworks.ProjTTNSum
— TypeProjTTNSum
ITensorNetworks.TreeTensorNetwork
— TypeTreeTensorNetwork{V} <: AbstractTreeTensorNetwork{V}
ITensorNetworks._DensityMartrixAlgGraph
— TypeThe struct stores data used in the density matrix algorithm. partition: The given tn partition outtree: the binary tree structure of the output ITensorNetwork root: root vertex of the bfstree for truncation innerindstosim: mapping each inner index of the tn represented by partition
to a sim index caches: all the cached density matrices
ITensorNetworks._DensityMatrixAlgCaches
— TypeThe struct contains cached density matrices and cached partial density matrices for each edge / set of edges in the tensor network.
Density matrix example: Consider a tensor network below, 1 /9 2 / /3 6 /| /4 5 7 8 / | | The density matrix for the edge NamedEdge(2, 3)
squares the subgraph with vertices 3, 4, 5 | 3 /| 4 5 | | 4 5 |/ 3 | The density matrix for the edge NamedEdge(5, 3)
squares the subgraph with vertices 1, 2, 3, 4, 6, 7, 8, 9 1 // 2 / // 3 6 9 /| /| 4 7 8 | | | | | 4 7 8 | |/ | / | 3 6 | | / | | / | 2 9 / |/ 1 The density matrix for the edge NamedEdge(4, 3)
squares the subgraph with vertices 1, 2, 3, 5, 6, 7, 8, 9 1 // 2 / // 3 6 9 /| /| 5 7 8 | | | | | 5 7 8 | |/ | / | 3 6 | | / | | / | 2 9 / |/ 1
Partial density matrix example: Consider a tensor network below, 1 /9 2 / /3 6 /| /4 5 7 8 / | | The partial density matrix for the Edge set Set([NamedEdge(2, 3), NamedEdge(5, 3)])
squares the subgraph with vertices 4, and contract with the tensor 3 | 3 / 4 - 4 - The partial density matrix for the Edge set Set([NamedEdge(4, 3), NamedEdge(5, 3)])
squares the subgraph with vertices 1, 2, 6, 7, 8, 9, and contract with the tensor 3 1 // 2 / // 3 6 9 /| /| 7 8 | | | | 7 8 | | / | 6 | / | | / | 2 9 / |/ 1 The density matrix for the Edge set Set([NamedEdge(4, 3), NamedEdge(2, 3)])
squares the subgraph with vertices 5. and contract with the tensor 3 | 3 / 5 - 5 -
ITensorNetworks._approx_itensornetwork_density_matrix!
— MethodApproximate a partition
into an output ITensorNetwork with the binary tree structure defined by out_tree
.
ITensorNetworks._approx_itensornetwork_ttn_svd!
— MethodApproximate a partition
into an output ITensorNetwork with the binary tree structure defined by out_tree
by first transforming the partition into a ttn, then truncating the ttn using a sequence of SVDs.
ITensorNetworks._binary_tree_structure
— MethodGiven a tn and outinds (a subset of noncommoninds of tn), get a DataGraph
with binary tree structure of outinds that will be used in the binary tree partition. If maximally_unbalanced=true, the binary tree will have a line/mps structure. The binary tree is recursively constructed from leaves to the root.
Example:
TODO
ITensorNetworks._contract_deltas
— MethodGiven an input tensor network tn
, remove redundent delta tensors in tn
and change inds accordingly to make the output tn
represent the same tensor network but with less delta tensors.
======== Example: julia> is = [Index(2, string(i)) for i in 1:6] julia> a = ITensor(is[1], is[2]) julia> b = ITensor(is[2], is[3]) julia> delta1 = delta(is[3], is[4]) julia> delta2 = delta(is[5], is[6]) julia> tn = ITensorNetwork([a, b, delta1, delta2]) julia> ITensorNetworks.contractdeltas(tn) ITensorNetwork{Int64} with 3 vertices: 3-element Vector{Int64}: 1 2 4
and 1 edge(s): 1 => 2
with vertex data: 3-element Dictionaries.Dictionary{Int64, Any} 1 │ ((dim=2|id=457|"1"), (dim=2|id=296|"2")) 2 │ ((dim=2|id=296|"2"), (dim=2|id=613|"4")) 4 │ ((dim=2|id=626|"6"), (dim=2|id=237|"5"))
ITensorNetworks._contract_deltas_ignore_leaf_partitions
— MethodGiven an input partition
, contract redundent delta tensors of non-leaf vertices in partition
without changing the tensor network value. root
is the root of the dfs_tree that defines the leaves. Note: for each vertex v
of partition
, the number of non-delta tensors in partition[v]
will not be changed. Note: only delta tensors of non-leaf vertices will be contracted. Note: this function assumes that all noncommoninds of the partition are in leaf partitions.
ITensorNetworks._delta_inds_disjointsets
— MethodGiven a list of delta tensors deltas
, return a DisjointSets
of all its indices such that each pair of indices adjacent to any delta tensor must be in the same disjoint set. If a disjoint set contains indices in rootinds
, then one of such indices in rootinds
must be the root of this set.
ITensorNetworks._densitymatrix_outinds_to_sim
— MethodReturns a dict that maps the partition's outinds that are adjacent to partition[root]
to siminds
ITensorNetworks._distance
— MethodSum of shortest path distances among all outinds.
ITensorNetworks._introot_union!
— MethodRewrite of the function DataStructures.root_union!(s::IntDisjointSet{T}, x::T, y::T) where {T<:Integer}
.
ITensorNetworks._maxweightoutinds_tn
— Methodcreate a tn with empty ITensors whose outinds weights are MAXWEIGHT The maxweighttn is constructed so that only commoninds of the tn will be considered in mincut.
ITensorNetworks._mincut
— MethodCalculate the mincut between two subsets of the uncontracted inds (sourceinds and terminalinds) of the input tn. Mincut of two inds list is defined as the mincut of two newly added vertices, each one neighboring to one inds subset.
ITensorNetworks._mincut_inds
— MethodFind a vector of indices within sourceindslist yielding the mincut of given tnpair. Args: tnpair: a pair of tns (tn1 => tn2), where tn2 is generated via _maxweightoutindstn(tn1) outtomaxweightind: a dict mapping each out ind in tn1 to out ind in tn2 sourceindslist: a list of vector of indices to be considered Note: For each sourceinds in sourceindslist, we consider its mincut within both tns (tn1, tn2) given in tnpair. The mincut in tn1 represents the rank upper bound when splitting sourceinds with other inds in outinds. The mincut in tn2 represents the rank upper bound when the weights of outinds are very large. The first mincut upper_bounds the number of non-zero singular values, while the second empirically reveals the singular value decay. We output the sourceinds where the first mincut value is the minimum, the secound mincut value is also the minimum under the condition that the first mincut is optimal, and the sourceinds have the lowest all-pair shortest path.
ITensorNetworks._mincut_partitions
— MethodCalculate the mincutpartitions between two subsets of the uncontracted inds (sourceinds and terminal_inds) of the input tn.
ITensorNetworks._mps_partition_inds_order
— MethodGiven a tn and outinds, returns a vector of indices representing MPS inds ordering.
ITensorNetworks._optcontract
— MethodContract of a vector of tensors, network
, with a contraction sequence generated via sa_bipartite
ITensorNetworks._partition
— MethodGiven an input tn and a rooted binary tree of indices, return a partition of tn with the same binary tree structure as indsbtree. Note: in the output partition, we add multiple delta tensors to the network so that the output graph is guaranteed to be the same binary tree as indsbtree. Note: in the output partition, we add multiple scalar tensors. These tensors are used to make the output partition connected, even if the input tn
is disconnected. Note: in the output partition, tensor vertex names will be changed. For a given input tensor with vertex name v
, its name in the output partition will be
(v, 1). Any delta tensor will have name
(v, 2), and any scalar tensor used to maintain the connectivity of the partition will have name
(v, 3). Note: for a given binary tree with n indices, the output partition will contain 2n-1 vertices, with each leaf vertex corresponding to a sub tn adjacent to one output index. Keeping these leaf vertices in the partition makes later
approxitensornetworkalgorithms more efficient. Note: name of vertices in the output partition are the same as the name of vertices in
indsbtree`.
ITensorNetworks._rem_leaf_vertices!
— MethodFor a given ITensorNetwork tn
and a root
vertex, remove leaf vertices in the directed tree with root root
without changing the tensor represented by tn. In particular, the tensor of each leaf vertex is contracted with the tensor of its parent vertex to keep the tensor unchanged.
ITensorNetworks._rem_vertex!
— MethodPerform truncation and remove root
vertex in the partition
and out_tree
of alg_graph
.
Example: Consider an alg_graph
whose
outtreeis shown below, 1 /9 2 / /3 6 /| /4 5 7 8 / | | when
root = 4, the output
outtreewill be 1 /9 2 / /3 6 /| /5 7 8 | | and the returned tensor
U` will be the projector at vertex 4 in the output tn.
ITensorNetworks._root_union!
— MethodRewrite of the function DataStructures.root_union!(s::DisjointSet{T}, x::T, y::T)
. The difference is that in the output of _root_union!
, x is guaranteed to be the root of y when setting left_root=true
, and y will be the root of x when setting left_root=false
. In DataStructures.root_union!
, the root value cannot be specified. A specified root is useful in functions such as _remove_deltas
, where when we union two indices into one disjointset, we want the index that is the outinds if the given tensor network to always be the root in the DisjointSets.
ITensorNetworks._sim
— MethodReplace the inds of partialdmtensor that are in keys of inds_to_siminds
to the corresponding value, and replace the inds that are in values of inds_to_siminds
to the corresponding key.
ITensorNetworks._update!
— MethodUpdate caches.e_to_dm[e]
and caches.es_to_pdm[es]
. caches: the caches of the density matrix algorithm. edge: the edge defining the density matrix children: the children vertices of dst(edge)
in the dfstree network: the tensor network at vertex dst(edge)
indsto_sim: a dict mapping inds to sim inds
ITensorNetworks.add
— MethodAdd two itensornetworks together by growing the bond dimension. The network structures need to be have the same vertex names, same site index on each vertex
ITensorNetworks.alternating_update
— Methodtdvp(Hs::Vector{MPO},init_state::MPS,t::Number; kwargs...)
tdvp(Hs::Vector{MPO},init_state::MPS,t::Number, sweeps::Sweeps; kwargs...)
Use the time dependent variational principle (TDVP) algorithm to compute exp(t*H)*init_state
using an efficient algorithm based on alternating optimization of the MPS tensors and local Krylov exponentiation of H.
This version of tdvp
accepts a representation of H as a Vector of MPOs, Hs = [H1,H2,H3,...] such that H is defined as H = H1+H2+H3+... Note that this sum of MPOs is not actually computed; rather the set of MPOs [H1,H2,H3,..] is efficiently looped over at each step of the algorithm when optimizing the MPS.
Returns:
state::MPS
- time-evolved MPS
ITensorNetworks.binary_tree_structure
— MethodGiven a tn
and outinds
(a subset of noncommoninds of tn
), outputs a directed binary tree DataGraph of outinds
defining the desired graph structure
ITensorNetworks.binary_tree_structure
— MethodOutputs a directed binary tree DataGraph defining the desired graph structure
ITensorNetworks.contract_approx
— FunctionApproximate a given ITensorNetwork tn
into an output ITensorNetwork with output_structure
. output_structure
outputs a directed binary tree DataGraph defining the desired graph structure.
ITensorNetworks.contract_approx
— MethodApproximate a binary_tree_partition
into an output ITensorNetwork with the same binary tree structure. root
is the root vertex of the pre-order depth-first-search traversal used to perform the truncations.
ITensorNetworks.contract_approx
— MethodApproximate a given ITensorNetwork tn
into an output ITensorNetwork with a binary tree structure. The binary tree structure is defined based on inds_btree
, which is a directed binary tree DataGraph of indices.
ITensorNetworks.contraction_sequence_to_digraph
— MethodTake a contraction sequence and return a directed graph.
ITensorNetworks.contraction_sequence_to_graph
— MethodTake a contractionsequence and return a graphical representation of it. The leaves of the graph represent the leaves of the sequence whilst the internalnodes of the graph define a tripartition of the graph and thus are named as an n = 3 element tuples, which each element specifying the keys involved. Edges connect parents/children within the contraction sequence.
ITensorNetworks.contraction_tree_leaf_bipartition
— MethodGet the vertex bi-partition that a given edge represents
ITensorNetworks.delta_network
— MethodRETURN A TENSOR NETWORK WITH COPY TENSORS ON EACH VERTEX. Note that passing a link_space will mean the indices of the resulting network don't match those of the input indsnetwork
ITensorNetworks.fidelity
— MethodCalculate the overlap of the gate acting on the previous p and q versus the new p and q in the presence of environments. This is the cost function that optimisepq will minimise
ITensorNetworks.findall_on_edges
— MethodFind all edges e
such that f(graph[e]) == true
ITensorNetworks.findall_on_vertices
— MethodFind all vertices v
such that f(graph[v]) == true
ITensorNetworks.findfirst_on_edges
— MethodFind the edge e
such that f(graph[e]) == true
ITensorNetworks.findfirst_on_vertices
— MethodFind the vertex v
such that f(graph[v]) == true
ITensorNetworks.gauge_error
— MethodFunction to measure the 'distance' of a state from the Vidal Gauge
ITensorNetworks.is_multi_edge
— MethodCheck if the edge of an itensornetwork has multiple indices
ITensorNetworks.optimise_p_q
— MethodDo Full Update Sweeping, Optimising the tensors p and q in the presence of the environments envs, Specifically this functions find the pcur and qcur which optimise envsgatepqdag(prime(pcur))*dag(prime(qcur))
ITensorNetworks.path_graph_structure
— MethodGiven a tn
and outinds
(a subset of noncommoninds of tn
), outputs a maximimally unbalanced directed binary tree DataGraph of outinds
defining the desired graph structure
ITensorNetworks.path_graph_structure
— MethodOutputs a maximimally unbalanced directed binary tree DataGraph defining the desired graph structure
ITensorNetworks.random_tensornetwork
— MethodBuild an ITensor network on a graph specified by the inds network s. Bonddim is given by linkspace and entries are randomized. The random distribution is based on the input argument distribution
.
ITensorNetworks.random_tensornetwork
— MethodBuild an ITensor network on a graph specified by the inds network s. Bonddim is given by linkspace and entries are randomised (normal distribution, mean 0 std 1)
ITensorNetworks.tdvp
— Methodtdvp(operator::TTN, t::Number, init_state::TTN; kwargs...)
Use the time dependent variational principle (TDVP) algorithm to approximately compute exp(operator*t)*init_state
using an efficient algorithm based on alternating optimization of the state tensors and local Krylov exponentiation of operator. The time parameter t
can be a real or complex number.
Returns:
state
- time-evolved state
Optional keyword arguments:
time_step::Number = t
- time step to use when evolving the state. Smaller time steps generally give more accurate results but can make the algorithm take more computational time to run.nsteps::Integer
- evolve by the requested total timet
by performingnsteps
of the TDVP algorithm. More steps can result in more accurate results but require more computational time to run. (Note that only one of thetime_step
ornsteps
parameters can be provided, not both.)outputlevel::Int = 1
- larger outputlevel values resulting in printing more information and 0 means no outputobserver
- object implementing the Observer interface which can perform measurements and stop earlywrite_when_maxdim_exceeds::Int
- when the allowed maxdim exceeds this value, begin saving tensors to disk to free memory in large calculations
ITensorNetworks.to_vec
— Methodto_vec(x)
Transform x
into a Vector
. Returns the vector and a closure which inverts the transformation.
Modeled after FiniteDifferences.to_vec
:
https://github.com/JuliaDiff/FiniteDifferences.jl/blob/main/src/to_vec.jl
ITensorNetworks.ttn
— Methodttn(os::OpSum, sites::IndsNetwork{<:Index}; kwargs...)
ttn(eltype::Type{<:Number}, os::OpSum, sites::IndsNetwork{<:Index}; kwargs...)
Convert an OpSum object os
to a TreeTensorNetwork, with indices given by sites
.
ITensorNetworks.ttn_svd
— Methodttn_svd(os::OpSum, sites::IndsNetwork, root_vertex, kwargs...)
Construct a TreeTensorNetwork from a symbolic OpSum representation of a Hamiltonian, compressing shared interaction channels.
ITensorNetworks.update
— MethodDo parallel updates between groups of edges of all message tensors Currently we send the full message tensor data struct to update for each edge_group. But really we only need the mts relevant to that group.
ITensorNetworks.update
— MethodDo a sequential update of the message tensors on edges
ITensorNetworks.update
— MethodMore generic interface for update, with default params
ITensorNetworks.update_factors
— MethodUpdate the tensornetwork inside the cache
ITensorNetworks.updated_message
— MethodCompute message tensor as product of incoming mts and local state
ITensorNetworks.vidal_gauge_isometry
— MethodFunction to construct the 'isometry' of a state in the Vidal Gauge on the given edge
ITensorNetworks.vidalitensornetwork_preserve_cache
— MethodUse an ITensorNetwork ψ, its bond tensors and belief propagation cache to put ψ into the vidal gauge, return the bond tensors and updated_ψ.
ITensors.apply
— MethodApply() function for an ITN in the Vidal Gauge. Hence the bond tensors are required. Gate does not necessarily need to be passed. Can supply an edge to do an identity update instead. Uses Simple Update procedure assuming gate is two-site
ITensors.apply
— MethodOverload of ITensors.apply
.
KrylovKit.linsolve
— Functionlinsolve(
A::ITensorNetworks.AbstractTreeTensorNetwork,
b::ITensorNetworks.AbstractTreeTensorNetwork,
x₀::ITensorNetworks.AbstractTreeTensorNetwork;
...
)
linsolve(
A::ITensorNetworks.AbstractTreeTensorNetwork,
b::ITensorNetworks.AbstractTreeTensorNetwork,
x₀::ITensorNetworks.AbstractTreeTensorNetwork,
a₀::Number;
...
)
linsolve(
A::ITensorNetworks.AbstractTreeTensorNetwork,
b::ITensorNetworks.AbstractTreeTensorNetwork,
x₀::ITensorNetworks.AbstractTreeTensorNetwork,
a₀::Number,
a₁::Number;
updater,
nsites,
nsweeps,
updater_kwargs,
kwargs...
)
Compute a solution x to the linear system:
(a₀ + a₁ * A)*x = b
using starting guess x₀. Leaving a₀, a₁ set to their default values solves the system A*x = b.
To adjust the balance between accuracy of solution and speed of the algorithm, it is recommed to first try adjusting the solver_tol
keyword argument descibed below.
Keyword arguments:
ishermitian::Bool=false
- should set to true if the MPO A is Hermitiansolver_krylovdim::Int=30
- max number of Krylov vectors to build on each solver iterationsolver_maxiter::Int=100
- max number outer iterations (restarts) to do in the solver stepsolver_tol::Float64=1E-14
- tolerance or error goal of the solver
Overload of KrylovKit.linsolve
.
NDTensors.contract
— MethodOverload of ITensors.contract
.
ITensorNetworks.ModelNetworks.ising_network
— MethodBUILD Z OF CLASSICAL ISING MODEL ON A GIVEN GRAPH AT INVERSE TEMP BETA H = -\sum{(v,v') \in edges}\sigma^{z}{v}\sigma^{z}_{v'} OPTIONAL ARGUMENT: h: EXTERNAL MAGNETIC FIELD szverts: A LIST OF VERTICES OVER WHICH TO APPLY A SZ. THE RESULTANT NETWORK CAN THEN BE CONTRACTED AND DIVIDED BY THE ACTUAL PARTITION FUNCTION TO GET THAT OBSERVABLE INDSNETWORK IS ASSUMED TO BE BUILT FROM A GRAPH (NO SITE INDS) AND OF LINK SPACE 2
ITensorNetworks.ModelNetworks.ising_network_state
— MethodBuild the wavefunction whose norm is equal to Z of the classical ising model s needs to have site indices in this case!
ITensorNetworks.ModelHamiltonians.heisenberg
— MethodRandom field J1-J2 Heisenberg model on a general graph
ITensorNetworks.ModelHamiltonians.heisenberg
— MethodRandom field J1-J2 Heisenberg model on a chain of length N
ITensorNetworks.ModelHamiltonians.hubbard
— Methodt-t' Hubbard Model g,i,v
ITensorNetworks.ModelHamiltonians.ising
— MethodNext-to-nearest-neighbor Ising model (ZZX) on a general graph
ITensorNetworks.ModelHamiltonians.ising
— MethodNext-to-nearest-neighbor Ising model (ZZX) on a chain of length N