Basics of Tree Tensor Networks#

Overview#

Tree tensor network state (TTNS) and tree tensor network operator (TTNO) are generalizations of MPS and MPO. In this notebook, we will show how to define tree tensor network states in Renormalizer. Overall the interfaces for tree tensor networks are very similar to that of the matrix product states. The biggest difference is how we define the model (tree structure), which we will describe in detail in the following.

Defining a Tree#

The tree structure in Renormalizer is specified through BasisTree, i.e., a tree of basis sets. This is similar to the list of basis sets for the definition of MPS.

The tree node of BasisTree is TreeNodeBasis. In the following we’ll construct a number of TreeNodeBasiss and connect them as a tree.

[1]:
from renormalizer import BasisHalfSpin
from renormalizer.tn import TreeNodeBasis, BasisTree
2024-08-24 09:07:24,399[INFO] Use NumPy as backend
2024-08-24 09:07:24,399[INFO] numpy random seed is 9012
2024-08-24 09:07:24,400[INFO] random seed is 1092
2024-08-24 09:07:24,411[INFO] Git Commit Hash: 78516d8ee697e5cfb88a1cf2c55ad3299ff3b640
2024-08-24 09:07:24,412[INFO] use 64 bits
[2]:
# some spin DOF for future usage
spins = [BasisHalfSpin(i) for i in range(5)]
spins
[2]:
[BasisHalfSpin(dof: 0, nbas: 2),
 BasisHalfSpin(dof: 1, nbas: 2),
 BasisHalfSpin(dof: 2, nbas: 2),
 BasisHalfSpin(dof: 3, nbas: 2),
 BasisHalfSpin(dof: 4, nbas: 2)]
[3]:
# construct tree nodes. One node for one basis set
nodes = [TreeNodeBasis(basis) for basis in spins]

Tree nodes are can connected together by the add_child function. We first pick a root and then add two children to the root

[4]:
root = nodes[0]
root.add_child(nodes[1:3])
[4]:
TreeNodeBasis(BasisHalfSpin(dof: 0, nbas: 2))

We can check they are already connected through .children and .parent attributes. Note that a node can have multiple children but can only have one parent.

[5]:
root.children
[5]:
[TreeNodeBasis(BasisHalfSpin(dof: 1, nbas: 2)),
 TreeNodeBasis(BasisHalfSpin(dof: 2, nbas: 2))]
[6]:
nodes[1].parent
[6]:
TreeNodeBasis(BasisHalfSpin(dof: 0, nbas: 2))

We then connect the rest of the nodes to one of the children

[7]:
root.children[0].add_child(nodes[3:])
[7]:
TreeNodeBasis(BasisHalfSpin(dof: 1, nbas: 2))

Now all of the nodes are connected, we can construct the tree by feeding the root node to the BasisTree class.

We can visualize the tree structure through the print function, which shows the degrees of freedom in the tree.

[8]:
tree = BasisTree(root)
tree.print()

             ┌[(3,)]
      ┌[(1,)]┤
      │      └[(4,)]
[(0,)]┤
      └[(2,)]

All nodes in a tree, including TTNS and TTNO, can be accessed via a list, according to pre-order traversal.

[9]:
tree.node_list
[9]:
[TreeNodeBasis(BasisHalfSpin(dof: 0, nbas: 2)),
 TreeNodeBasis(BasisHalfSpin(dof: 1, nbas: 2)),
 TreeNodeBasis(BasisHalfSpin(dof: 3, nbas: 2)),
 TreeNodeBasis(BasisHalfSpin(dof: 4, nbas: 2)),
 TreeNodeBasis(BasisHalfSpin(dof: 2, nbas: 2))]

The BasisTree object can be considered as immutable. If a new tree is desired by modifying the existing tree, the BasisTree instance should be recreated with the corresponding root node.

Advanced Tree Structure#

In the example above, the number of BasisSets in each node is one. Renormalizer permits any number of BasisSets for each node.

First of all, if no BasisSet is provided, the node will become purely virtual and will not associated with any physical degree of freedom. A dummy basis set is attached to the node for programming convenience. The dimension of the Hilbert space of this dummy degree of freedom is 1 and the only allowed operator is identity.

[10]:
node1 = TreeNodeBasis()
node1
[10]:
TreeNodeBasis(BasisDummy(dof: ('Virtual DOF', 0), nbas: 1))

In addition, each tree node can have multiple associated BasisSet.

[11]:
node2 = TreeNodeBasis([BasisHalfSpin(i) for i in range(3)])
node2
[11]:
TreeNodeBasis(BasisHalfSpin(dof: 0, nbas: 2), BasisHalfSpin(dof: 1, nbas: 2), BasisHalfSpin(dof: 2, nbas: 2))

Since building the tree structure can be tedious, Rernomalizer has a number of builtin tree structure that can be used out of box. The most simple example is the MPS (linear) tree strucure.

[12]:
BasisTree.linear([BasisHalfSpin(i) for i in range(5)]).print()

[(0,)]─[(1,)]─[(2,)]─[(3,)]─[(4,)]

Another example is simple binary tree

[13]:
BasisTree.binary([BasisHalfSpin(i) for i in range(7)]).print()

             ┌[(3,)]
      ┌[(1,)]┤
      │      └[(4,)]
[(0,)]┤
      │      ┌[(5,)]
      └[(2,)]┤
             └[(6,)]

Another frequent tree structure is the “MCTDH-style” tree, or hierarchical Tucker format.

The feature of this type of tree is that all physical degrees of freedom are attached to the leaf nodes. Also, each leaf node typically has more than one physical degrees of freedom.

The following shows this type of tree with a binary structure and a ternary structure

[14]:
BasisTree.binary_mctdh([BasisHalfSpin(i) for i in range(8)]).print()

                                                   ┌[(0,), (1,)]
                         ┌[(('MCTDH virtual', 1),)]┤
                         │                         └[(2,), (3,)]
[(('MCTDH virtual', 0),)]┤
                         │                         ┌[(4,), (5,)]
                         └[(('MCTDH virtual', 2),)]┤
                                                   └[(6,), (7,)]
[15]:
BasisTree.ternary_mctdh([BasisHalfSpin(i) for i in range(9)]).print()

                         ┌[(0,), (1,), (2,)]
[(('MCTDH virtual', 0),)]┼[(3,), (4,), (5,)]
                         └[(6,), (7,), (8,)]

Note that the internal nodes (non-leaf nodes) is not associated with any physical degree of freedom and a virtual degree of freedom is attached.

TTNS and TTNO#

Having defined the tree structure using BasisTree. We can now construct TTNS and TTNO accordingly. The interfaces are very similar to that of MPS, but BasisTree take over the roles of Model. With this modification it is more convenient to perform TTNO/TTNS algebra when the idea of the Hamiltonian is not well defined, such as in the case of quantum circuit simulation.

We recommend reading previous tutorials if the following code seems to be too complicated.

TTNS#

We first construct TTNS based on a complete binary spin tree with 7 sites in total.

[16]:
from renormalizer.tn import TTNS, TTNO
basis = BasisTree.binary([BasisHalfSpin(i) for i in range(7)])
[17]:
ttns = TTNS.random(basis, qntot=0, m_max=5)
ttns
[17]:
<class 'renormalizer.tn.tree.TTNS'> with 7 nodes

Each node in TTNS is a TreeNodeTensor whose core data is the numerical tensor. The tensor indices are arranged as follows: indices to the children, physical indices (if any) and the index to the parent.

[18]:
ttns.node_list
[18]:
[TreeNodeTensor((5, 5, 2, 1),float64),
 TreeNodeTensor((2, 2, 2, 5),float64),
 TreeNodeTensor((2, 2),float64),
 TreeNodeTensor((2, 2),float64),
 TreeNodeTensor((2, 2, 2, 5),float64),
 TreeNodeTensor((2, 2),float64),
 TreeNodeTensor((2, 2),float64)]
[19]:
ttns.node_list[2].tensor
[19]:
array([[-0.99496648, -0.10020826],
       [-0.10020826,  0.99496648]])
[20]:
ttns.calc_1site_rdm()
2024-08-24 09:07:24,602[DEBUG] # of operator terms: 1
2024-08-24 09:07:24,605[DEBUG] Input operator terms: 1
2024-08-24 09:07:24,608[DEBUG] After combination of the same terms: 1
[20]:
{0: array([[ 0.36219816, -0.05108033],
        [-0.05108033,  0.63780184]]),
 1: array([[0.52219779, 0.04792124],
        [0.04792124, 0.47780221]]),
 2: array([[0.51527382, 0.01204695],
        [0.01204695, 0.48472618]]),
 3: array([[ 0.45914641, -0.03186007],
        [-0.03186007,  0.54085359]]),
 4: array([[0.53935311, 0.17021523],
        [0.17021523, 0.46064689]]),
 5: array([[ 0.39225728, -0.05994727],
        [-0.05994727,  0.60774272]]),
 6: array([[ 0.67911799, -0.05841049],
        [-0.05841049,  0.32088201]])}

TTNS by default will construct a hartree product state given the basis tree

[21]:
TTNS(basis).bond_dims
[21]:
[1, 1, 1, 1, 1, 1, 1]

The bond dimension is the dimension of the parent index for each node in pre-order traversal.

TTNO#

We next turn to TTNO. We will boroow the testing Holstein Hamiltonian for MPS to illustrate TTNO construction.

[22]:
from renormalizer.tests.parameter import holstein_model
[23]:
holstein_model.basis
[23]:
[BasisSimpleElectron(dof: 0, nbas: 2, qn: [[0], [1]]),
 BasisSHO(dof: (0, 0), x0: 0.0, omega: 0.0004852952677876329, nbas: 4),
 BasisSHO(dof: (0, 1), x0: 0.0, omega: 0.007087607302666907, nbas: 4),
 BasisSimpleElectron(dof: 1, nbas: 2, qn: [[0], [1]]),
 BasisSHO(dof: (1, 0), x0: 0.0, omega: 0.0004852952677876329, nbas: 4),
 BasisSHO(dof: (1, 1), x0: 0.0, omega: 0.007087607302666907, nbas: 4),
 BasisSimpleElectron(dof: 2, nbas: 2, qn: [[0], [1]]),
 BasisSHO(dof: (2, 0), x0: 0.0, omega: 0.0004852952677876329, nbas: 4),
 BasisSHO(dof: (2, 1), x0: 0.0, omega: 0.007087607302666907, nbas: 4)]
[24]:
holstein_model.ham_terms
[24]:
[Op('a^\\dagger a', [0, 0], 0.10016074648883302, [[1], [-1]]),
 Op('a^\\dagger a', [0, 1], -0.003674932217565499, [[1], [-1]]),
 Op('a^\\dagger a', [0, 2], -0.007349864435130998, [[1], [-1]]),
 Op('a^\\dagger a', [1, 0], -0.003674932217565499, [[1], [-1]]),
 Op('a^\\dagger a', [1, 1], 0.10016074648883302, [[1], [-1]]),
 Op('a^\\dagger a', [1, 2], -0.011024796652696497, [[1], [-1]]),
 Op('a^\\dagger a', [2, 0], -0.007349864435130998, [[1], [-1]]),
 Op('a^\\dagger a', [2, 1], -0.011024796652696497, [[1], [-1]]),
 Op('a^\\dagger a', [2, 2], 0.10016074648883302, [[1], [-1]]),
 Op('p^2', [(0, 0)], 0.5),
 Op('x^2', [(0, 0)], 1.1775574846853516e-07),
 Op('p^2', [(0, 1)], 0.5),
 Op('x^2', [(0, 1)], 2.5117088638408635e-05),
 Op('p^2', [(1, 0)], 0.5),
 Op('x^2', [(1, 0)], 1.1775574846853516e-07),
 Op('p^2', [(1, 1)], 0.5),
 Op('x^2', [(1, 1)], 2.5117088638408635e-05),
 Op('p^2', [(2, 0)], 0.5),
 Op('x^2', [(2, 0)], 1.1775574846853516e-07),
 Op('p^2', [(2, 1)], 0.5),
 Op('x^2', [(2, 1)], 2.5117088638408635e-05),
 Op('a^\\dagger a x', [0, 0, (0, 0)], -7.097609983192488e-06, [[1], [-1], [0]]),
 Op('a^\\dagger a x', [0, 0, (0, 1)], -0.00044069941383179025, [[1], [-1], [0]]),
 Op('a^\\dagger a x', [1, 1, (1, 0)], -7.097609983192488e-06, [[1], [-1], [0]]),
 Op('a^\\dagger a x', [1, 1, (1, 1)], -0.00044069941383179025, [[1], [-1], [0]]),
 Op('a^\\dagger a x', [2, 2, (2, 0)], -7.097609983192488e-06, [[1], [-1], [0]]),
 Op('a^\\dagger a x', [2, 2, (2, 1)], -0.00044069941383179025, [[1], [-1], [0]])]

The basis tree is construct manually to reflect the interaction in the Hamiltonian

[25]:
nodes = []
for i in range(3):
    node1 = TreeNodeBasis(holstein_model.basis[3*i])
    node2 = TreeNodeBasis(holstein_model.basis[3*i+1:3*i+3])
    node1.add_child(node2)
    nodes.append(node1)
nodes[1].add_child([nodes[0], nodes[2]])

basis = BasisTree(nodes[1])
basis.print()

      ┌[((1, 0),), ((1, 1),)]
[(1,)]┼[(0,)]─[((0, 0),), ((0, 1),)]
      └[(2,)]─[((2, 0),), ((2, 1),)]
[26]:
ttno = TTNO(basis, holstein_model.ham_terms)
2024-08-24 09:07:24,697[DEBUG] # of operator terms: 27
2024-08-24 09:07:24,697[DEBUG] Input operator terms: 27
2024-08-24 09:07:24,698[DEBUG] After combination of the same terms: 27
[27]:
ttno.bond_dims
[27]:
[1, 3, 4, 3, 4, 3]
[28]:
ttno.node_list
[28]:
[TreeNodeTensor((3, 4, 4, 2, 2, 1),float64),
 TreeNodeTensor((4, 4, 4, 4, 3),float64),
 TreeNodeTensor((3, 2, 2, 4),float64),
 TreeNodeTensor((4, 4, 4, 4, 3),float64),
 TreeNodeTensor((3, 2, 2, 4),float64),
 TreeNodeTensor((4, 4, 4, 4, 3),float64)]

TTNO/TTNS manipulation#

We demonstrate TTNO/TTNS manipulation using the Holstein model

[29]:
ttns = TTNS(basis)
[30]:
ttns + ttns
[30]:
<class 'renormalizer.tn.tree.TTNS'> with 6 nodes
[31]:
ttno @ ttns
[31]:
<class 'renormalizer.tn.tree.TTNS'> with 6 nodes
[32]:
ttns.expectation(ttno)
[32]:
0.011359353855681813