{
"cells": [
{
"attachments": {},
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"(tensor-network-design)=\n",
"\n",
"# Design\n",
"\n",
"The section contains some notes on the design of quimb's tensor network\n",
"functionality that should be helpful for a) understanding how certain functions\n",
"work better, and b) implementing more advanced functions or contributing to quimb."
]
},
{
"cell_type": "code",
"execution_count": 1,
"metadata": {},
"outputs": [],
"source": [
"%config InlineBackend.figure_formats = ['svg']\n",
"\n",
"# the main tensor and tensor network functionality\n",
"import quimb.tensor as qtn\n",
"\n",
"# backend agnostic array functions\n",
"from autoray import do"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"## `Tensor`\n",
"\n",
"The {class}`~quimb.tensor.tensor_core.Tensor` object, unsurprisingly,\n",
"is the core object in `quimb.tensor`. It has three main attributes.\n",
"\n",
"### `Tensor.data`\n",
"\n",
"This is **any n-dimensional array looking object** (i.e. something with a\n",
"`.shape` attribute). Whenever array functions are performed on a\n",
"tensors `.data`, [autoray](https://github.com/jcmgray/autoray)\n",
"is used to dispatch to the correct backend functions, so `quimb`\n",
"itself does no explicit conversion of array backends."
]
},
{
"cell_type": "code",
"execution_count": 2,
"metadata": {},
"outputs": [],
"source": [
"data = do('random.uniform', size=(3, 4, 5), like='cupy')"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"Having said that, if the supplied `data` has no `.shape` attribute\n",
"{func}`numpy.asarray` will be called on it.\n",
"\n",
"### `Tensor.inds`\n",
"\n",
"This is an **ordered tuple of index names** - one for each of the dimensions\n",
"of `.data`. Functions that would normally require 'axes' to be specified\n",
"instead can just use names, and these are propagated through all operations.\n",
"For example, operations such as transposing the underyling array can just\n",
"be done via reordering the `inds` of the tensor."
]
},
{
"cell_type": "code",
"execution_count": 3,
"metadata": {},
"outputs": [],
"source": [
"inds = ('a', 'b', 'c')"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"Moreover, the names of indices explicitly define the geometry of bonds\n",
"(edges) in tensor networks at this local tensor level - there are no extra\n",
"'index' objects and this information is thus completely local.\n",
"\n",
"### `Tensor.tags`\n",
"\n",
"These are an **ordered set** (quimb has a quick `dict` based implementation\n",
"\\- {class}`~quimb.utils.oset`) of **an arbitrary number of extra\n",
"identifiers** for the tensor. The main use of these is that *within a\n",
"tensor network* you can use any number of simultenous labelling schemes to\n",
"identify your tensors with. You probably don't need to use any if just using\n",
"raw tensors."
]
},
{
"cell_type": "code",
"execution_count": 4,
"metadata": {},
"outputs": [],
"source": [
"tags = qtn.oset(['hello', 'world'])"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
":::{note}\n",
"{class}`~quimb.utils.oset` is used in many places across `quimb.tensor`\n",
"**in order to make all operations deterministic**, unlike {class}`set`. You can\n",
"supply any (ideally itself ordered) iterable to `tags` and it will\n",
"be converted for you.\n",
":::\n",
"\n",
"## Creating & Contracting Tensors\n",
"\n",
"We can use our `data`, `inds` and (optional) `tags` to create a tensor:"
]
},
{
"cell_type": "code",
"execution_count": 5,
"metadata": {},
"outputs": [
{
"data": {
"text/html": [
"Tensor (shape=(3 , 4 , 5 ), inds=[a , b , c ], tags={hello , world }), backend=cupy , dtype=float64 , data=array([[[0.97641966, 0.6964992 , 0.14297803, 0.4075414 , 0.38340721],\n",
" [0.11856711, 0.89633028, 0.59211389, 0.547883 , 0.5934418 ],\n",
" [0.75794605, 0.89411517, 0.3868147 , 0.0090452 , 0.29680125],\n",
" [0.91689399, 0.84223811, 0.46408656, 0.58087305, 0.89949236]],\n",
"\n",
" [[0.74613888, 0.77655533, 0.05940234, 0.23631638, 0.82518739],\n",
" [0.57981722, 0.54969852, 0.79932289, 0.77402098, 0.86315628],\n",
" [0.46177098, 0.23686848, 0.35011062, 0.5614713 , 0.13329136],\n",
" [0.082935 , 0.51894518, 0.50928339, 0.95217395, 0.67710415]],\n",
"\n",
" [[0.22488959, 0.53207111, 0.74554138, 0.95702639, 0.78858658],\n",
" [0.51136781, 0.55735051, 0.77407297, 0.20138574, 0.72253493],\n",
" [0.73572433, 0.78192978, 0.6545592 , 0.38579206, 0.46297037],\n",
" [0.52019734, 0.28735351, 0.57257375, 0.69800421, 0.21275464]]]) "
],
"text/plain": [
"Tensor(shape=(3, 4, 5), inds=('a', 'b', 'c'), tags=oset(['hello', 'world']))"
]
},
"execution_count": 5,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"T = qtn.Tensor(data=data, inds=inds, tags=tags)\n",
"T"
]
},
{
"cell_type": "code",
"execution_count": 6,
"metadata": {},
"outputs": [
{
"data": {
"image/svg+xml": [
"\n",
"\n",
"\n",
" \n",
" \n",
" \n",
" \n",
" 2023-11-29T16:12:21.820935 \n",
" image/svg+xml \n",
" \n",
" \n",
" Matplotlib v3.7.2, https://matplotlib.org/ \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n"
],
"text/plain": [
""
]
},
"metadata": {
"needs_background": "light"
},
"output_type": "display_data"
}
],
"source": [
"T.draw(color='hello')"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"Just on their own (i.e. not within a tensor network) just the labelled dimension\n",
"aspect of tensors is already very convenient.\n",
"Here's {meth}`~quimb.tensor.tensor_core.Tensor.fuse`,\n",
"for example, which would usually require some non-trivial combination of axis finding,\n",
"permuting, reshaping etc:"
]
},
{
"cell_type": "code",
"execution_count": 7,
"metadata": {},
"outputs": [
{
"data": {
"text/html": [
"Tensor (shape=(15 , 4 ), inds=[ac , b ], tags={hello , world }), backend=cupy , dtype=float64 , data=array([[0.97641966, 0.11856711, 0.75794605, 0.91689399],\n",
" [0.6964992 , 0.89633028, 0.89411517, 0.84223811],\n",
" [0.14297803, 0.59211389, 0.3868147 , 0.46408656],\n",
" [0.4075414 , 0.547883 , 0.0090452 , 0.58087305],\n",
" [0.38340721, 0.5934418 , 0.29680125, 0.89949236],\n",
" [0.74613888, 0.57981722, 0.46177098, 0.082935 ],\n",
" [0.77655533, 0.54969852, 0.23686848, 0.51894518],\n",
" [0.05940234, 0.79932289, 0.35011062, 0.50928339],\n",
" [0.23631638, 0.77402098, 0.5614713 , 0.95217395],\n",
" [0.82518739, 0.86315628, 0.13329136, 0.67710415],\n",
" [0.22488959, 0.51136781, 0.73572433, 0.52019734],\n",
" [0.53207111, 0.55735051, 0.78192978, 0.28735351],\n",
" [0.74554138, 0.77407297, 0.6545592 , 0.57257375],\n",
" [0.95702639, 0.20138574, 0.38579206, 0.69800421],\n",
" [0.78858658, 0.72253493, 0.46297037, 0.21275464]]) "
],
"text/plain": [
"Tensor(shape=(15, 4), inds=('ac', 'b'), tags=oset(['hello', 'world']))"
]
},
"execution_count": 7,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"T1 = T.fuse({'ac': ('a', 'c')})\n",
"T1"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"More importantly, the indices of a tensor also exactly define the subscripts\n",
"as would appear in a\n",
"[implicit Einstein summation](https://en.wikipedia.org/wiki/Einstein_notation),\n",
"such that if we create the following tensor:"
]
},
{
"cell_type": "code",
"execution_count": 8,
"metadata": {},
"outputs": [
{
"data": {
"text/html": [
"Tensor (shape=(6 , 15 ), inds=[d , ac ], tags={T2 }), backend=cupy , dtype=float64 , data=array([[0.61133688, 0.19963178, 0.54512192, 0.45136319, 0.95200169,\n",
" 0.68349322, 0.72054855, 0.35749849, 0.01521056, 0.76230481,\n",
" 0.95791951, 0.55055866, 0.84917076, 0.38681014, 0.795198 ],\n",
" [0.71160786, 0.90332059, 0.76197671, 0.53927703, 0.42098982,\n",
" 0.82002418, 0.90061322, 0.48990883, 0.75245645, 0.54850696,\n",
" 0.15722231, 0.14761825, 0.38377117, 0.71195033, 0.22025585],\n",
" [0.22271323, 0.6002882 , 0.29158765, 0.47496975, 0.89326597,\n",
" 0.94999516, 0.37679517, 0.8857649 , 0.61252652, 0.66570427,\n",
" 0.43343891, 0.10384877, 0.94062542, 0.72933167, 0.29090906],\n",
" [0.77569595, 0.91219689, 0.67441759, 0.0982548 , 0.86256351,\n",
" 0.43845473, 0.18648777, 0.40573815, 0.11182709, 0.03758457,\n",
" 0.66733917, 0.22062413, 0.37580556, 0.02836813, 0.46582853],\n",
" [0.62867367, 0.81370598, 0.10292322, 0.9020288 , 0.36236472,\n",
" 0.29828873, 0.07717501, 0.36800433, 0.2745711 , 0.24406048,\n",
" 0.73428318, 0.35486658, 0.56554947, 0.53758826, 0.45742711],\n",
" [0.41337104, 0.88415806, 0.50730996, 0.77998229, 0.97972492,\n",
" 0.4472324 , 0.92755498, 0.04104936, 0.48683115, 0.97383488,\n",
" 0.97269154, 0.96569825, 0.08422208, 0.4742093 , 0.76365892]]) "
],
"text/plain": [
"Tensor(shape=(6, 15), inds=('d', 'ac'), tags=oset(['T2']))"
]
},
"execution_count": 8,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"T2 = qtn.Tensor(\n",
" # use whatever backend function T1's data is using\n",
" data=do('random.uniform', size=(6, 15), like=T1.data),\n",
" # 'd' is new, 'ac' is shared with T1\n",
" inds=['d', 'ac'],\n",
" tags=['T2'])\n",
"T2"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"The fact that both tensors have a `'ac'` index means that when combined\n",
"the two corresponding dimensions are either implicitly or explicitly summed\n",
"over (contracted) like so:\n",
"\n",
"```{math}\n",
"T_{3}[b, d] = \\sum_{ac} T_1[ac, b] T_2[d, ac]\n",
"```"
]
},
{
"cell_type": "code",
"execution_count": 9,
"metadata": {},
"outputs": [
{
"data": {
"image/svg+xml": [
"\n",
"\n",
"\n",
" \n",
" \n",
" \n",
" \n",
" 2023-11-29T16:12:21.959380 \n",
" image/svg+xml \n",
" \n",
" \n",
" Matplotlib v3.7.2, https://matplotlib.org/ \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n"
],
"text/plain": [
""
]
},
"metadata": {
"needs_background": "light"
},
"output_type": "display_data"
}
],
"source": [
"T1.add_tag('T1')\n",
"(T1 | T2).draw(color=['T1', 'T2'])"
]
},
{
"cell_type": "code",
"execution_count": 10,
"metadata": {},
"outputs": [
{
"data": {
"text/html": [
"Tensor (shape=(4 , 6 ), inds=[b , d ], tags={hello , world , T1 , T2 }), backend=cupy , dtype=float64 , data=array([[5.22498485, 5.03995703, 4.74283749, 3.35558005, 3.92328952,\n",
" 5.61195934],\n",
" [5.24083022, 5.07076672, 5.40682846, 3.77210843, 3.95661664,\n",
" 5.91813226],\n",
" [4.06943567, 3.88964851, 3.79027758, 3.51389366, 3.43819655,\n",
" 4.50403753],\n",
" [4.82537913, 5.21932489, 5.04876845, 3.84178214, 4.11865728,\n",
" 5.68091703]]) "
],
"text/plain": [
"Tensor(shape=(4, 6), inds=('b', 'd'), tags=oset(['hello', 'world', 'T1', 'T2']))"
]
},
"execution_count": 10,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"T3 = T1 @ T2 \n",
"# == qtn.tensor_contract(T1, T2)\n",
"# == T1.contract(T2)\n",
"T3"
]
},
{
"cell_type": "code",
"execution_count": 11,
"metadata": {},
"outputs": [
{
"data": {
"image/svg+xml": [
"\n",
"\n",
"\n",
" \n",
" \n",
" \n",
" \n",
" 2023-11-29T16:12:22.072938 \n",
" image/svg+xml \n",
" \n",
" \n",
" Matplotlib v3.7.2, https://matplotlib.org/ \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n"
],
"text/plain": [
""
]
},
"metadata": {
"needs_background": "light"
},
"output_type": "display_data"
}
],
"source": [
"T3.draw(color=['T1'])"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"The `@` operator eagerly contracts two tensors according to this\n",
"summation (via [`qtn.tensor_contract`](quimb.tensor.tensor_core.tensor_contract)),\n",
"thereby removing the `'ac'` index. As you can see, contracting\n",
"two tensors also merges their tags.\n",
"\n",
":::{hint}\n",
"In `quimb`, **the names of the indices entirely define the geometry\n",
"of the tensor network**. If you want two tensor dimensions to be contracted\n",
"just name them the same thing. Alternatively use\n",
"[`qtn.connect`](quimb.tensor.tensor_core.connect) to name them for you. The actual\n",
"name is not so important, in as much as its unique and shared, you can\n",
"generate new names automatically using\n",
"[`qtn.rand_uuid`](quimb.tensor.tensor_core.rand_uuid).\n",
"\n",
"The fact that the specific name of internal indices has no effect on the\n",
"overall oject represented is made use of by\n",
"[`TensorNetwork`](quimb.tensor.tensor_core.TensorNetwork) -\n",
"see - {ref}`indexmangling`.\n",
":::\n",
"\n",
"## Copying & Modifying Tensors\n",
"\n",
"You can copy a tensor with [`Tensor.copy`](quimb.tensor.tensor_core.Tensor.copy), (note\n",
"by default this **doesn't** copy the `.data` - see below), and update it with\n",
"[`Tensor.modify`](quimb.tensor.tensor_core.Tensor.modify). Routing all changes to tensors\n",
"through [`Tensor.modify`](quimb.tensor.tensor_core.Tensor.modify) allows various checks, and\n",
"for tensors to automatically inform any tensor networks they belong to, for example.\n",
"\n",
"An example of a typical tensor method that demonstrates the above is the following:\n"
]
},
{
"cell_type": "code",
"execution_count": 12,
"metadata": {},
"outputs": [
{
"data": {
"text/html": [
"Tensor (shape=(15 ), inds=[ac ], tags={hello , world , T1 }), backend=cupy , dtype=float64 , data=array([2.76982681, 3.32918276, 1.58599319, 1.54534265, 2.17314261,\n",
" 1.87066208, 2.08206752, 1.71811924, 2.5239826 , 2.49873919,\n",
" 1.99217906, 2.15870491, 2.7467473 , 2.2422084 , 2.18684652]) "
],
"text/plain": [
"Tensor(shape=(15,), inds=('ac',), tags=oset(['hello', 'world', 'T1']))"
]
},
"execution_count": 12,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"def sum_reduce(self, ind, inplace=False):\n",
" \"\"\"Sum the array axis corresponding to ``ind``.\n",
" \"\"\"\n",
" # easily enables copy vs inplace version\n",
" t = self if inplace else self.copy()\n",
" \n",
" # get the axis `ind` corresponds to\n",
" axis = t.inds.index(ind)\n",
" \n",
" # use autoray.do to be backend agnostic\n",
" new_data = do('sum', t.data, axis=axis)\n",
" new_inds = t.inds[:axis] + t.inds[axis + 1:]\n",
" \n",
" # update tensor using modify\n",
" t.modify(data=new_data, inds=new_inds)\n",
" \n",
" # return `t` even if inplace to allow method chaining\n",
" return t\n",
"\n",
"sum_reduce(T1, 'b')"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"A convention that `quimb` adopts is for many common methods to have inplace\n",
"versions specifically named with a trailing underscore, like:\n",
"\n",
"- [`Tensor.sum_reduce`](quimb.tensor.tensor_core.Tensor.sum_reduce)\n",
"- [`Tensor.sum_reduce_`](quimb.tensor.tensor_core.Tensor.sum_reduce_)"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
":::{warning}\n",
"All array functions in `quimb` as assumed to be **pure**, i.e. return\n",
"a new array and not modify the original.\n",
"This allows tensor copies by default *to not copy their* `.data`, which saves\n",
"memory and is more efficient with backends that e.g. construct an explicit\n",
"computational graph. You should only inplace modify the actual arrays if\n",
"you know what you are doing, or after calling `.copy(deep=True)` first.\n",
":::"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"## `TensorNetwork`\n",
"\n",
"As described above, there is already an implicit tensor network\n",
"geometry introduced by simply matching up tensor index names -\n",
"and indeed you could just keep track of the tensors yourself.\n",
"\n",
"However, its generally useful to collect tensors into a\n",
"{class}`~quimb.tensor.tensor_core.TensorNetwork` object,\n",
"which enables $\\mathscr{O}(1)$ access to any tensors based on\n",
"their `inds` or `tags`, as well as encapsulating all sorts of\n",
"functionality the depends on the network relations of the tensors."
]
},
{
"cell_type": "code",
"execution_count": 13,
"metadata": {},
"outputs": [
{
"data": {
"text/html": [
"TensorNetwork (tensors=3, indices=6)Tensor (shape=(2 , 3 , 3 ), inds=[a , x , y ], tags={A }), backend=numpy , dtype=float64 , data=array([[[ 1.74058538, -0.09727624, -1.51295378],\n",
" [ 2.89531464, -0.57855413, 1.20852059],\n",
" [-0.76578221, 0.90925054, 0.4563162 ]],\n",
"\n",
" [[ 1.22858738, -0.41212638, -0.83124623],\n",
" [ 0.81559135, -0.81969156, -1.11280435],\n",
" [-1.715578 , 0.46736194, 1.31020563]]])Tensor (shape=(2 , 3 , 3 ), inds=[b , y , z ], tags={B }), backend=numpy , dtype=float64 , data=array([[[-0.77826738, -1.24821725, -0.62410868],\n",
" [-1.14399318, -0.16202123, -0.15981997],\n",
" [-1.31423348, -0.89974453, -0.46182192]],\n",
"\n",
" [[-0.71201596, 1.62508318, -1.30336751],\n",
" [ 0.17445857, 1.12671226, -0.58820474],\n",
" [-1.34346682, 1.01526495, -0.05090115]]])Tensor (shape=(2 , 3 , 3 ), inds=[c , z , x ], tags={C }), backend=numpy , dtype=float64 , data=array([[[ 0.14340117, -0.99054769, 0.30215627],\n",
" [ 0.01621206, 1.88997099, -0.47993354],\n",
" [-1.2574597 , -1.6274538 , 0.11786352]],\n",
"\n",
" [[-0.72554233, -0.39697657, -0.70370705],\n",
" [-0.70023463, 0.6874081 , 0.61294832],\n",
" [ 0.92082505, 0.36431587, 0.3457875 ]]]) "
],
"text/plain": [
"TensorNetwork(tensors=3, indices=6)"
]
},
"execution_count": 13,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"ta = qtn.rand_tensor((2, 3, 3), ['a', 'x', 'y'], tags='A')\n",
"tb = qtn.rand_tensor((2, 3, 3), ['b', 'y', 'z'], tags='B')\n",
"tc = qtn.rand_tensor((2, 3, 3), ['c', 'z', 'x'], tags='C')\n",
"\n",
"tn = qtn.TensorNetwork([ta, tb, tc])\n",
"tn"
]
},
{
"cell_type": "code",
"execution_count": 14,
"metadata": {},
"outputs": [
{
"data": {
"image/svg+xml": [
"\n",
"\n",
"\n",
" \n",
" \n",
" \n",
" \n",
" 2023-11-29T16:12:22.219877 \n",
" image/svg+xml \n",
" \n",
" \n",
" Matplotlib v3.7.2, https://matplotlib.org/ \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n"
],
"text/plain": [
""
]
},
"metadata": {
"needs_background": "light"
},
"output_type": "display_data"
}
],
"source": [
"tn.draw(color=['A', 'B', 'C'])"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
":::{note}\n",
"By default TNs take copies of the input tensors,\n",
"(but again, not the actual data array).\n",
"This can be turned off and it is\n",
"possible e.g. to create multiple TNs that view\n",
"the same tensors - see {ref}`virtualtns`.\n",
":::\n",
"\n",
"Tensor networks also have three key attributes...\n",
"\n",
"### `TensorNetwork.tensor_map`\n",
"\n",
"The key storage of a {class}`~quimb.tensor.tensor_core.TensorNetwork`\n",
"is the `tensor_map` attribute. This is a mapping of 'tids' - unique\n",
"integers reprensenting nodes in the network - to the actual tensor at\n",
"that node, which need not be unique."
]
},
{
"cell_type": "code",
"execution_count": 15,
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"{0: Tensor(shape=(2, 3, 3), inds=('a', 'x', 'y'), tags=oset(['A'])),\n",
" 1: Tensor(shape=(2, 3, 3), inds=('b', 'y', 'z'), tags=oset(['B'])),\n",
" 2: Tensor(shape=(2, 3, 3), inds=('c', 'z', 'x'), tags=oset(['C']))}"
]
},
"execution_count": 15,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"tn.tensor_map"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"This means that the same tensor can be in the same or different tensor networks multiple times.\n",
"\n",
":::{hint}\n",
"`'tids'` are essentially vertex labels in the (hyper) graph of the TN.\n",
":::\n",
"\n",
"### `TensorNetwork.ind_map`\n",
"\n",
"This is a mapping of every index in the tensor network to the ordered\n",
"set of 'tids' of the tensors which have that index."
]
},
{
"cell_type": "code",
"execution_count": 16,
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"{'a': oset([0]),\n",
" 'x': oset([0, 2]),\n",
" 'y': oset([0, 1]),\n",
" 'b': oset([1]),\n",
" 'z': oset([1, 2]),\n",
" 'c': oset([2])}"
]
},
"execution_count": 16,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"tn.ind_map"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"For example if we wanted all tensors with the ``'z'`` index, we could call:"
]
},
{
"cell_type": "code",
"execution_count": 17,
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"[Tensor(shape=(2, 3, 3), inds=('b', 'y', 'z'), tags=oset(['B'])),\n",
" Tensor(shape=(2, 3, 3), inds=('c', 'z', 'x'), tags=oset(['C']))]"
]
},
"execution_count": 17,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"[tn.tensor_map[tid] for tid in tn.ind_map['z']]"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
":::{hint}\n",
"`ind_map[ix]` is essentially the set of (hyper) graph vertices the\n",
"(hyper) edge `ix` is incident to.\n",
":::\n",
"\n",
"### `TensorNetwork.tag_map`\n",
"\n",
"Similarly, this is a mapping of every tag in the tensor network to the ordered\n",
"set of 'tids' of the tensors which have that tag."
]
},
{
"cell_type": "code",
"execution_count": 18,
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"{'A': oset([0]), 'B': oset([1]), 'C': oset([2])}"
]
},
"execution_count": 18,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"tn.tag_map"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"If we wanted all tensors tagged with ``'B'`` we could call:"
]
},
{
"cell_type": "code",
"execution_count": 19,
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"[Tensor(shape=(2, 3, 3), inds=('b', 'y', 'z'), tags=oset(['B']))]"
]
},
"execution_count": 19,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"[tn.tensor_map[tid] for tid in tn.tag_map['B']]"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"## Selecting Tensors Based on Tags - the `which` kwarg\n",
"\n",
"`tags` are the main high level method for labelling and\n",
"accessing tensors. Many functions which take a `tags` and\n",
"`which` arg, e.g. [`TensorNetwork.select`](quimb.tensor.tensor_core.TensorNetwork.select)\n",
"internally call\n",
"{meth}`~quimb.tensor.tensor_core.TensorNetwork._get_tids_from_tags`.\n",
"When supplied to this, the `which` kwarg takes four options:\n",
"\n",
"- `'all'` return tids for tensors which match **all** the tags\n",
"- `'!all'` return tids for tensors which *don't* match **all** the tags\n",
"- `'any'` return tids for tensors which match **any** of the tags\n",
"- `'!any'` return tids for tensors which *don't* match **any** of the tags"
]
},
{
"cell_type": "code",
"execution_count": 20,
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"oset([0, 1])"
]
},
"execution_count": 20,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"tn._get_tids_from_tags(['A', 'B'], which='any')"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"## Copying & Modifying Tensor Networks\n",
"\n",
"Tensor networks also have a\n",
"[`TensorNetwork.copy`](quimb.tensor.tensor_core.TensorNetwork.copy) method which\n",
"is used heavily. As with creation, by default this copies both\n",
"the TN and the tensors it contains (but again, not their data).\n",
"\n",
"In terms of modifying tensor networks, `quimb` takes the general\n",
"approach that you **directly modify tensors in a TN,\n",
"which then automatically tell the TN, and possibly other\n",
"TNs, about any relevant changes**."
]
},
{
"cell_type": "code",
"execution_count": 21,
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"{'x': oset([0, 2]),\n",
" 'y': oset([0, 1]),\n",
" 'b': oset([1, 0]),\n",
" 'z': oset([1, 2]),\n",
" 'c': oset([2])}"
]
},
"execution_count": 21,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"# get the tensor identified by tag 'A'\n",
"t = tn['A'] \n",
"\n",
"# inplace rename an index\n",
"t.reindex_({'a': 'b'})\n",
"\n",
"# changes are automatically reflected in ``tn``\n",
"tn.ind_map"
]
},
{
"cell_type": "code",
"execution_count": 22,
"metadata": {},
"outputs": [
{
"data": {
"image/svg+xml": [
"\n",
"\n",
"\n",
" \n",
" \n",
" \n",
" \n",
" 2023-11-29T16:12:22.388794 \n",
" image/svg+xml \n",
" \n",
" \n",
" Matplotlib v3.7.2, https://matplotlib.org/ \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n"
],
"text/plain": [
""
]
},
"metadata": {
"needs_background": "light"
},
"output_type": "display_data"
}
],
"source": [
"# 'b' should now be traced over\n",
"tn.draw(color=['A', 'B', 'C'])"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"This functionality is enabled by tensors having 'owners' - see {ref}`tensorowners`.\n",
"\n",
"The following function demonstrates the anatomy of a typical\n",
"tensor network function. It performs the (probably not useful task) of:\n",
"\n",
"1. Getting all the tensors corresponding to a particular `tag`\n",
"2. Getting all the tensors neighboring to these\n",
"3. Performing a backend-agnostic `tanh` on the data of all these tensors"
]
},
{
"cell_type": "code",
"execution_count": 23,
"metadata": {},
"outputs": [],
"source": [
"def tanh_neighbors(self, tag, inplace=False):\n",
" \"\"\"Example function that demonstrates the anatomy of a \n",
" typical tensor network method, using ``tensor_map``, \n",
" ``ind_map``, ``tag_map``, and ``modify``.\n",
" \"\"\"\n",
" # easily handle both copy and non copy methods\n",
" tn = self if inplace else self.copy()\n",
" \n",
" # get the unique tensor identifiers from our tag\n",
" tids = tn.tag_map[tag]\n",
" \n",
" # collect neighbors in here -> oset for determinism\n",
" neighbors = qtn.oset()\n",
" \n",
" # for every tid in our tagged region\n",
" for tid in tids:\n",
" \n",
" # get the tensor\n",
" t = tn.tensor_map[tid]\n",
" \n",
" # for each of its inds\n",
" for ix in t.inds:\n",
" \n",
" # add all tensors with that index to neighbors\n",
" neighbors |= tn.ind_map[ix]\n",
" \n",
" # now apply tanh to our expanded region\n",
" for tid in neighbors:\n",
" # always use modify and autoray to update tensors\n",
" t = tn.tensor_map[tid]\n",
" t.modify(data=do('tanh', t.data))\n",
" \n",
" # return tn even if self (i.e. inplace) for method chaining\n",
" return tn\n",
"\n",
"qtn.TensorNetwork.tanh_neighbors = tanh_neighbors\n",
"\n",
"# create trailing underscore inplace version\n",
"import functools\n",
"qtn.TensorNetwork.tanh_neighbors_ = \\\n",
" functools.partialmethod(tanh_neighbors, inplace=True)"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"In use:"
]
},
{
"cell_type": "code",
"execution_count": 24,
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"1.7405853798277557"
]
},
"execution_count": 24,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"tn['A'].data.item(0)"
]
},
{
"cell_type": "code",
"execution_count": 25,
"metadata": {},
"outputs": [
{
"data": {
"text/html": [
"TensorNetwork (tensors=3, indices=5)Tensor (shape=(2 , 3 , 3 ), inds=[b , x , y ], tags={A }), backend=numpy , dtype=float64 , data=array([[[ 0.94029449, -0.09697056, -0.90746183],\n",
" [ 0.99390651, -0.52161375, 0.83623517],\n",
" [-0.64447019, 0.72077234, 0.42707706]],\n",
"\n",
" [[ 0.84216909, -0.39027667, -0.6811446 ],\n",
" [ 0.67266318, -0.67490196, -0.80505145],\n",
" [-0.93732835, 0.43606535, 0.86432743]]])Tensor (shape=(2 , 3 , 3 ), inds=[b , y , z ], tags={B }), backend=numpy , dtype=float64 , data=array([[[-0.6517111 , -0.84778297, -0.55398225],\n",
" [-0.81575435, -0.16061823, -0.158473 ],\n",
" [-0.86534269, -0.71617345, -0.43156796]],\n",
"\n",
" [[-0.61193943, 0.92535818, -0.86258756],\n",
" [ 0.17270992, 0.80989104, -0.52860321],\n",
" [-0.87250242, 0.76793088, -0.05085724]]])Tensor (shape=(2 , 3 , 3 ), inds=[c , z , x ], tags={C }), backend=numpy , dtype=float64 , data=array([[[ 0.14242623, -0.75759576, 0.29328465],\n",
" [ 0.01621064, 0.95537059, -0.44619038],\n",
" [-0.85036226, -0.92569812, 0.11732076]],\n",
"\n",
" [[-0.62033062, -0.37735903, -0.60671552],\n",
" [-0.60451669, 0.59631434, 0.54619917],\n",
" [ 0.72628749, 0.34900989, 0.33263441]]]) "
],
"text/plain": [
"TensorNetwork(tensors=3, indices=5)"
]
},
"execution_count": 25,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"tn.tanh_neighbors_('B')"
]
},
{
"cell_type": "code",
"execution_count": 26,
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"0.9402944918012002"
]
},
"execution_count": 26,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"# since A is connected to B, its data should have been tanh'd\n",
"tn['A'].data.item(0)"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"## Combining TNs\n",
"\n",
"(indexmangling)=\n",
"\n",
"### Index Mangling\n",
"\n",
"While it is very convenient generally to have the tensor network geometry\n",
"defined 'locally' just by `Tensor.inds`, **there is one main drawback**.\n",
"\n",
"1. Imagine we have two tensor networks, `tn1` and `tn2`,\n",
" both with the same outer indices - e.g. `('k0', 'k1', 'k2', 'k3', ...).\n",
"2. When we combine these these networks we expect the new TN to automatically\n",
" represent the overlap of these networks - i.e. with the `'k{}'` indices\n",
" contracted.\n",
"3. **However**, if `tn2` has the same internal indices as `tn1` (e.g. it\n",
" actually is `tn1`, or is a copy of or is derived from `tn1` etc.),\n",
" then these indices will now clash and appear four times in the new TN.\n",
"\n",
"The solution `quimb` adopts is that when you combine two or more tensor\n",
"networks, **inner indices in the latter will be mangled** - 'inner' being defined\n",
"as those appearing $\\geq 2$ times. This is fairly natural since renaming\n",
"internal indices has no effect on the overall TN object represented, but it does\n",
"mean **you should only rely on outer index names being preserved**.\n",
"\n",
":::{hint}\n",
"In general in `quimb`, you should keep track of the external indices of\n",
"a TN and the tags describing the internal tensor structure. If and when you\n",
"need explicit index names you can retrieve them from the tensors with, e.g.:\n",
"\n",
"- [`qtn.bonds`](quimb.tensor.tensor_core.bonds)\n",
"- [`qtn.group_inds`](quimb.tensor.tensor_core.group_inds)\n",
"- [`Tensor.filter_bonds`](quimb.tensor.tensor_core.Tensor.filter_bonds)\n",
":::\n",
"\n",
"(virtualtns)=\n",
"\n",
"### Virtual TNs\n",
"\n",
"When you create a {class}`~quimb.tensor.tensor_core.TensorNetwork`, from some\n",
"collection of tensors and/or other tensor networks, the tensors are added via:\n",
"\n",
"- [`TensorNetwork.add_tensor`](quimb.tensor.tensor_core.TensorNetwork.add_tensor)\n",
"- [`TensorNetwork.add_tensor_network`](quimb.tensor.tensor_core.TensorNetwork.add_tensor_network)\n",
"\n",
"By default the tensors (**but not the data - see above**), are copied, so that\n",
"they only appear in the new TN. This behavior corresponds to the\n",
"`virtual=False` option and overloaded operators:\n",
"\n",
"- `&`\n",
"- `&=`\n",
"\n",
"used to combine tensors and tensor networks."
]
},
{
"cell_type": "code",
"execution_count": 27,
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"False"
]
},
"execution_count": 27,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"tx = qtn.rand_tensor([3, 4], ['a', 'b'], tags='X')\n",
"ty = qtn.rand_tensor([4, 5], ['b', 'c'], tags='Y')\n",
"\n",
"tn = tx & ty\n",
"# == qtn.TensorNetwork([tx, ty], virtual=False)\n",
"\n",
"tx is tn['X']"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"Any changes to `tx` *won't* affect `tn`.\n",
"\n",
"If you want to tensor networks to 'view' existing tensors, either to\n",
"have the tensors appear in multiple networks, or simply because you\n",
"know a copy is not needed, you can use the `virtual=True` option.\n",
"This corresponds to the overloaded operators:\n",
"\n",
"- `|`\n",
"- `|=`\n",
"\n",
"used to combine tensors and tensor networks."
]
},
{
"cell_type": "code",
"execution_count": 28,
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"True"
]
},
"execution_count": 28,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"tn = tx | ty\n",
"# == qtn.TensorNetwork([tx, ty], virtual=True)\n",
"\n",
"tx is tn['X']"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"Now any changes to ``tx`` *will* affect ``tn``.\n",
"The ``virtual`` kwarg can also be supplied to ``copy``."
]
},
{
"cell_type": "code",
"execution_count": 29,
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"True"
]
},
"execution_count": 29,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"tx is tn.copy(virtual=True)['X']"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"(tensorowners)=\n",
"\n",
"## Tensor 'Owners'\n",
"\n",
"The piece of technology that enables tensors to tell, possibly mutiple,\n",
"tensor networks about changes to their `inds` and `tags` is **owners**.\n",
"You probably won't need to interact with this, but it might be useful\n",
"to understand, what happens behind the scenes.\n",
"\n",
"When a tensor is added to a tensor network, the tensor itself stores\n",
"a [weakref](https://docs.python.org/3/library/weakref.html) to the tensor\n",
"network. If the tensor's indices or tags are then changed using\n",
"[`Tensor.modify`](quimb.tensor.tensor_core.Tensor.modify), it can tell each tensor\n",
"network it has been added to to update their `ind_map` and `tag_map`\n",
"correctly.\n"
]
},
{
"cell_type": "code",
"execution_count": 30,
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"{8749283625669: (,\n",
" 0),\n",
" 8749283625672: (,\n",
" 0)}"
]
},
"execution_count": 30,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"tx = qtn.rand_tensor([3, 4], ['a', 'b'], tags='X')\n",
"ty = qtn.rand_tensor([4, 5], ['b', 'c'], tags='Y')\n",
"tz = qtn.rand_tensor([4, 5], ['b', 'c'], tags='Z')\n",
"\n",
"tn_xy = tx | ty\n",
"tn_xz = tx | tz\n",
"\n",
"tx.owners"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"A `weakref` is used so as not to prevent garbage collection of TNs, with stale\n",
"weakrefs being simply cleared out.\n",
"\n",
":::{note}\n",
"Since {class}`~quimb.tensor.tensor_core.TensorNetwork` by default copies\n",
"tensors before adding them, and a freshly copied tensor has no owners,\n",
"most tensors are usually only 'owned' by a single tensor network.\n",
":::"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"## Structured (1D, 2D, ...) Tensor Networks\n",
"\n",
"So far we have only talked about the design of generic tensor networks\n",
"of any geometry, and noted that the external indices of a TN, plus the\n",
"tags of the tensors inside it, are the things that essentially define it.\n",
"\n",
"Specific, structured, tensor networks such as MPS or PEPS are\n",
"implemented as subclasses of\n",
"{class}`~quimb.tensor.tensor_1d.TensorNetwork`, with extra properties\n",
"that describe, for example, a 'format' for how the physical indices and\n",
"tensors are tagged as a function of coordinate. Methods that rely on\n",
"assuming such structure can then be invoked.\n",
"\n",
"As an example consider\n",
"{class}`~quimb.tensor.tensor_1d.MatrixProductState`,\n",
"which is a subclass of\n",
"{class}`~quimb.tensor.tensor_1d.TensorNetwork`\n",
"but also has methods [mixed in](https://en.wikipedia.org/wiki/Mixin)\n",
"from the following:\n",
"\n",
"- {class}`~quimb.tensor.tensor_1d.TensorNetwork1D` -\n",
" methods that apply to any TN with a 1D tagging structure, for example,\n",
" [`TensorNetwork1D.site_tag`](quimb.tensor.tensor_1d.TensorNetwork1D.site_tag).\n",
"- {class}`~quimb.tensor.tensor_1d.TensorNetwork1DFlat` -\n",
" methods that apply to any 1D TN with a single tensor per site, for example,\n",
" [`TensorNetwork1DFlat.compress`](quimb.tensor.tensor_1d.TensorNetwork1DFlat.compress)\n",
"- {class}`~quimb.tensor.tensor_1d.TensorNetwork1DVector`.\n",
" methods that apply to any 1D TN with a single physical index per site, for example,\n",
" [`TensorNetwork1DVector.site_ind`](quimb.tensor.tensor_1d.TensorNetwork1DVector.site_ind).\n",
"\n",
"The MPS class then has the following extra properties that are required:"
]
},
{
"cell_type": "code",
"execution_count": 31,
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"('_site_tag_id', '_site_ind_id', 'cyclic', '_L')"
]
},
"execution_count": 31,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"qtn.MatrixProductState._EXTRA_PROPS"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"- `'_site_tag_id'` being the string formatter that converts site coordinate to site tag\n",
"- `'_site_ind_id'` being the string formatter that converts site coordinate to site ind\n",
"- `'cyclic'` describing whether the MPS has periodic or open boundary conditions\n",
"- `'_L'` describing the number of sites.\n",
"\n",
"These are usually generated automatically, from defaults or context:"
]
},
{
"cell_type": "code",
"execution_count": 32,
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"('I{}', 'k{}', True, 100)"
]
},
"execution_count": 32,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"mps = qtn.MPS_rand_state(L=100, bond_dim=20, cyclic=True)\n",
"mps.site_tag_id, mps.site_ind_id, mps.cyclic, mps.L"
]
},
{
"cell_type": "code",
"execution_count": 33,
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"('I{}', 'b{}', False, 4)"
]
},
"execution_count": 33,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"mps = qtn.MatrixProductState(\n",
" [\n",
" do('random.normal', size=(7, 2), like='numpy'),\n",
" do('random.normal', size=(7, 7, 2), like='numpy'),\n",
" do('random.normal', size=(7, 7, 2), like='numpy'),\n",
" do('random.normal', size=(7, 2), like='numpy'),\n",
" ], \n",
" site_ind_id='b{}'\n",
")\n",
"mps.site_tag_id, mps.site_ind_id, mps.cyclic, mps.L"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"### Converting Structured Tensor Networks\n",
"\n",
"Sometimes you might have created or modified *another*\n",
"{class}`~quimb.tensor.tensor_core.TensorNetwork` such that\n",
"you know it matches the structure of some more specific\n",
"geometry and now want to use relevant methods, for example,\n",
"if you have constructed a 2D tensor network and now want to use\n",
"[`TensorNetwork2D.contract_boundary`](quimb.tensor.tensor_2d.TensorNetwork2D.contract_boundary).\n",
"The following describes how to convert between generic and\n",
"structured tensor networks without going via the raw constructors,\n",
"that might be fiddly and inefficient.\n",
"\n",
"Imagine we have the following tensors and generic TN:"
]
},
{
"cell_type": "code",
"execution_count": 34,
"metadata": {},
"outputs": [],
"source": [
"t00 = qtn.rand_tensor([2, 4, 4], inds=['k0,0', '00-01', '00-10'], tags=['I0,0', 'X0', 'Y0'])\n",
"t01 = qtn.rand_tensor([2, 4, 4], inds=['k0,1', '00-01', '01-11'], tags=['I0,1', 'X0', 'Y0'])\n",
"t10 = qtn.rand_tensor([2, 4, 4], inds=['k1,0', '00-10', '10-11'], tags=['I1,0', 'X1', 'Y1'])\n",
"t11 = qtn.rand_tensor([2, 4, 4], inds=['k1,1', '01-11', '10-11'], tags=['I1,1', 'X1', 'Y1'])"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"These make up a little PEPS-like 2x2 TN:"
]
},
{
"cell_type": "code",
"execution_count": 35,
"metadata": {},
"outputs": [
{
"data": {
"text/html": [
"TensorNetwork (tensors=4, indices=8)Tensor (shape=(2 , 4 , 4 ), inds=[k0,0 , 00-01 , 00-10 ], tags={I0,0 , X0 , Y0 }), backend=numpy , dtype=float64 , data=array([[[-1.10079799, 0.0596712 , 0.11051359, 0.14294902],\n",
" [ 0.48908166, 0.42848178, -0.75141165, -0.44749067],\n",
" [ 0.86966155, 0.31911253, 1.68113226, 2.62923484],\n",
" [ 1.26233162, -0.24575065, 0.50711996, -0.3204049 ]],\n",
"\n",
" [[ 1.67505259, 0.35965073, 1.97026062, -1.42291639],\n",
" [ 0.05697117, 0.63406443, 0.4456802 , 0.74561541],\n",
" [ 1.86274148, -0.98858507, 0.35782756, -0.10051891],\n",
" [-0.11145768, -1.30228361, -0.93935626, -1.40401747]]])Tensor (shape=(2 , 4 , 4 ), inds=[k0,1 , 00-01 , 01-11 ], tags={I0,1 , X0 , Y0 }), backend=numpy , dtype=float64 , data=array([[[ 0.38860772, -0.23785043, 0.53180618, 1.00064276],\n",
" [ 1.81362309, 2.16673389, -1.128748 , 0.24049856],\n",
" [ 0.02505905, 0.86045316, 0.02844933, 0.28908474],\n",
" [-0.35201407, 0.94671534, -1.69440475, 0.15656306]],\n",
"\n",
" [[ 0.11541973, 0.14393308, -1.02624579, -0.18799492],\n",
" [ 1.50641591, 0.97460991, -0.75742392, 0.22991593],\n",
" [ 0.13252203, -0.92264156, 0.87362977, 1.13883566],\n",
" [-0.90282265, -1.60227634, 1.1316849 , -0.41218239]]])Tensor (shape=(2 , 4 , 4 ), inds=[k1,0 , 00-10 , 10-11 ], tags={I1,0 , X1 , Y1 }), backend=numpy , dtype=float64 , data=array([[[ 0.70519837, 0.6616975 , 1.0114127 , -1.3656219 ],\n",
" [-0.02186246, 0.15600819, -3.09894756, -1.00428208],\n",
" [-0.23836214, -0.74322461, -1.58058123, 1.10484127],\n",
" [-1.84649077, -0.03061097, -0.22689998, 0.13139946]],\n",
"\n",
" [[-1.41904566, 0.16417185, -1.35531029, -0.68284241],\n",
" [ 0.16804185, 1.2119743 , -1.04827067, 0.03572145],\n",
" [-0.3913261 , -0.33323323, -0.04251186, 0.69972619],\n",
" [ 0.03534986, -0.81910189, 0.78524742, 0.71433832]]])Tensor (shape=(2 , 4 , 4 ), inds=[k1,1 , 01-11 , 10-11 ], tags={I1,1 , X1 , Y1 }), backend=numpy , dtype=float64 , data=array([[[-0.33062103, 1.48883878, 1.74355326, -0.1238563 ],\n",
" [-1.4736216 , -0.8324853 , -0.86660347, -1.14581125],\n",
" [ 0.69396446, -0.25231889, 1.18294182, 0.37714141],\n",
" [ 1.54063522, 0.40692399, -1.41923073, 0.09172671]],\n",
"\n",
" [[-0.89979094, -1.6911207 , 0.15734748, -2.79835174],\n",
" [-0.07682657, 0.07359371, 0.39672031, -0.66140958],\n",
" [ 0.97479304, -1.40525561, -2.38343859, 0.89768701],\n",
" [-0.93541728, -1.48125354, 0.73259751, -0.84065252]]]) "
],
"text/plain": [
"TensorNetwork(tensors=4, indices=8)"
]
},
"execution_count": 35,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"tn = (t00 | t01 | t10 | t11)\n",
"tn"
]
},
{
"cell_type": "code",
"execution_count": 36,
"metadata": {},
"outputs": [
{
"data": {
"image/svg+xml": [
"\n",
"\n",
"\n",
" \n",
" \n",
" \n",
" \n",
" 2023-11-29T16:12:22.789856 \n",
" image/svg+xml \n",
" \n",
" \n",
" Matplotlib v3.7.2, https://matplotlib.org/ \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n"
],
"text/plain": [
""
]
},
"metadata": {
"needs_background": "light"
},
"output_type": "display_data"
}
],
"source": [
"tn.draw(color=['I0,0', 'I0,1', 'I1,0', 'I1,1'])"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"However it doesn't have any of the special {class}`~quimb.tensor.tensor_2d.PEPS`\n",
"methods yet. Rather than retrieving raw arrays and calling `PEPS(arrays)`,\n",
"we can use the following functions to directly instantiate a PEPS:\n",
"\n",
":::{note}\n",
"Note the following only work because **we have already correctly tagged and\n",
"labelled the tensors**, thus 'promising' that the PEPS structure exists.\n",
":::\n",
"\n",
"### `from_TN`"
]
},
{
"cell_type": "code",
"execution_count": 37,
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"('_site_tag_id', '_x_tag_id', '_y_tag_id', '_Lx', '_Ly', '_site_ind_id')"
]
},
"execution_count": 37,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"# we need to tell the constructor what values we are using for the following:\n",
"qtn.PEPS._EXTRA_PROPS"
]
},
{
"cell_type": "code",
"execution_count": 38,
"metadata": {},
"outputs": [
{
"data": {
"text/html": [
"PEPS (tensors=4, indices=8, Lx=2, Ly=2, max_bond=4)Tensor (shape=(2 , 4 , 4 ), inds=[k0,0 , 00-01 , 00-10 ], tags={I0,0 , X0 , Y0 }), backend=numpy , dtype=float64 , data=array([[[-1.10079799, 0.0596712 , 0.11051359, 0.14294902],\n",
" [ 0.48908166, 0.42848178, -0.75141165, -0.44749067],\n",
" [ 0.86966155, 0.31911253, 1.68113226, 2.62923484],\n",
" [ 1.26233162, -0.24575065, 0.50711996, -0.3204049 ]],\n",
"\n",
" [[ 1.67505259, 0.35965073, 1.97026062, -1.42291639],\n",
" [ 0.05697117, 0.63406443, 0.4456802 , 0.74561541],\n",
" [ 1.86274148, -0.98858507, 0.35782756, -0.10051891],\n",
" [-0.11145768, -1.30228361, -0.93935626, -1.40401747]]])Tensor (shape=(2 , 4 , 4 ), inds=[k0,1 , 00-01 , 01-11 ], tags={I0,1 , X0 , Y0 }), backend=numpy , dtype=float64 , data=array([[[ 0.38860772, -0.23785043, 0.53180618, 1.00064276],\n",
" [ 1.81362309, 2.16673389, -1.128748 , 0.24049856],\n",
" [ 0.02505905, 0.86045316, 0.02844933, 0.28908474],\n",
" [-0.35201407, 0.94671534, -1.69440475, 0.15656306]],\n",
"\n",
" [[ 0.11541973, 0.14393308, -1.02624579, -0.18799492],\n",
" [ 1.50641591, 0.97460991, -0.75742392, 0.22991593],\n",
" [ 0.13252203, -0.92264156, 0.87362977, 1.13883566],\n",
" [-0.90282265, -1.60227634, 1.1316849 , -0.41218239]]])Tensor (shape=(2 , 4 , 4 ), inds=[k1,0 , 00-10 , 10-11 ], tags={I1,0 , X1 , Y1 }), backend=numpy , dtype=float64 , data=array([[[ 0.70519837, 0.6616975 , 1.0114127 , -1.3656219 ],\n",
" [-0.02186246, 0.15600819, -3.09894756, -1.00428208],\n",
" [-0.23836214, -0.74322461, -1.58058123, 1.10484127],\n",
" [-1.84649077, -0.03061097, -0.22689998, 0.13139946]],\n",
"\n",
" [[-1.41904566, 0.16417185, -1.35531029, -0.68284241],\n",
" [ 0.16804185, 1.2119743 , -1.04827067, 0.03572145],\n",
" [-0.3913261 , -0.33323323, -0.04251186, 0.69972619],\n",
" [ 0.03534986, -0.81910189, 0.78524742, 0.71433832]]])Tensor (shape=(2 , 4 , 4 ), inds=[k1,1 , 01-11 , 10-11 ], tags={I1,1 , X1 , Y1 }), backend=numpy , dtype=float64 , data=array([[[-0.33062103, 1.48883878, 1.74355326, -0.1238563 ],\n",
" [-1.4736216 , -0.8324853 , -0.86660347, -1.14581125],\n",
" [ 0.69396446, -0.25231889, 1.18294182, 0.37714141],\n",
" [ 1.54063522, 0.40692399, -1.41923073, 0.09172671]],\n",
"\n",
" [[-0.89979094, -1.6911207 , 0.15734748, -2.79835174],\n",
" [-0.07682657, 0.07359371, 0.39672031, -0.66140958],\n",
" [ 0.97479304, -1.40525561, -2.38343859, 0.89768701],\n",
" [-0.93541728, -1.48125354, 0.73259751, -0.84065252]]]) "
],
"text/plain": [
"PEPS(tensors=4, indices=8, Lx=2, Ly=2, max_bond=4)"
]
},
"execution_count": 38,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"peps = qtn.PEPS.from_TN(\n",
" tn, \n",
" Lx=2, \n",
" Ly=2,\n",
" site_tag_id='I{},{}', \n",
" site_ind_id='k{},{}',\n",
" x_tag_id='X{}',\n",
" y_tag_id='Y{}',\n",
")\n",
"peps"
]
},
{
"cell_type": "code",
"execution_count": 39,
"metadata": {},
"outputs": [
{
"name": "stdout",
"output_type": "stream",
"text": [
" 4 \n",
" ●━━━━●\n",
"╱┃4 ╱┃4 \n",
" ┃ 4 ┃ \n",
" ●━━━━●\n",
"╱ ╱ \n"
]
}
],
"source": [
"peps.show()"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"### `TensorNetwork.view_as`\n",
"\n",
"We can also use the method version,\n",
"[`TensorNetwork.view_as`](quimb.tensor.tensor_core.TensorNetwork.view_as),\n",
"which enables an inplace option:"
]
},
{
"cell_type": "code",
"execution_count": 40,
"metadata": {},
"outputs": [
{
"data": {
"text/html": [
"TensorNetwork2DFlat (tensors=4, indices=8, Lx=2, Ly=2, max_bond=4)Tensor (shape=(2 , 4 , 4 ), inds=[k0,0 , 00-01 , 00-10 ], tags={I0,0 , X0 , Y0 }), backend=numpy , dtype=float64 , data=array([[[-1.10079799, 0.0596712 , 0.11051359, 0.14294902],\n",
" [ 0.48908166, 0.42848178, -0.75141165, -0.44749067],\n",
" [ 0.86966155, 0.31911253, 1.68113226, 2.62923484],\n",
" [ 1.26233162, -0.24575065, 0.50711996, -0.3204049 ]],\n",
"\n",
" [[ 1.67505259, 0.35965073, 1.97026062, -1.42291639],\n",
" [ 0.05697117, 0.63406443, 0.4456802 , 0.74561541],\n",
" [ 1.86274148, -0.98858507, 0.35782756, -0.10051891],\n",
" [-0.11145768, -1.30228361, -0.93935626, -1.40401747]]])Tensor (shape=(2 , 4 , 4 ), inds=[k0,1 , 00-01 , 01-11 ], tags={I0,1 , X0 , Y0 }), backend=numpy , dtype=float64 , data=array([[[ 0.38860772, -0.23785043, 0.53180618, 1.00064276],\n",
" [ 1.81362309, 2.16673389, -1.128748 , 0.24049856],\n",
" [ 0.02505905, 0.86045316, 0.02844933, 0.28908474],\n",
" [-0.35201407, 0.94671534, -1.69440475, 0.15656306]],\n",
"\n",
" [[ 0.11541973, 0.14393308, -1.02624579, -0.18799492],\n",
" [ 1.50641591, 0.97460991, -0.75742392, 0.22991593],\n",
" [ 0.13252203, -0.92264156, 0.87362977, 1.13883566],\n",
" [-0.90282265, -1.60227634, 1.1316849 , -0.41218239]]])Tensor (shape=(2 , 4 , 4 ), inds=[k1,0 , 00-10 , 10-11 ], tags={I1,0 , X1 , Y1 }), backend=numpy , dtype=float64 , data=array([[[ 0.70519837, 0.6616975 , 1.0114127 , -1.3656219 ],\n",
" [-0.02186246, 0.15600819, -3.09894756, -1.00428208],\n",
" [-0.23836214, -0.74322461, -1.58058123, 1.10484127],\n",
" [-1.84649077, -0.03061097, -0.22689998, 0.13139946]],\n",
"\n",
" [[-1.41904566, 0.16417185, -1.35531029, -0.68284241],\n",
" [ 0.16804185, 1.2119743 , -1.04827067, 0.03572145],\n",
" [-0.3913261 , -0.33323323, -0.04251186, 0.69972619],\n",
" [ 0.03534986, -0.81910189, 0.78524742, 0.71433832]]])Tensor (shape=(2 , 4 , 4 ), inds=[k1,1 , 01-11 , 10-11 ], tags={I1,1 , X1 , Y1 }), backend=numpy , dtype=float64 , data=array([[[-0.33062103, 1.48883878, 1.74355326, -0.1238563 ],\n",
" [-1.4736216 , -0.8324853 , -0.86660347, -1.14581125],\n",
" [ 0.69396446, -0.25231889, 1.18294182, 0.37714141],\n",
" [ 1.54063522, 0.40692399, -1.41923073, 0.09172671]],\n",
"\n",
" [[-0.89979094, -1.6911207 , 0.15734748, -2.79835174],\n",
" [-0.07682657, 0.07359371, 0.39672031, -0.66140958],\n",
" [ 0.97479304, -1.40525561, -2.38343859, 0.89768701],\n",
" [-0.93541728, -1.48125354, 0.73259751, -0.84065252]]]) "
],
"text/plain": [
"TensorNetwork2DFlat(tensors=4, indices=8, Lx=2, Ly=2, max_bond=4)"
]
},
"execution_count": 40,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"tn.view_as_(\n",
" qtn.tensor_2d.TensorNetwork2DFlat,\n",
" Lx=2, \n",
" Ly=2,\n",
" site_tag_id='I{},{}',\n",
" x_tag_id='X{}',\n",
" y_tag_id='Y{}',\n",
")"
]
},
{
"cell_type": "code",
"execution_count": 41,
"metadata": {},
"outputs": [
{
"name": "stdout",
"output_type": "stream",
"text": [
" 4 \n",
" ●━━━━●\n",
" ┃4 ┃4 \n",
" ┃ 4 ┃ \n",
" ●━━━━●\n",
" \n"
]
}
],
"source": [
"# `tn` is now a 'FLat 2D TN'\n",
"tn.show()"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"### `TensorNetwork.view_like`\n",
"\n",
"Finally, and often most conveniently, if you have an\n",
"existing structured TN with extra attributes you want to match,\n",
"you can call\n",
"[`TensorNetwork.view_like`]quimb.tensor.tensor_core.TensorNetwork.view_like), which\n",
"defaults to picking up all the necessary attributes from whatever\n",
"you supply to the `like=` kwarg:"
]
},
{
"cell_type": "code",
"execution_count": 42,
"metadata": {},
"outputs": [
{
"data": {
"text/html": [
"PEPS (tensors=4, indices=8, Lx=2, Ly=2, max_bond=4)Tensor (shape=(2 , 4 , 4 ), inds=[k0,0 , 00-01 , 00-10 ], tags={I0,0 , X0 , Y0 }), backend=numpy , dtype=float64 , data=array([[[-1.10079799, 0.0596712 , 0.11051359, 0.14294902],\n",
" [ 0.48908166, 0.42848178, -0.75141165, -0.44749067],\n",
" [ 0.86966155, 0.31911253, 1.68113226, 2.62923484],\n",
" [ 1.26233162, -0.24575065, 0.50711996, -0.3204049 ]],\n",
"\n",
" [[ 1.67505259, 0.35965073, 1.97026062, -1.42291639],\n",
" [ 0.05697117, 0.63406443, 0.4456802 , 0.74561541],\n",
" [ 1.86274148, -0.98858507, 0.35782756, -0.10051891],\n",
" [-0.11145768, -1.30228361, -0.93935626, -1.40401747]]])Tensor (shape=(2 , 4 , 4 ), inds=[k0,1 , 00-01 , 01-11 ], tags={I0,1 , X0 , Y0 }), backend=numpy , dtype=float64 , data=array([[[ 0.38860772, -0.23785043, 0.53180618, 1.00064276],\n",
" [ 1.81362309, 2.16673389, -1.128748 , 0.24049856],\n",
" [ 0.02505905, 0.86045316, 0.02844933, 0.28908474],\n",
" [-0.35201407, 0.94671534, -1.69440475, 0.15656306]],\n",
"\n",
" [[ 0.11541973, 0.14393308, -1.02624579, -0.18799492],\n",
" [ 1.50641591, 0.97460991, -0.75742392, 0.22991593],\n",
" [ 0.13252203, -0.92264156, 0.87362977, 1.13883566],\n",
" [-0.90282265, -1.60227634, 1.1316849 , -0.41218239]]])Tensor (shape=(2 , 4 , 4 ), inds=[k1,0 , 00-10 , 10-11 ], tags={I1,0 , X1 , Y1 }), backend=numpy , dtype=float64 , data=array([[[ 0.70519837, 0.6616975 , 1.0114127 , -1.3656219 ],\n",
" [-0.02186246, 0.15600819, -3.09894756, -1.00428208],\n",
" [-0.23836214, -0.74322461, -1.58058123, 1.10484127],\n",
" [-1.84649077, -0.03061097, -0.22689998, 0.13139946]],\n",
"\n",
" [[-1.41904566, 0.16417185, -1.35531029, -0.68284241],\n",
" [ 0.16804185, 1.2119743 , -1.04827067, 0.03572145],\n",
" [-0.3913261 , -0.33323323, -0.04251186, 0.69972619],\n",
" [ 0.03534986, -0.81910189, 0.78524742, 0.71433832]]])Tensor (shape=(2 , 4 , 4 ), inds=[k1,1 , 01-11 , 10-11 ], tags={I1,1 , X1 , Y1 }), backend=numpy , dtype=float64 , data=array([[[-0.33062103, 1.48883878, 1.74355326, -0.1238563 ],\n",
" [-1.4736216 , -0.8324853 , -0.86660347, -1.14581125],\n",
" [ 0.69396446, -0.25231889, 1.18294182, 0.37714141],\n",
" [ 1.54063522, 0.40692399, -1.41923073, 0.09172671]],\n",
"\n",
" [[-0.89979094, -1.6911207 , 0.15734748, -2.79835174],\n",
" [-0.07682657, 0.07359371, 0.39672031, -0.66140958],\n",
" [ 0.97479304, -1.40525561, -2.38343859, 0.89768701],\n",
" [-0.93541728, -1.48125354, 0.73259751, -0.84065252]]]) "
],
"text/plain": [
"PEPS(tensors=4, indices=8, Lx=2, Ly=2, max_bond=4)"
]
},
"execution_count": 42,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"other_peps = qtn.PEPS.rand(Lx=2, Ly=2, bond_dim=5)\n",
"\n",
"tn.view_like_(other_peps)"
]
},
{
"cell_type": "code",
"execution_count": 43,
"metadata": {},
"outputs": [
{
"name": "stdout",
"output_type": "stream",
"text": [
" 9 \n",
" ●━━━━●\n",
"╱┃9 ╱┃9 \n",
" ┃ 9 ┃ \n",
" ●━━━━●\n",
"╱ ╱ \n"
]
}
],
"source": [
"# `tn` is now a PEPS with its special methods\n",
"tn.expand_bond_dimension(9).show()"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"This is useful if you perform modifications to a specific type of TN,\n",
"such that it loses its structured identity, but then you want to recast\n",
"it as original type, knowing that all the tags and indices are still correct.\n",
"\n",
"### Compatible Subclasses\n",
"\n",
"One final feature to note regarding specific tensor network\n",
"subclasses is that when you combine certain TNs with the\n",
"`&` or `|` operators, if they are both 'compatible' versions\n",
"of an inherited structured TN, they keep that structure:"
]
},
{
"cell_type": "code",
"execution_count": 44,
"metadata": {},
"outputs": [
{
"data": {
"text/html": [
"TensorNetwork1D (tensors=20, indices=28, L=10, max_bond=7)Tensor (shape=(7 , 2 ), inds=[_5ea058AAADz , k0 ], tags={I0 }), backend=numpy , dtype=float64 , data=array([[ 0.10667958, 0.43460376],\n",
" [ 0.18356072, 0.43281806],\n",
" [-0.47654754, -0.56836425],\n",
" [ 0.15192568, -0.23423402],\n",
" [-0.01854748, 0.31846122],\n",
" [-0.72145024, 0.27601746],\n",
" [-0.04482127, -0.32750076]])Tensor (shape=(7 , 7 , 2 ), inds=[_5ea058AAADz , _5ea058AAAEA , k1 ], tags={I1 }), backend=numpy , dtype=float64 , data=array([[[ 0.36463648, 0.30896422],\n",
" [ 0.11027168, 0.21970176],\n",
" [ 0.30721991, 0.15809952],\n",
" [-0.42111672, 0.05444968],\n",
" [ 0.65572221, 0.31595011],\n",
" [ 0.14696424, -0.22293619],\n",
" [-0.3376481 , -0.10456457]],\n",
"\n",
" [[ 0.31075154, 0.28757696],\n",
" [ 0.15469386, 0.38839358],\n",
" [-0.08841303, 0.1111061 ],\n",
" [-0.26989522, -0.06952075],\n",
" [ 0.02361441, 0.10043134],\n",
" [ 0.52958256, 0.33113841],\n",
" [-0.3940553 , 0.1510596 ]],\n",
"\n",
" [[ 0.44370206, 0.41030918],\n",
" [ 0.1495414 , -0.44699449],\n",
" [-0.29230541, -0.09843075],\n",
" [ 0.68732758, -0.01888702],\n",
" [-0.22780574, 0.10113459],\n",
" [-0.09463108, 0.19805312],\n",
" [-0.18875575, 0.16284019]],\n",
"\n",
" [[-0.25523676, -0.30228712],\n",
" [ 0.50935002, -0.20146142],\n",
" [-0.13966482, -0.45600358],\n",
" [ 0.08364221, -0.01145466],\n",
" [ 0.13890749, -0.20123609],\n",
" [ 0.21198093, 0.44293016],\n",
" [ 0.09268257, 0.43444058]],\n",
"\n",
" [[ 0.28231305, -0.08302645],\n",
" [-0.0208292 , 0.48394883],\n",
" [-0.09448289, -0.56449763],\n",
" [-0.23068247, 0.39692065],\n",
" [-0.12987429, 0.41848041],\n",
" [-0.27097851, 0.05013033],\n",
" [ 0.0061537 , -0.21103555]],\n",
"\n",
" [[-0.08354841, 0.02570106],\n",
" [ 0.04664423, 0.05903268],\n",
" [-0.25907206, 0.54869689],\n",
" [ 0.36122059, -0.74456619],\n",
" [-0.49386319, 0.03206811],\n",
" [-0.24166144, -0.42894566],\n",
" [ 0.06146449, 0.28889836]],\n",
"\n",
" [[ 0.2979051 , 0.21847142],\n",
" [-0.32912169, 0.34359257],\n",
" [ 0.71212591, 0.25432835],\n",
" [ 0.27494262, -0.08021264],\n",
" [ 0.03841404, 0.37004461],\n",
" [-0.06217767, 0.08544313],\n",
" [-0.28950812, 0.69772085]]])Tensor (shape=(7 , 7 , 2 ), inds=[_5ea058AAAEA , _5ea058AAAEB , k2 ], tags={I2 }), backend=numpy , dtype=float64 , data=array([[[-0.0674622 , 0.09531434],\n",
" [ 0.15364688, 0.15622314],\n",
" [-0.21357801, -0.41162665],\n",
" [ 0.6631182 , -0.28943395],\n",
" [ 0.26401233, 0.06741417],\n",
" [-0.0411162 , -0.0489503 ],\n",
" [ 0.02282935, -0.18976925]],\n",
"\n",
" [[ 0.20116235, 0.65447471],\n",
" [ 0.32566958, 0.20191934],\n",
" [ 0.38621883, 0.09371959],\n",
" [ 0.19572957, -0.43756167],\n",
" [ 0.41129239, -0.22185977],\n",
" [-0.03347971, 0.30013001],\n",
" [-0.24000569, 0.52581296]],\n",
"\n",
" [[ 0.42691706, 0.23970259],\n",
" [ 0.00223142, -0.97076676],\n",
" [ 0.45769233, -0.13484798],\n",
" [ 0.26981586, -0.13052823],\n",
" [ 0.59690455, 0.29096862],\n",
" [-0.14098378, -0.13961376],\n",
" [ 0.29224598, 0.26579672]],\n",
"\n",
" [[ 0.14108392, -0.33582179],\n",
" [-0.39195789, -0.06199785],\n",
" [ 0.13640521, -0.06081178],\n",
" [ 0.03779749, -0.51120059],\n",
" [ 0.26986231, -0.08686723],\n",
" [-0.16874409, 0.27254808],\n",
" [ 0.08551108, 0.10663054]],\n",
"\n",
" [[ 0.12609989, -0.29698382],\n",
" [ 0.27114845, -0.13791775],\n",
" [-0.42853545, -0.22433165],\n",
" [ 0.56818702, -0.19075512],\n",
" [ 0.74381924, -0.73046484],\n",
" [ 0.20199693, 0.26213131],\n",
" [-0.03560184, 0.29315886]],\n",
"\n",
" [[ 0.00679922, -0.29847277],\n",
" [ 0.01639847, 0.04672636],\n",
" [ 0.3687735 , -0.10892058],\n",
" [-0.055516 , -0.02295742],\n",
" [-0.11262614, 0.2034438 ],\n",
" [ 0.20534674, -0.01656091],\n",
" [ 0.80461272, 0.37958948]],\n",
"\n",
" [[-0.04886204, 0.02001407],\n",
" [ 0.19129469, 0.22252211],\n",
" [ 0.2164089 , 0.10759276],\n",
" [ 0.5321792 , 0.5016342 ],\n",
" [ 0.08243591, -0.23207251],\n",
" [ 0.21088515, 0.55510654],\n",
" [-0.25027732, 0.34689043]]])Tensor (shape=(7 , 7 , 2 ), inds=[_5ea058AAAEB , _5ea058AAAEC , k3 ], tags={I3 }), backend=numpy , dtype=float64 , data=array([[[ 0.01877289, -0.46872294],\n",
" [-0.32273201, -0.52393566],\n",
" [ 0.14766874, 0.22990295],\n",
" [ 0.46738969, -0.2962717 ],\n",
" [-0.32613964, -0.09734032],\n",
" [ 0.09629676, -0.72959512],\n",
" [-0.50785321, 0.09920142]],\n",
"\n",
" [[ 0.00571944, -0.02300505],\n",
" [ 0.11230096, -0.00875314],\n",
" [ 0.21619887, 0.23794723],\n",
" [-0.02395151, 0.06949099],\n",
" [-0.38053327, -0.35882451],\n",
" [-0.10114924, 0.19781418],\n",
" [-0.32650307, 0.50972773]],\n",
"\n",
" [[-0.17568399, -0.24117432],\n",
" [-0.02346306, 0.30784853],\n",
" [ 0.28284282, -0.16375829],\n",
" [ 0.44676813, -0.11073052],\n",
" [-0.47283812, -0.43061851],\n",
" [-0.03967218, -0.02517368],\n",
" [ 0.35029486, 0.38697832]],\n",
"\n",
" [[-0.31207125, 0.34888349],\n",
" [-0.06298621, 0.14020835],\n",
" [-0.53112243, 0.12062568],\n",
" [-0.24809299, -0.35613623],\n",
" [-0.67188807, -0.01809549],\n",
" [-0.47233916, 0.01759363],\n",
" [ 0.01025288, -0.00743861]],\n",
"\n",
" [[-0.96222181, 0.44479526],\n",
" [-0.40429563, 0.41391342],\n",
" [-0.23238099, 0.84558789],\n",
" [ 0.2833067 , 0.20216819],\n",
" [-0.57650274, -0.13932612],\n",
" [-0.09260767, -0.1607125 ],\n",
" [ 0.03761778, 0.30235942]],\n",
"\n",
" [[-0.40370173, 0.10895752],\n",
" [ 0.11319523, 0.3211773 ],\n",
" [ 0.22661239, -0.23341429],\n",
" [ 0.60734721, -0.08282081],\n",
" [-0.07995627, -0.40095008],\n",
" [ 0.20150293, -0.04646896],\n",
" [-0.26213453, -0.24730749]],\n",
"\n",
" [[ 0.01166313, 0.29810887],\n",
" [ 0.24602542, 0.23117637],\n",
" [ 0.39975465, -0.36875793],\n",
" [-0.03667836, 0.32735391],\n",
" [-0.06952987, -0.01338072],\n",
" [ 0.23284358, -0.03186916],\n",
" [ 0.09296511, 0.35095442]]])Tensor (shape=(7 , 7 , 2 ), inds=[_5ea058AAAEC , _5ea058AAAED , k4 ], tags={I4 }), backend=numpy , dtype=float64 , data=array([[[-4.37366253e-01, -9.58257257e-02],\n",
" [ 5.94729822e-03, 5.86142358e-01],\n",
" [-1.32582101e-01, -1.05507291e-01],\n",
" [ 4.06790712e-01, 5.06505157e-02],\n",
" [ 4.14753982e-01, -2.90224835e-01],\n",
" [-5.05053403e-01, -9.29292347e-02],\n",
" [ 3.10814503e-01, -5.22636538e-02]],\n",
"\n",
" [[ 2.86446393e-01, -3.14892027e-01],\n",
" [-3.28840177e-02, -1.63118851e-01],\n",
" [ 3.02865561e-01, 2.62736053e-01],\n",
" [ 3.71641055e-01, -1.33745452e-01],\n",
" [ 1.51440502e-01, -4.28394966e-01],\n",
" [-7.49159663e-01, 1.27625388e-02],\n",
" [-2.55485760e-01, -1.20265476e-01]],\n",
"\n",
" [[ 9.48451482e-02, -1.71938086e-01],\n",
" [ 4.53733950e-01, -4.03677015e-01],\n",
" [ 1.56149253e-01, 2.24257627e-01],\n",
" [ 6.04253433e-01, 1.49966817e-01],\n",
" [-1.73911953e-01, 3.34109428e-01],\n",
" [ 2.04143045e-01, 6.44356234e-01],\n",
" [ 2.51640379e-01, 6.95786954e-01]],\n",
"\n",
" [[-2.59426186e-02, 1.81668832e-01],\n",
" [-6.90665156e-02, 8.12671695e-02],\n",
" [ 6.18857017e-04, -1.96899058e-01],\n",
" [-4.47952354e-02, -1.46856986e-01],\n",
" [-7.24173341e-01, 2.14029251e-01],\n",
" [-3.37541552e-02, -2.37108557e-02],\n",
" [ 2.80703333e-02, -5.51674220e-02]],\n",
"\n",
" [[ 3.88436394e-01, -2.98399687e-01],\n",
" [-8.02585897e-01, 7.21096401e-02],\n",
" [ 1.64654421e-01, -8.05473521e-01],\n",
" [ 3.90798150e-01, 4.81143158e-01],\n",
" [-1.04691886e-01, 2.24809706e-02],\n",
" [ 7.68435437e-02, 5.50124520e-01],\n",
" [-2.31065940e-01, -9.36867669e-04]],\n",
"\n",
" [[-4.72971863e-02, 2.63528032e-01],\n",
" [ 3.91571653e-01, 4.88336618e-01],\n",
" [ 5.75985494e-02, -1.61371278e-01],\n",
" [ 6.78926864e-01, -7.88697575e-02],\n",
" [ 7.19814058e-02, -2.47750798e-01],\n",
" [ 1.83073701e-01, 1.47386788e-01],\n",
" [-5.26582144e-01, 3.71396017e-01]],\n",
"\n",
" [[ 2.50881906e-01, 2.18197787e-01],\n",
" [-6.01454192e-01, 7.24352750e-02],\n",
" [ 1.69706572e-01, 9.29186700e-03],\n",
" [-2.74899397e-01, -2.17504562e-01],\n",
" [ 2.36616703e-01, 4.00599110e-01],\n",
" [ 1.15157890e-01, 3.28059035e-01],\n",
" [ 4.64251333e-01, -4.86459666e-01]]])Tensor (shape=(7 , 7 , 2 ), inds=[_5ea058AAAED , _5ea058AAAEE , k5 ], tags={I5 }), backend=numpy , dtype=float64 , data=array([[[-0.37319043, -0.24048094],\n",
" [ 0.26421341, -0.2955586 ],\n",
" [ 0.14359197, -0.24042257],\n",
" [-0.48558971, -0.15783224],\n",
" [ 0.48457079, -0.45503191],\n",
" [ 0.27530444, 0.37429223],\n",
" [-0.20017181, -0.47650099]],\n",
"\n",
" [[ 0.26580465, 0.50552621],\n",
" [ 0.05109448, 0.28938821],\n",
" [ 0.05460396, -0.06116214],\n",
" [-0.04545442, -0.38893047],\n",
" [ 0.33268127, 0.63356642],\n",
" [-0.00518775, -0.1182908 ],\n",
" [-0.15145152, 0.38530031]],\n",
"\n",
" [[-0.63662418, 0.22903477],\n",
" [-0.12927536, 0.18247706],\n",
" [ 0.1124482 , 0.17705573],\n",
" [ 0.79587811, 0.15708001],\n",
" [ 0.30540201, 0.58779977],\n",
" [ 0.20052386, 0.34651597],\n",
" [-0.41739019, 0.40943056]],\n",
"\n",
" [[-0.34403997, -0.12042728],\n",
" [ 0.30757856, 0.13171954],\n",
" [-0.01963703, 0.09991927],\n",
" [-0.44709559, 0.52676178],\n",
" [-0.13306523, 0.3764855 ],\n",
" [-0.21889239, -0.19986471],\n",
" [-0.13660681, 0.24563983]],\n",
"\n",
" [[ 0.14774158, 0.04077004],\n",
" [ 0.42138934, -0.30846518],\n",
" [ 0.32909284, 0.09477495],\n",
" [ 0.09579384, -0.01579478],\n",
" [-0.16357703, 0.58785853],\n",
" [-0.2296893 , -0.27447461],\n",
" [-0.44506946, -0.04805672]],\n",
"\n",
" [[-0.22614202, -0.06629081],\n",
" [ 0.1520737 , -0.30768033],\n",
" [-0.40959512, -0.21091417],\n",
" [ 0.03813645, -0.46538803],\n",
" [-0.06145405, -0.05055765],\n",
" [-0.10123858, 0.27761483],\n",
" [-0.17135357, 0.43209868]],\n",
"\n",
" [[ 0.44679691, -0.50490497],\n",
" [ 0.11513539, 0.15390761],\n",
" [ 0.0393863 , -0.1639887 ],\n",
" [-0.14713693, -0.32184128],\n",
" [-0.25192071, 0.12723526],\n",
" [ 0.29058953, -0.16738317],\n",
" [-0.98102126, -0.04164749]]])Tensor (shape=(7 , 7 , 2 ), inds=[_5ea058AAAEE , _5ea058AAAEF , k6 ], tags={I6 }), backend=numpy , dtype=float64 , data=array([[[ 0.33646236, -0.23158709],\n",
" [-0.0990975 , -0.01770377],\n",
" [ 0.08546441, -0.02207693],\n",
" [ 0.30466761, 0.1849371 ],\n",
" [-0.01761516, -0.10857565],\n",
" [ 0.31819795, 0.43944061],\n",
" [ 0.20171031, -0.23998167]],\n",
"\n",
" [[ 0.00993387, 0.75810682],\n",
" [-0.54994766, -0.01721529],\n",
" [ 0.25201891, 0.12435247],\n",
" [ 0.21839765, 0.23557524],\n",
" [-0.22855849, -0.08689068],\n",
" [-0.00521593, -0.12853708],\n",
" [-0.25892937, 0.33326892]],\n",
"\n",
" [[-0.30996272, -0.01771899],\n",
" [-0.39106169, -0.09125562],\n",
" [ 0.09152637, 0.01754536],\n",
" [-0.16667507, 0.30241393],\n",
" [ 0.3920397 , 0.2828755 ],\n",
" [ 0.16605798, -0.24619574],\n",
" [-0.07172969, 0.31317368]],\n",
"\n",
" [[ 0.39665261, -0.59550305],\n",
" [-0.25453357, -0.17268722],\n",
" [-0.05346267, -0.08504999],\n",
" [-0.02641703, 0.40507495],\n",
" [ 0.41651265, 0.02810667],\n",
" [-0.6307981 , -0.09185413],\n",
" [-0.19168373, -0.15503258]],\n",
"\n",
" [[-0.56857567, -0.00905714],\n",
" [-0.0718782 , -0.21073597],\n",
" [-0.16097361, 0.12340444],\n",
" [ 0.7915139 , 0.36682742],\n",
" [ 0.29904958, -0.05178732],\n",
" [-0.15951871, 0.24854104],\n",
" [-0.05548262, -0.4102568 ]],\n",
"\n",
" [[ 0.14383098, 0.3309253 ],\n",
" [ 0.1629967 , -0.09343837],\n",
" [-0.58783769, 0.29034635],\n",
" [-0.10327785, -0.19325075],\n",
" [ 0.0157258 , -0.06917776],\n",
" [ 0.14086102, 0.46545015],\n",
" [ 0.59108785, 0.24374455]],\n",
"\n",
" [[-0.40337958, 0.39045244],\n",
" [-0.12640999, 0.28484834],\n",
" [ 0.10607389, -0.05495159],\n",
" [-0.11195234, -0.24112136],\n",
" [-0.26000908, 0.58962589],\n",
" [ 0.44485924, 0.30720018],\n",
" [ 0.29958977, -0.28277345]]])Tensor (shape=(7 , 7 , 2 ), inds=[_5ea058AAAEF , _5ea058AAAEG , k7 ], tags={I7 }), backend=numpy , dtype=float64 , data=array([[[-0.65800368, -0.12871699],\n",
" [-0.2267694 , -0.40363445],\n",
" [-0.08374168, 0.13119852],\n",
" [-0.03592643, -0.55884454],\n",
" [ 0.44337707, 0.18501541],\n",
" [-0.23693387, 0.71610163],\n",
" [-0.31863494, -0.49842103]],\n",
"\n",
" [[-0.4650215 , -0.0588458 ],\n",
" [-0.22587947, -0.53966384],\n",
" [ 0.28018098, -0.10022684],\n",
" [ 0.15699422, -0.13864303],\n",
" [ 0.05603301, 0.24891115],\n",
" [-0.39931207, 0.14955892],\n",
" [-0.72154545, -0.41683622]],\n",
"\n",
" [[ 0.09052714, -0.01655503],\n",
" [-0.21155056, 0.11277397],\n",
" [ 0.12414092, -0.43916799],\n",
" [ 0.13666432, 0.23105145],\n",
" [-0.13002448, -0.10921213],\n",
" [-0.12279022, 0.39396002],\n",
" [ 0.27803708, -0.00229879]],\n",
"\n",
" [[-0.34842423, 0.33007101],\n",
" [ 0.29752476, 0.15363679],\n",
" [ 0.57046213, 0.34905194],\n",
" [ 0.31231561, 0.0089798 ],\n",
" [-0.10953436, -0.84333714],\n",
" [ 0.04308136, 0.03969103],\n",
" [-0.1342694 , 0.05800314]],\n",
"\n",
" [[-0.1680112 , 0.34590386],\n",
" [ 0.39473769, 0.27517079],\n",
" [-0.36191761, -0.03798922],\n",
" [-0.40307874, -0.06952252],\n",
" [ 0.2025674 , -0.20939284],\n",
" [ 0.2422281 , -0.01955989],\n",
" [ 0.03105964, -0.22304362]],\n",
"\n",
" [[ 0.28067431, 0.33849262],\n",
" [ 0.35385787, 0.0497614 ],\n",
" [-0.20989213, -0.20203254],\n",
" [ 0.16103388, 0.00209478],\n",
" [ 0.08297031, 0.26294267],\n",
" [ 0.49397951, -0.24814357],\n",
" [-0.29752701, -0.40784807]],\n",
"\n",
" [[ 0.15512403, 1.11251792],\n",
" [ 0.31413473, 0.08692281],\n",
" [-0.27255984, -0.14534735],\n",
" [-0.17142208, -0.28043824],\n",
" [ 0.09315027, 0.08817523],\n",
" [-0.58168943, -0.70566156],\n",
" [-0.19742525, -0.17252979]]])Tensor (shape=(7 , 7 , 2 ), inds=[_5ea058AAAEG , _5ea058AAAEH , k8 ], tags={I8 }), backend=numpy , dtype=float64 , data=array([[[ 0.40793018, 0.73757475],\n",
" [-0.01552022, 0.7196414 ],\n",
" [ 0.44014854, 0.07947017],\n",
" [ 0.2728351 , -0.51344033],\n",
" [ 0.36840917, 0.04201604],\n",
" [ 0.35868278, 0.27465024],\n",
" [-0.13424817, -0.06131688]],\n",
"\n",
" [[-0.47464014, -0.03673834],\n",
" [ 0.10086346, 0.20311108],\n",
" [-0.37488543, -0.02355421],\n",
" [ 0.10076884, 0.37445736],\n",
" [-0.28314453, 0.13011809],\n",
" [ 0.44963574, 0.68127882],\n",
" [-0.33825994, -0.05233583]],\n",
"\n",
" [[ 0.45425561, -0.12516105],\n",
" [ 0.25616401, 0.09762562],\n",
" [-0.48146508, -0.10824366],\n",
" [ 0.24154111, 0.22869483],\n",
" [-0.01044511, -0.35609876],\n",
" [ 0.57737104, 0.01781104],\n",
" [-0.25961963, 0.2295816 ]],\n",
"\n",
" [[-0.50145923, -0.06361463],\n",
" [-0.2705376 , -0.24379656],\n",
" [ 0.07847862, 0.16573037],\n",
" [ 0.2647418 , 0.25473303],\n",
" [-0.35044035, -0.52064109],\n",
" [ 0.40498533, 0.31632263],\n",
" [ 0.21653174, 0.35346459]],\n",
"\n",
" [[-0.28976109, 0.03497829],\n",
" [-0.10360018, -0.00098986],\n",
" [-0.02633051, -0.14786901],\n",
" [-0.18282857, 0.0332661 ],\n",
" [-0.117132 , -0.63453923],\n",
" [ 0.0634816 , -0.18597299],\n",
" [-0.26441395, -0.26023838]],\n",
"\n",
" [[ 0.50188394, 0.43775824],\n",
" [-0.1708556 , 0.04291476],\n",
" [-0.02036392, -0.07761347],\n",
" [-0.2146938 , -0.21185677],\n",
" [-0.46667896, 0.35474908],\n",
" [ 0.36145081, -0.11671075],\n",
" [ 0.12656567, -0.13309296]],\n",
"\n",
" [[ 0.06034238, 0.15806676],\n",
" [-0.06910878, -0.31254089],\n",
" [ 0.06831289, 0.05236893],\n",
" [ 0.50363766, -0.00090138],\n",
" [ 0.08649658, -0.26779918],\n",
" [-0.168796 , 0.0617152 ],\n",
" [ 0.44480359, -0.02352649]]])Tensor (shape=(7 , 2 ), inds=[_5ea058AAAEH , k9 ], tags={I9 }), backend=numpy , dtype=float64 , data=array([[ 0.00077166, 0.17227245],\n",
" [-0.13695738, 0.05373427],\n",
" [-0.05799634, 0.03574974],\n",
" [-0.14151641, 0.20650976],\n",
" [-0.05642301, 0.07264482],\n",
" [-0.08117216, -0.17874181],\n",
" [-0.20722386, -0.00054867]])Tensor (shape=(7 , 2 ), inds=[_5ea058AAAER , k0 ], tags={I0 }), backend=numpy , dtype=float64 , data=array([[ 0.47395716, -0.5498284 ],\n",
" [-0.04183531, -0.08086085],\n",
" [-0.04009412, -0.24181029],\n",
" [-0.85635818, 0.19401176],\n",
" [ 0.37738836, -0.00809324],\n",
" [-0.17171806, -0.58638134],\n",
" [ 0.21554743, -0.088293 ]])Tensor (shape=(7 , 7 , 2 ), inds=[_5ea058AAAER , _5ea058AAAES , k1 ], tags={I1 }), backend=numpy , dtype=float64 , data=array([[[ 0.37936801, -0.33524963],\n",
" [-0.00340863, -0.0716447 ],\n",
" [-0.40136615, -0.23344966],\n",
" [-0.28717451, -0.43441655],\n",
" [-0.14551537, -0.47234646],\n",
" [ 0.06015349, -0.00376255],\n",
" [-0.17055977, 0.33431599]],\n",
"\n",
" [[ 0.18749911, 0.44551025],\n",
" [-0.27165046, -0.26065426],\n",
" [-0.09622248, -0.02663148],\n",
" [-0.26903189, 0.2954287 ],\n",
" [-0.32594954, -0.03835274],\n",
" [-0.47602812, 0.31669006],\n",
" [ 0.1240176 , -0.24603839]],\n",
"\n",
" [[-0.32539047, -0.20313063],\n",
" [-0.18039085, -0.23169417],\n",
" [ 0.64880353, 0.12534983],\n",
" [ 0.4086111 , -0.26397866],\n",
" [-0.14896794, -0.63225808],\n",
" [-0.22194728, 1.0425243 ],\n",
" [ 0.3497811 , 0.05070243]],\n",
"\n",
" [[-0.11564027, 0.6999554 ],\n",
" [ 0.11100782, 0.09015069],\n",
" [ 0.3424783 , 0.11267941],\n",
" [-0.32007561, -0.67625251],\n",
" [-0.14753574, -0.68019184],\n",
" [ 0.5979037 , 0.26525484],\n",
" [-0.18916273, -0.07389329]],\n",
"\n",
" [[ 0.02814397, 0.15646152],\n",
" [-0.83556811, 0.19369412],\n",
" [ 0.0467395 , 0.05109171],\n",
" [ 0.14227155, 0.26713256],\n",
" [-0.21474369, -0.12171365],\n",
" [ 0.11172592, -0.20680093],\n",
" [ 0.30195661, -0.31662739]],\n",
"\n",
" [[ 0.08180827, -0.24563849],\n",
" [-0.29335326, -0.3103717 ],\n",
" [-0.1112403 , -0.11045994],\n",
" [ 0.10624333, -0.37794902],\n",
" [-0.43856719, 0.01241654],\n",
" [ 0.51704076, 0.2858201 ],\n",
" [ 0.77424708, 0.31145638]],\n",
"\n",
" [[-0.24639427, -0.01445314],\n",
" [-0.04828148, -0.13911124],\n",
" [ 0.50270079, 0.16621341],\n",
" [-0.02136756, -0.12510867],\n",
" [ 0.17637523, -0.05745199],\n",
" [ 0.05684481, -0.25533401],\n",
" [-0.57380813, -0.1062739 ]]])Tensor (shape=(7 , 7 , 2 ), inds=[_5ea058AAAES , _5ea058AAAET , k2 ], tags={I2 }), backend=numpy , dtype=float64 , data=array([[[ 0.18617003, -0.20597041],\n",
" [-0.0406378 , 0.41378189],\n",
" [ 0.12622731, -0.32095298],\n",
" [ 0.13739043, 0.11821688],\n",
" [-0.57861509, 0.10228486],\n",
" [ 0.21299202, -0.17573596],\n",
" [ 0.20379633, 0.14619856]],\n",
"\n",
" [[ 0.25206414, -0.14748162],\n",
" [ 0.31547876, -0.31596766],\n",
" [-0.33494531, -0.61102241],\n",
" [ 0.33621982, 0.65216703],\n",
" [ 0.76699593, 0.19788671],\n",
" [ 0.39924671, -0.27471272],\n",
" [-0.19833743, 0.1583553 ]],\n",
"\n",
" [[-0.19148593, 0.31003284],\n",
" [-0.31336887, 0.3961306 ],\n",
" [-0.14011545, 0.2264034 ],\n",
" [-0.44780403, 0.44858553],\n",
" [ 0.13063928, -0.12151794],\n",
" [-0.17714548, 0.48442646],\n",
" [ 0.28540035, 0.2926135 ]],\n",
"\n",
" [[ 0.2075224 , -0.41541152],\n",
" [-0.28158466, -0.41812664],\n",
" [-0.15470568, -0.55597806],\n",
" [-0.21738657, -0.00121673],\n",
" [ 0.08745507, -0.65638892],\n",
" [ 0.33671299, -0.04528144],\n",
" [-0.47511456, 0.23474704]],\n",
"\n",
" [[ 0.61249883, 0.25040958],\n",
" [-0.02364388, -0.13047417],\n",
" [-0.30103411, -0.68340891],\n",
" [ 0.25342073, 0.62377141],\n",
" [ 0.1443986 , -0.29887648],\n",
" [ 0.43886196, 0.33722572],\n",
" [ 0.34918237, -0.13933612]],\n",
"\n",
" [[ 0.3287949 , -0.23236067],\n",
" [-0.42462116, -0.05373624],\n",
" [-0.36562411, 0.36594392],\n",
" [ 0.25887158, 0.04388555],\n",
" [-0.20741747, -0.25791231],\n",
" [-0.15105595, -0.21845654],\n",
" [-0.08856615, -0.16183803]],\n",
"\n",
" [[ 0.03535534, -0.10461482],\n",
" [ 0.01657496, 0.27218973],\n",
" [-0.06205513, -0.51618432],\n",
" [-0.26280241, -0.25926144],\n",
" [ 0.1929324 , -0.22571108],\n",
" [-0.22882741, -0.26216967],\n",
" [ 0.34696286, 0.16661178]]])Tensor (shape=(7 , 7 , 2 ), inds=[_5ea058AAAET , _5ea058AAAEU , k3 ], tags={I3 }), backend=numpy , dtype=float64 , data=array([[[-0.61429138, 0.15948804],\n",
" [ 0.11135708, 0.1793775 ],\n",
" [ 0.02364871, 0.20209696],\n",
" [ 0.33468274, -0.17069658],\n",
" [ 0.22944445, 0.6038976 ],\n",
" [-0.17303623, 0.13690422],\n",
" [ 0.14040089, -0.15102526]],\n",
"\n",
" [[-0.58961828, -0.15462562],\n",
" [ 0.1529977 , 0.01058179],\n",
" [-0.21566089, -0.20922437],\n",
" [-0.30491512, 0.42069114],\n",
" [-0.30145244, 0.08673978],\n",
" [-0.02473495, 0.07622847],\n",
" [ 0.13491501, 0.21787682]],\n",
"\n",
" [[-0.48947543, -0.18206527],\n",
" [-0.91481232, 0.42034758],\n",
" [-0.19826054, -0.13786642],\n",
" [ 0.11704214, 0.27013931],\n",
" [ 0.34059784, -0.20993431],\n",
" [-0.01852917, -0.20168444],\n",
" [-0.16445711, -0.13836657]],\n",
"\n",
" [[-0.00167299, 0.20938732],\n",
" [ 0.45506169, 0.49781732],\n",
" [ 0.31355304, -0.45655236],\n",
" [ 0.31657332, -0.41331866],\n",
" [ 0.38170817, 0.58938958],\n",
" [-0.065245 , -0.22745516],\n",
" [ 0.12223286, -0.34800888]],\n",
"\n",
" [[-0.36322088, 0.10435939],\n",
" [-0.03852791, 0.15730841],\n",
" [-0.37585379, 0.19827358],\n",
" [-0.22372705, 0.11302068],\n",
" [ 0.34064066, 0.38286242],\n",
" [ 0.04932938, 0.04347563],\n",
" [-0.23293536, -0.09850435]],\n",
"\n",
" [[ 0.12204388, -0.17168604],\n",
" [ 0.02803196, -0.16610741],\n",
" [-0.73502209, -0.05792304],\n",
" [-0.22800673, 0.08277428],\n",
" [ 0.56930488, 0.17278178],\n",
" [ 0.22730077, 0.03863085],\n",
" [ 0.16645379, 0.55584833]],\n",
"\n",
" [[ 0.27870564, -0.1417862 ],\n",
" [-0.0016648 , -0.30653317],\n",
" [-0.65416888, -0.18487159],\n",
" [-0.16745256, 0.41772158],\n",
" [-0.05795884, -0.47801255],\n",
" [ 0.24239019, -0.23091666],\n",
" [-0.36323672, 0.18109307]]])Tensor (shape=(7 , 7 , 2 ), inds=[_5ea058AAAEU , _5ea058AAAEV , k4 ], tags={I4 }), backend=numpy , dtype=float64 , data=array([[[ 0.78731107, 0.21004647],\n",
" [ 0.52741506, 0.43994766],\n",
" [ 0.47708148, 0.01481966],\n",
" [ 0.30408889, 0.51773454],\n",
" [ 0.15537385, 0.30171181],\n",
" [-0.48200314, -0.41248418],\n",
" [ 0.65291394, -0.33286303]],\n",
"\n",
" [[-0.11552407, -0.3085487 ],\n",
" [-0.47943249, -0.09614137],\n",
" [ 0.50069197, 0.22250983],\n",
" [-0.54797332, -0.18607586],\n",
" [ 0.38843852, 0.00265275],\n",
" [-0.11764076, 0.31642717],\n",
" [ 0.17207371, -0.85904319]],\n",
"\n",
" [[ 0.19963637, 0.18677644],\n",
" [ 0.18561212, -0.1731359 ],\n",
" [ 0.46093437, -0.13234629],\n",
" [-0.00962448, 0.39038536],\n",
" [ 0.48169 , -0.09446604],\n",
" [-0.37886595, 0.27405837],\n",
" [-0.15396321, 0.10296841]],\n",
"\n",
" [[ 0.29467235, 0.49771329],\n",
" [ 0.24879117, 0.79315666],\n",
" [ 0.0874336 , 0.07108681],\n",
" [-0.15998795, -0.21133163],\n",
" [-0.21104931, 0.36667451],\n",
" [-0.36025673, -0.06683243],\n",
" [-0.09955251, -0.10086256]],\n",
"\n",
" [[-0.14837517, 0.14994806],\n",
" [ 0.14236737, -0.06020421],\n",
" [ 0.25331839, 0.03563144],\n",
" [ 0.09036226, -0.40487843],\n",
" [-0.20131316, 0.01037212],\n",
" [-0.38364995, 0.3672459 ],\n",
" [ 0.18759547, -0.29545806]],\n",
"\n",
" [[ 0.5049154 , -0.22958862],\n",
" [-0.1597303 , -0.0887779 ],\n",
" [ 0.28666437, -0.59732172],\n",
" [-0.11373574, -0.14742915],\n",
" [ 0.026962 , 0.15264385],\n",
" [-0.30225079, -0.28475259],\n",
" [-0.09374053, -0.16452018]],\n",
"\n",
" [[-0.46108207, -0.29431788],\n",
" [ 0.50594233, -0.11105374],\n",
" [-0.25389343, -0.08446522],\n",
" [ 0.08727862, 0.24602351],\n",
" [-0.06263019, -0.02132593],\n",
" [ 0.68480669, -0.03278971],\n",
" [ 0.1320918 , 0.00572759]]])Tensor (shape=(7 , 7 , 2 ), inds=[_5ea058AAAEV , _5ea058AAAEW , k5 ], tags={I5 }), backend=numpy , dtype=float64 , data=array([[[ 0.51364797, -0.05129313],\n",
" [ 0.70325937, 0.02386915],\n",
" [ 0.19114355, 0.23573459],\n",
" [-0.63226395, -0.62909883],\n",
" [-0.22932703, 0.09187213],\n",
" [-0.41967212, -0.02675523],\n",
" [-0.29576288, -0.2168272 ]],\n",
"\n",
" [[-0.45869116, -0.01779355],\n",
" [-0.04927913, 0.01325764],\n",
" [-0.03131577, 0.49120286],\n",
" [ 0.19169085, -0.14233112],\n",
" [-0.34506342, -0.17719807],\n",
" [-0.41280227, 0.43697613],\n",
" [ 0.28949735, -0.48692094]],\n",
"\n",
" [[ 0.35971934, 0.085507 ],\n",
" [-0.17927927, 0.16269724],\n",
" [ 0.26949982, 0.24089774],\n",
" [ 0.04723488, -0.39356784],\n",
" [ 0.18284115, 0.20338727],\n",
" [ 0.05865349, -0.60845679],\n",
" [ 0.23929981, -0.1059546 ]],\n",
"\n",
" [[-0.35244058, -0.08371223],\n",
" [-0.00801144, 0.44759691],\n",
" [-0.59055561, 0.39347605],\n",
" [ 0.0668733 , 0.22434765],\n",
" [-0.01507118, -0.3575622 ],\n",
" [ 0.1431948 , 0.51901979],\n",
" [-0.44053523, 0.3831614 ]],\n",
"\n",
" [[-0.30060194, 0.0931552 ],\n",
" [-0.14095805, 0.17144856],\n",
" [ 0.52046393, -0.21912552],\n",
" [-0.17362052, -0.8160126 ],\n",
" [-0.44542313, 0.59560646],\n",
" [ 0.18510886, 0.20923147],\n",
" [-0.31669808, 0.25502686]],\n",
"\n",
" [[-0.18306577, -0.20746615],\n",
" [ 0.21120363, 0.36578189],\n",
" [ 0.26170018, -0.0115031 ],\n",
" [ 0.26300372, -0.55516702],\n",
" [-0.01230597, -0.04557538],\n",
" [-0.04333503, 0.22461592],\n",
" [-0.17056938, -0.31543471]],\n",
"\n",
" [[ 0.22401276, 0.02752952],\n",
" [ 0.29523743, 0.38584941],\n",
" [ 0.35414564, -0.0815388 ],\n",
" [-0.03486628, 0.10008203],\n",
" [-0.60031227, -0.43144617],\n",
" [ 0.38524697, -0.09778537],\n",
" [ 0.18044684, 0.56252633]]])Tensor (shape=(7 , 7 , 2 ), inds=[_5ea058AAAEW , _5ea058AAAEX , k6 ], tags={I6 }), backend=numpy , dtype=float64 , data=array([[[-0.08416254, -0.0579689 ],\n",
" [-0.26615288, -0.09438523],\n",
" [ 0.55961413, -0.12349958],\n",
" [-0.38895929, -0.01365775],\n",
" [-0.29974783, 0.50690448],\n",
" [-0.17886188, -0.14730136],\n",
" [-0.23385431, -0.026047 ]],\n",
"\n",
" [[ 0.51090189, 0.35497452],\n",
" [ 0.18397158, 0.39974937],\n",
" [ 0.16056089, 0.20187852],\n",
" [-0.0386433 , 0.42650666],\n",
" [-0.13475823, -0.27543706],\n",
" [ 0.15131959, 0.02987652],\n",
" [ 0.07220134, -0.56702063]],\n",
"\n",
" [[-0.39611837, 0.24513859],\n",
" [-0.20262517, 0.56382248],\n",
" [-0.49713261, -0.10003621],\n",
" [-0.18759944, -0.04383945],\n",
" [-0.06965295, 0.28456265],\n",
" [ 0.27040769, -0.39774807],\n",
" [ 0.00507573, -0.23567965]],\n",
"\n",
" [[ 0.00394807, -0.56985675],\n",
" [ 0.09711578, -0.05933225],\n",
" [-0.21970798, -0.04312317],\n",
" [ 0.2233734 , 0.44252441],\n",
" [-0.08452594, -0.1732618 ],\n",
" [ 0.42290529, -0.12223644],\n",
" [ 0.32487263, 0.00483605]],\n",
"\n",
" [[ 0.03777547, -0.08431554],\n",
" [ 0.03941421, 0.59094796],\n",
" [ 0.18934892, -0.36456678],\n",
" [ 0.15828623, 0.29132 ],\n",
" [ 0.50695886, 0.12483893],\n",
" [ 0.37394813, 0.0600786 ],\n",
" [-0.05031781, -0.07473872]],\n",
"\n",
" [[ 0.44370141, 0.13628731],\n",
" [-0.15349218, -0.26250782],\n",
" [ 0.37029134, 0.41577255],\n",
" [ 0.0776724 , -0.06259455],\n",
" [-0.25977525, -0.49513776],\n",
" [ 0.00827536, 0.67746206],\n",
" [ 0.21727701, -0.24625169]],\n",
"\n",
" [[ 0.16761199, 0.02220682],\n",
" [-0.09161784, -0.02836317],\n",
" [-0.66472973, -0.15587899],\n",
" [ 0.58247699, 0.2643906 ],\n",
" [-0.01236895, 0.21687767],\n",
" [-0.25117893, 0.24708616],\n",
" [ 0.71025628, 0.10578184]]])Tensor (shape=(7 , 7 , 2 ), inds=[_5ea058AAAEX , _5ea058AAAEY , k7 ], tags={I7 }), backend=numpy , dtype=float64 , data=array([[[-0.45844823, 0.19821437],\n",
" [-0.50070115, 0.17094377],\n",
" [-0.67319816, -0.07647333],\n",
" [ 0.03047687, 0.02406828],\n",
" [-0.62353805, 0.10828566],\n",
" [-0.09422612, 0.48445926],\n",
" [ 0.36529513, 0.40705443]],\n",
"\n",
" [[-0.04076505, -0.03983574],\n",
" [ 0.37114129, 0.02757671],\n",
" [ 0.05976666, -0.01888297],\n",
" [ 0.16732854, -0.19853326],\n",
" [-0.37967514, 0.03365298],\n",
" [-0.49694625, -0.66504023],\n",
" [ 0.70764006, 0.17127784]],\n",
"\n",
" [[ 0.16653726, -0.18978569],\n",
" [ 0.20946685, 0.03558293],\n",
" [-0.19452723, -0.10302914],\n",
" [-0.53704991, 0.24556788],\n",
" [ 0.30772713, 0.16027986],\n",
" [ 0.85879266, 0.19407688],\n",
" [-0.53258289, 0.05346453]],\n",
"\n",
" [[ 0.40837715, -0.29533384],\n",
" [-0.47594543, 0.21333849],\n",
" [-0.1050386 , -0.05675576],\n",
" [ 0.22193326, -0.41177346],\n",
" [-0.17850445, 0.23798179],\n",
" [-0.36880912, -0.0121511 ],\n",
" [-0.12499594, -0.11726721]],\n",
"\n",
" [[-0.02470055, -0.05922319],\n",
" [-0.05955615, -0.12797192],\n",
" [-0.15261432, 0.35397322],\n",
" [ 0.16668137, 0.07480621],\n",
" [-0.00882383, -0.02154783],\n",
" [-0.34210173, -0.16146632],\n",
" [-0.4732369 , -0.66080432]],\n",
"\n",
" [[-0.13566065, -0.04801924],\n",
" [-0.26504744, 0.30688566],\n",
" [ 0.0988583 , 0.1705232 ],\n",
" [ 0.10661639, 0.47083813],\n",
" [ 0.11864655, 0.04496947],\n",
" [ 0.26307149, -0.69095818],\n",
" [-0.47671824, 0.04503324]],\n",
"\n",
" [[-0.25374946, 0.19250075],\n",
" [ 0.47149054, 0.21589819],\n",
" [-0.36058823, 0.22073824],\n",
" [ 0.2440586 , -0.18564929],\n",
" [ 0.14559286, 0.10587352],\n",
" [-0.06487019, 0.15311495],\n",
" [-0.21291694, -0.45274407]]])Tensor (shape=(7 , 7 , 2 ), inds=[_5ea058AAAEY , _5ea058AAAEZ , k8 ], tags={I8 }), backend=numpy , dtype=float64 , data=array([[[ 0.24456271, -0.12176326],\n",
" [-0.40440777, -0.3713822 ],\n",
" [-0.09944227, -0.01234654],\n",
" [-0.75625724, 0.51645599],\n",
" [ 0.18954741, -0.35448119],\n",
" [ 0.33923677, -0.31658588],\n",
" [-0.62842783, -0.03817052]],\n",
"\n",
" [[-0.37128572, 0.12138385],\n",
" [ 0.22812109, 0.64689099],\n",
" [ 0.22224472, -0.17281071],\n",
" [ 0.35256084, -0.36033104],\n",
" [-0.02964014, -0.11188638],\n",
" [ 0.42826982, -0.39809753],\n",
" [ 0.16748572, -0.15635467]],\n",
"\n",
" [[ 0.31508692, 0.05269544],\n",
" [-0.05012618, -0.33509913],\n",
" [ 0.02628795, 0.27904 ],\n",
" [ 0.55128451, -0.49773928],\n",
" [-0.28395441, 0.16990847],\n",
" [-0.00684591, -0.18512134],\n",
" [-0.1172634 , -0.05983457]],\n",
"\n",
" [[ 0.1411343 , -0.69501957],\n",
" [ 0.04219888, 0.06855935],\n",
" [-0.41487571, 0.06884335],\n",
" [ 0.40358933, -0.49723671],\n",
" [ 0.11392077, -0.45377457],\n",
" [-0.30389272, -0.0832783 ],\n",
" [-0.30842307, 0.00122607]],\n",
"\n",
" [[-0.04746192, -0.12647542],\n",
" [-0.3377732 , -0.00438588],\n",
" [-0.02402572, 0.3029833 ],\n",
" [ 0.13637175, -0.23409117],\n",
" [ 0.36040963, 0.59882607],\n",
" [ 0.14990102, -0.03343051],\n",
" [ 0.10663151, 0.37216732]],\n",
"\n",
" [[ 0.14140998, 0.20763977],\n",
" [ 0.23547196, -0.27723627],\n",
" [-0.21712285, -0.40360787],\n",
" [-0.05885538, -0.05081516],\n",
" [-0.00518586, 0.18117953],\n",
" [-0.54896781, 0.34085905],\n",
" [ 0.23078799, 0.24785405]],\n",
"\n",
" [[-0.43450751, -0.01882217],\n",
" [ 0.40313653, -0.11985614],\n",
" [-0.13927992, -0.21346225],\n",
" [ 0.17328442, 0.44976657],\n",
" [ 0.46067941, -0.2740861 ],\n",
" [ 0.21052843, -0.09491916],\n",
" [-0.35267016, -0.15636201]]])Tensor (shape=(7 , 2 ), inds=[_5ea058AAAEZ , k9 ], tags={I9 }), backend=numpy , dtype=float64 , data=array([[ 0.16846815, 0.11286849],\n",
" [ 0.1167139 , 0.02678773],\n",
" [-0.17737708, -0.06013245],\n",
" [ 0.22978633, 0.09973961],\n",
" [-0.11355793, 0.14300286],\n",
" [ 0.20092254, -0.30521252],\n",
" [ 0.01617367, -0.04109249]]) "
],
"text/plain": [
"TensorNetwork1D(tensors=20, indices=28, L=10, max_bond=7)"
]
},
"execution_count": 44,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"# so if we combine two MPS, which are both TensorNetwork1D\n",
"mps_a = qtn.MPS_rand_state(10, 7)\n",
"mps_b = qtn.MPS_rand_state(10, 7)\n",
"\n",
"# we get a TensorNetwork1D\n",
"mps_a | mps_b"
]
},
{
"cell_type": "code",
"execution_count": 45,
"metadata": {},
"outputs": [
{
"data": {
"text/html": [
"TensorNetwork2D (tensors=200, indices=460, Lx=10, Ly=10, max_bond=7)Tensor (shape=(7 , 7 , 2 ), inds=[_5ea058AAAEj , _5ea058AAAEk , k0,0 ], tags={I0,0 , X0 , Y0 }), backend=numpy , dtype=float64 , data=array([[[-0.18945712, 0.04733577],\n",
" [ 0.11395101, 0.25020231],\n",
" [-0.27039268, 0.33628383],\n",
" [ 0.22535564, -0.18309576],\n",
" [-0.00395838, 0.01527437],\n",
" [-0.10016121, -0.18655138],\n",
" [-0.0214583 , -0.27461605]],\n",
"\n",
" [[ 0.0644115 , -0.271691 ],\n",
" [ 0.1673902 , 0.01083355],\n",
" [ 0.15142456, -0.00455654],\n",
" [-0.04770221, 0.03159731],\n",
" [ 0.00251603, -0.34378662],\n",
" [-0.24126223, -0.00207497],\n",
" [-0.04420883, -0.13983361]],\n",
"\n",
" [[ 0.50898858, 0.08123947],\n",
" [ 0.0193988 , -0.01183026],\n",
" [ 0.0199097 , -0.18194533],\n",
" [-0.16061977, -0.08799774],\n",
" [ 0.14128153, 0.1963858 ],\n",
" [-0.16963567, -0.0523913 ],\n",
" [-0.1930432 , -0.37175708]],\n",
"\n",
" [[ 0.08123199, -0.19606467],\n",
" [-0.12662249, 0.17052092],\n",
" [ 0.15104046, 0.07922823],\n",
" [-0.0964309 , -0.09962798],\n",
" [ 0.04545444, 0.04060369],\n",
" [ 0.14834529, 0.02477108],\n",
" [ 0.23559065, 0.23742286]],\n",
"\n",
" [[ 0.04460572, -0.01688814],\n",
" [-0.06462653, -0.10666646],\n",
" [ 0.1485165 , 0.06126698],\n",
" [ 0.24760386, 0.26277296],\n",
" [ 0.43826333, 0.30841239],\n",
" [ 0.00659603, 0.08740246],\n",
" [-0.05959811, 0.23595577]],\n",
"\n",
" [[ 0.01975747, -0.07942149],\n",
" [ 0.01864525, 0.33346222],\n",
" [ 0.00312937, 0.40834082],\n",
" [-0.00192261, -0.35181204],\n",
" [ 0.24373095, 0.24885491],\n",
" [-0.07843895, 0.15079061],\n",
" [ 0.14989711, 0.03337135]],\n",
"\n",
" [[ 0.22640118, 0.03562621],\n",
" [ 0.10399615, -0.08323904],\n",
" [ 0.25796734, -0.02694456],\n",
" [ 0.04500069, -0.13867799],\n",
" [-0.14243468, 0.11007008],\n",
" [ 0.11964526, 0.1383061 ],\n",
" [ 0.11297899, -0.16352407]]])Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAEl , _5ea058AAAEm , _5ea058AAAEk , k0,1 ], tags={I0,1 , X0 , Y1 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAEn , _5ea058AAAEo , _5ea058AAAEm , k0,2 ], tags={I0,2 , X0 , Y2 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAEp , _5ea058AAAEq , _5ea058AAAEo , k0,3 ], tags={I0,3 , X0 , Y3 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAEr , _5ea058AAAEs , _5ea058AAAEq , k0,4 ], tags={I0,4 , X0 , Y4 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAEt , _5ea058AAAEu , _5ea058AAAEs , k0,5 ], tags={I0,5 , X0 , Y5 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAEv , _5ea058AAAEw , _5ea058AAAEu , k0,6 ], tags={I0,6 , X0 , Y6 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAEx , _5ea058AAAEy , _5ea058AAAEw , k0,7 ], tags={I0,7 , X0 , Y7 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAEz , _5ea058AAAFA , _5ea058AAAEy , k0,8 ], tags={I0,8 , X0 , Y8 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 2 ), inds=[_5ea058AAAFB , _5ea058AAAFA , k0,9 ], tags={I0,9 , X0 , Y9 }), backend=numpy , dtype=float64 , data=array([[[-0.12157857, 0.24161 ],\n",
" [ 0.41460762, -0.04849642],\n",
" [-0.0332495 , -0.01508166],\n",
" [-0.04748053, -0.22953934],\n",
" [ 0.26398502, 0.23571598],\n",
" [ 0.09435663, 0.31150425],\n",
" [ 0.28419376, 0.23574838]],\n",
"\n",
" [[ 0.21980332, 0.24134622],\n",
" [ 0.16117748, -0.10802498],\n",
" [ 0.51308873, -0.37409859],\n",
" [-0.19663492, 0.04822279],\n",
" [-0.0327776 , 0.12235521],\n",
" [ 0.03224348, 0.28776351],\n",
" [ 0.20878521, -0.08138795]],\n",
"\n",
" [[ 0.03165915, -0.12784148],\n",
" [-0.02866177, -0.01234716],\n",
" [ 0.23303737, -0.03226831],\n",
" [-0.0441102 , -0.16063379],\n",
" [ 0.16779209, -0.09493676],\n",
" [ 0.11110184, -0.13365511],\n",
" [-0.15680939, 0.19855325]],\n",
"\n",
" [[ 0.02178109, -0.27466432],\n",
" [-0.003194 , 0.21448861],\n",
" [ 0.03580574, 0.24084592],\n",
" [ 0.18790656, -0.36077226],\n",
" [ 0.07536634, -0.03645931],\n",
" [ 0.11691943, 0.25219886],\n",
" [ 0.12164802, 0.11102753]],\n",
"\n",
" [[ 0.13786039, 0.10199639],\n",
" [ 0.30040412, -0.36796953],\n",
" [ 0.08601876, -0.00957148],\n",
" [ 0.11099973, -0.15113817],\n",
" [ 0.22138195, 0.12668388],\n",
" [-0.11017555, 0.00325563],\n",
" [ 0.18984643, -0.08675919]],\n",
"\n",
" [[ 0.24661761, 0.18865507],\n",
" [-0.11735149, 0.3072796 ],\n",
" [ 0.13419097, -0.23414407],\n",
" [ 0.0371796 , 0.07015309],\n",
" [-0.34676174, 0.12863882],\n",
" [ 0.02824097, -0.17235414],\n",
" [-0.23060837, 0.15910656]],\n",
"\n",
" [[ 0.17518182, -0.09890955],\n",
" [ 0.16654388, -0.03732843],\n",
" [-0.10735401, 0.06406169],\n",
" [ 0.17888228, 0.11146607],\n",
" [ 0.11045499, -0.05856184],\n",
" [ 0.27557675, 0.09962027],\n",
" [-0.17436828, -0.16866183]]])Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAFC , _5ea058AAAFD , _5ea058AAAEj , k1,0 ], tags={I1,0 , X1 , Y0 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAFE , _5ea058AAAFF , _5ea058AAAEl , _5ea058AAAFD , k1,1 ], tags={I1,1 , X1 , Y1 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAFG , _5ea058AAAFH , _5ea058AAAEn , _5ea058AAAFF , k1,2 ], tags={I1,2 , X1 , Y2 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAFI , _5ea058AAAFJ , _5ea058AAAEp , _5ea058AAAFH , k1,3 ], tags={I1,3 , X1 , Y3 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAFK , _5ea058AAAFL , _5ea058AAAEr , _5ea058AAAFJ , k1,4 ], tags={I1,4 , X1 , Y4 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAFM , _5ea058AAAFN , _5ea058AAAEt , _5ea058AAAFL , k1,5 ], tags={I1,5 , X1 , Y5 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAFO , _5ea058AAAFP , _5ea058AAAEv , _5ea058AAAFN , k1,6 ], tags={I1,6 , X1 , Y6 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAFQ , _5ea058AAAFR , _5ea058AAAEx , _5ea058AAAFP , k1,7 ], tags={I1,7 , X1 , Y7 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAFS , _5ea058AAAFT , _5ea058AAAEz , _5ea058AAAFR , k1,8 ], tags={I1,8 , X1 , Y8 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAFU , _5ea058AAAFB , _5ea058AAAFT , k1,9 ], tags={I1,9 , X1 , Y9 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAFV , _5ea058AAAFW , _5ea058AAAFC , k2,0 ], tags={I2,0 , X2 , Y0 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAFX , _5ea058AAAFY , _5ea058AAAFE , _5ea058AAAFW , k2,1 ], tags={I2,1 , X2 , Y1 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAFZ , _5ea058AAAFa , _5ea058AAAFG , _5ea058AAAFY , k2,2 ], tags={I2,2 , X2 , Y2 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAFb , _5ea058AAAFc , _5ea058AAAFI , _5ea058AAAFa , k2,3 ], tags={I2,3 , X2 , Y3 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAFd , _5ea058AAAFe , _5ea058AAAFK , _5ea058AAAFc , k2,4 ], tags={I2,4 , X2 , Y4 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAFf , _5ea058AAAFg , _5ea058AAAFM , _5ea058AAAFe , k2,5 ], tags={I2,5 , X2 , Y5 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAFh , _5ea058AAAFi , _5ea058AAAFO , _5ea058AAAFg , k2,6 ], tags={I2,6 , X2 , Y6 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAFj , _5ea058AAAFk , _5ea058AAAFQ , _5ea058AAAFi , k2,7 ], tags={I2,7 , X2 , Y7 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAFl , _5ea058AAAFm , _5ea058AAAFS , _5ea058AAAFk , k2,8 ], tags={I2,8 , X2 , Y8 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAFn , _5ea058AAAFU , _5ea058AAAFm , k2,9 ], tags={I2,9 , X2 , Y9 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAFo , _5ea058AAAFp , _5ea058AAAFV , k3,0 ], tags={I3,0 , X3 , Y0 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAFq , _5ea058AAAFr , _5ea058AAAFX , _5ea058AAAFp , k3,1 ], tags={I3,1 , X3 , Y1 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAFs , _5ea058AAAFt , _5ea058AAAFZ , _5ea058AAAFr , k3,2 ], tags={I3,2 , X3 , Y2 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAFu , _5ea058AAAFv , _5ea058AAAFb , _5ea058AAAFt , k3,3 ], tags={I3,3 , X3 , Y3 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAFw , _5ea058AAAFx , _5ea058AAAFd , _5ea058AAAFv , k3,4 ], tags={I3,4 , X3 , Y4 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAFy , _5ea058AAAFz , _5ea058AAAFf , _5ea058AAAFx , k3,5 ], tags={I3,5 , X3 , Y5 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGA , _5ea058AAAGB , _5ea058AAAFh , _5ea058AAAFz , k3,6 ], tags={I3,6 , X3 , Y6 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGC , _5ea058AAAGD , _5ea058AAAFj , _5ea058AAAGB , k3,7 ], tags={I3,7 , X3 , Y7 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGE , _5ea058AAAGF , _5ea058AAAFl , _5ea058AAAGD , k3,8 ], tags={I3,8 , X3 , Y8 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAGG , _5ea058AAAFn , _5ea058AAAGF , k3,9 ], tags={I3,9 , X3 , Y9 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAGH , _5ea058AAAGI , _5ea058AAAFo , k4,0 ], tags={I4,0 , X4 , Y0 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGJ , _5ea058AAAGK , _5ea058AAAFq , _5ea058AAAGI , k4,1 ], tags={I4,1 , X4 , Y1 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGL , _5ea058AAAGM , _5ea058AAAFs , _5ea058AAAGK , k4,2 ], tags={I4,2 , X4 , Y2 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGN , _5ea058AAAGO , _5ea058AAAFu , _5ea058AAAGM , k4,3 ], tags={I4,3 , X4 , Y3 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGP , _5ea058AAAGQ , _5ea058AAAFw , _5ea058AAAGO , k4,4 ], tags={I4,4 , X4 , Y4 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGR , _5ea058AAAGS , _5ea058AAAFy , _5ea058AAAGQ , k4,5 ], tags={I4,5 , X4 , Y5 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGT , _5ea058AAAGU , _5ea058AAAGA , _5ea058AAAGS , k4,6 ], tags={I4,6 , X4 , Y6 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGV , _5ea058AAAGW , _5ea058AAAGC , _5ea058AAAGU , k4,7 ], tags={I4,7 , X4 , Y7 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGX , _5ea058AAAGY , _5ea058AAAGE , _5ea058AAAGW , k4,8 ], tags={I4,8 , X4 , Y8 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAGZ , _5ea058AAAGG , _5ea058AAAGY , k4,9 ], tags={I4,9 , X4 , Y9 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAGa , _5ea058AAAGb , _5ea058AAAGH , k5,0 ], tags={I5,0 , X5 , Y0 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGc , _5ea058AAAGd , _5ea058AAAGJ , _5ea058AAAGb , k5,1 ], tags={I5,1 , X5 , Y1 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGe , _5ea058AAAGf , _5ea058AAAGL , _5ea058AAAGd , k5,2 ], tags={I5,2 , X5 , Y2 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGg , _5ea058AAAGh , _5ea058AAAGN , _5ea058AAAGf , k5,3 ], tags={I5,3 , X5 , Y3 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGi , _5ea058AAAGj , _5ea058AAAGP , _5ea058AAAGh , k5,4 ], tags={I5,4 , X5 , Y4 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGk , _5ea058AAAGl , _5ea058AAAGR , _5ea058AAAGj , k5,5 ], tags={I5,5 , X5 , Y5 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGm , _5ea058AAAGn , _5ea058AAAGT , _5ea058AAAGl , k5,6 ], tags={I5,6 , X5 , Y6 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGo , _5ea058AAAGp , _5ea058AAAGV , _5ea058AAAGn , k5,7 ], tags={I5,7 , X5 , Y7 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGq , _5ea058AAAGr , _5ea058AAAGX , _5ea058AAAGp , k5,8 ], tags={I5,8 , X5 , Y8 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAGs , _5ea058AAAGZ , _5ea058AAAGr , k5,9 ], tags={I5,9 , X5 , Y9 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAGt , _5ea058AAAGu , _5ea058AAAGa , k6,0 ], tags={I6,0 , X6 , Y0 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGv , _5ea058AAAGw , _5ea058AAAGc , _5ea058AAAGu , k6,1 ], tags={I6,1 , X6 , Y1 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGx , _5ea058AAAGy , _5ea058AAAGe , _5ea058AAAGw , k6,2 ], tags={I6,2 , X6 , Y2 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAGz , _5ea058AAAHA , _5ea058AAAGg , _5ea058AAAGy , k6,3 ], tags={I6,3 , X6 , Y3 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAHB , _5ea058AAAHC , _5ea058AAAGi , _5ea058AAAHA , k6,4 ], tags={I6,4 , X6 , Y4 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAHD , _5ea058AAAHE , _5ea058AAAGk , _5ea058AAAHC , k6,5 ], tags={I6,5 , X6 , Y5 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAHF , _5ea058AAAHG , _5ea058AAAGm , _5ea058AAAHE , k6,6 ], tags={I6,6 , X6 , Y6 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAHH , _5ea058AAAHI , _5ea058AAAGo , _5ea058AAAHG , k6,7 ], tags={I6,7 , X6 , Y7 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAHJ , _5ea058AAAHK , _5ea058AAAGq , _5ea058AAAHI , k6,8 ], tags={I6,8 , X6 , Y8 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAHL , _5ea058AAAGs , _5ea058AAAHK , k6,9 ], tags={I6,9 , X6 , Y9 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAHM , _5ea058AAAHN , _5ea058AAAGt , k7,0 ], tags={I7,0 , X7 , Y0 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAHO , _5ea058AAAHP , _5ea058AAAGv , _5ea058AAAHN , k7,1 ], tags={I7,1 , X7 , Y1 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAHQ , _5ea058AAAHR , _5ea058AAAGx , _5ea058AAAHP , k7,2 ], tags={I7,2 , X7 , Y2 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAHS , _5ea058AAAHT , _5ea058AAAGz , _5ea058AAAHR , k7,3 ], tags={I7,3 , X7 , Y3 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAHU , _5ea058AAAHV , _5ea058AAAHB , _5ea058AAAHT , k7,4 ], tags={I7,4 , X7 , Y4 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAHW , _5ea058AAAHX , _5ea058AAAHD , _5ea058AAAHV , k7,5 ], tags={I7,5 , X7 , Y5 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAHY , _5ea058AAAHZ , _5ea058AAAHF , _5ea058AAAHX , k7,6 ], tags={I7,6 , X7 , Y6 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAHa , _5ea058AAAHb , _5ea058AAAHH , _5ea058AAAHZ , k7,7 ], tags={I7,7 , X7 , Y7 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAHc , _5ea058AAAHd , _5ea058AAAHJ , _5ea058AAAHb , k7,8 ], tags={I7,8 , X7 , Y8 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAHe , _5ea058AAAHL , _5ea058AAAHd , k7,9 ], tags={I7,9 , X7 , Y9 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAHf , _5ea058AAAHg , _5ea058AAAHM , k8,0 ], tags={I8,0 , X8 , Y0 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAHh , _5ea058AAAHi , _5ea058AAAHO , _5ea058AAAHg , k8,1 ], tags={I8,1 , X8 , Y1 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAHj , _5ea058AAAHk , _5ea058AAAHQ , _5ea058AAAHi , k8,2 ], tags={I8,2 , X8 , Y2 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAHl , _5ea058AAAHm , _5ea058AAAHS , _5ea058AAAHk , k8,3 ], tags={I8,3 , X8 , Y3 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAHn , _5ea058AAAHo , _5ea058AAAHU , _5ea058AAAHm , k8,4 ], tags={I8,4 , X8 , Y4 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAHp , _5ea058AAAHq , _5ea058AAAHW , _5ea058AAAHo , k8,5 ], tags={I8,5 , X8 , Y5 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAHr , _5ea058AAAHs , _5ea058AAAHY , _5ea058AAAHq , k8,6 ], tags={I8,6 , X8 , Y6 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAHt , _5ea058AAAHu , _5ea058AAAHa , _5ea058AAAHs , k8,7 ], tags={I8,7 , X8 , Y7 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 7 , 2 ), inds=[_5ea058AAAHv , _5ea058AAAHw , _5ea058AAAHc , _5ea058AAAHu , k8,8 ], tags={I8,8 , X8 , Y8 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAHx , _5ea058AAAHe , _5ea058AAAHw , k8,9 ], tags={I8,9 , X8 , Y9 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 2 ), inds=[_5ea058AAAHy , _5ea058AAAHf , k9,0 ], tags={I9,0 , X9 , Y0 }), backend=numpy , dtype=float64 , data=array([[[ 0.41095166, -0.25378787],\n",
" [-0.35543473, -0.0435038 ],\n",
" [ 0.16586621, -0.11328361],\n",
" [ 0.22856655, -0.06458924],\n",
" [ 0.07386491, -0.13653966],\n",
" [ 0.06506377, 0.07501375],\n",
" [-0.11399144, -0.1518549 ]],\n",
"\n",
" [[ 0.15119689, -0.16860582],\n",
" [ 0.13829012, 0.03127736],\n",
" [-0.08314063, 0.07096206],\n",
" [ 0.16171264, 0.16351128],\n",
" [-0.02334824, -0.32312151],\n",
" [-0.04392508, 0.17779242],\n",
" [-0.10118641, 0.15466673]],\n",
"\n",
" [[-0.20247927, -0.20889139],\n",
" [-0.12521677, -0.1455537 ],\n",
" [-0.06994869, -0.09362385],\n",
" [-0.08286166, 0.30008448],\n",
" [ 0.23543103, 0.17068613],\n",
" [ 0.23195817, 0.06582497],\n",
" [-0.01197463, -0.27970348]],\n",
"\n",
" [[ 0.23034017, -0.31215142],\n",
" [ 0.13855897, -0.26866738],\n",
" [ 0.22677439, 0.02264464],\n",
" [-0.18822234, -0.21160647],\n",
" [ 0.12650527, -0.41123326],\n",
" [ 0.02079681, 0.03439444],\n",
" [-0.09803844, 0.09402988]],\n",
"\n",
" [[-0.25783507, -0.11126753],\n",
" [-0.04608879, -0.13867866],\n",
" [-0.32301813, 0.08634624],\n",
" [ 0.2792831 , 0.01240703],\n",
" [-0.40054515, 0.04032573],\n",
" [-0.05173055, 0.12809015],\n",
" [-0.04522922, -0.12201716]],\n",
"\n",
" [[-0.23811474, 0.21585966],\n",
" [-0.10873568, -0.11782258],\n",
" [-0.12358463, 0.05955648],\n",
" [ 0.18989939, -0.08107662],\n",
" [-0.154403 , 0.01750408],\n",
" [ 0.2021266 , 0.00810342],\n",
" [ 0.05971326, -0.00138899]],\n",
"\n",
" [[-0.07307709, 0.09017588],\n",
" [ 0.24079384, 0.00536199],\n",
" [ 0.2296317 , 0.06026196],\n",
" [ 0.12816055, 0.03475 ],\n",
" [-0.11037976, 0.30677594],\n",
" [ 0.09008634, 0.03713451],\n",
" [ 0.11890037, -0.36858366]]])Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAHz , _5ea058AAAHh , _5ea058AAAHy , k9,1 ], tags={I9,1 , X9 , Y1 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAIA , _5ea058AAAHj , _5ea058AAAHz , k9,2 ], tags={I9,2 , X9 , Y2 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAIB , _5ea058AAAHl , _5ea058AAAIA , k9,3 ], tags={I9,3 , X9 , Y3 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAIC , _5ea058AAAHn , _5ea058AAAIB , k9,4 ], tags={I9,4 , X9 , Y4 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAID , _5ea058AAAHp , _5ea058AAAIC , k9,5 ], tags={I9,5 , X9 , Y5 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAIE , _5ea058AAAHr , _5ea058AAAID , k9,6 ], tags={I9,6 , X9 , Y6 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAIF , _5ea058AAAHt , _5ea058AAAIE , k9,7 ], tags={I9,7 , X9 , Y7 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 7 , 2 ), inds=[_5ea058AAAIG , _5ea058AAAHv , _5ea058AAAIF , k9,8 ], tags={I9,8 , X9 , Y8 }), backend=numpy , dtype=float64 , data=...Tensor (shape=(7 , 7 , 2 ), inds=[_5ea058AAAHx , _5ea058AAAIG , k9,9 ], tags={I9,9 , X9 , Y9 }), backend=numpy , dtype=float64 , data=array([[[-0.24573353, -0.14727387],\n",
" [-0.05657491, -0.28271307],\n",
" [ 0.04672743, -0.04831343],\n",
" [ 0.13421271, 0.12574853],\n",
" [ 0.10747916, 0.00695116],\n",
" [-0.00882143, 0.10984055],\n",
" [-0.08258746, -0.03237999]],\n",
"\n",
" [[-0.30437505, -0.119507 ],\n",
" [-0.40839816, -0.15348508],\n",
" [-0.07300267, -0.14723577],\n",
" [-0.02308698, 0.19687783],\n",
" [ 0.08437634, 0.10334722],\n",
" [-0.23326305, 0.21261936],\n",
" [-0.29865867, 0.05015615]],\n",
"\n",
" [[ 0.22007881, -0.03879863],\n",
" [-0.00212981, 0.02821524],\n",
" [ 0.12168066, 0.02650084],\n",
" [-0.17446343, 0.19189974],\n",
" [ 0.45476891, 0.14171278],\n",
" [-0.03252183, 0.20127262],\n",
" [-0.12840357, -0.04993555]],\n",
"\n",
" [[-0.06608241, 0.33742109],\n",
" [-0.37238338, -0.18034993],\n",
" [-0.10258263, -0.45572748],\n",
" [-0.09355153, 0.42465491],\n",
" [ 0.02882308, -0.04461744],\n",
" [ 0.15523016, -0.02903826],\n",
" [ 0.29048334, -0.1578078 ]],\n",
"\n",
" [[-0.08972963, 0.07892939],\n",
" [-0.34735261, 0.16064733],\n",
" [ 0.17720014, -0.00359091],\n",
" [ 0.02845276, -0.18494751],\n",
" [ 0.06469324, -0.18093769],\n",
" [-0.04889939, 0.08259947],\n",
" [-0.02221498, -0.17020292]],\n",
"\n",
" [[ 0.16730442, 0.15091053],\n",
" [-0.04144885, -0.0344569 ],\n",
" [-0.05661704, 0.41863354],\n",
" [-0.07263036, 0.10514986],\n",
" [-0.04628764, 0.34267562],\n",
" [ 0.1686543 , -0.19185386],\n",
" [-0.27959411, 0.02731203]],\n",
"\n",
" [[-0.13132225, -0.1291075 ],\n",
" [ 0.19778585, 0.214208 ],\n",
" [ 0.10681338, 0.10651628],\n",
" [-0.15393502, -0.03712542],\n",
" [-0.00672591, 0.15715481],\n",
" [-0.03222432, -0.19906711],\n",
" [-0.10373137, -0.03450763]]])...
"
],
"text/plain": [
"TensorNetwork2D(tensors=200, indices=460, Lx=10, Ly=10, max_bond=7)"
]
},
"execution_count": 45,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"# and if we combine two PEPS, which are both TensorNetwork2D\n",
"peps_a = qtn.PEPS.rand(10, 10, 7)\n",
"peps_b = qtn.PEPS.rand(10, 10, 7)\n",
"\n",
"# we get a TensorNetwork2D\n",
"peps_a | peps_b"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"One key thing this allows is for structured contraction schemes\n",
"to remain avaiable once norm- or expectation- TNs have been formed\n",
"for example, without calling\n",
"{meth}`~quimb.tensor.tensor_core.TensorNetwork.view_like`."
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"## Standard vs. Hyper Indices & Tensor Networks\n",
"\n",
"- In standard summation expressions and tensor networks, if an index appears\n",
" twice it is contracted, and if it appears once, it is an 'outer' index that\n",
" will be retained. Such a setup is assumed to be the default when contracting\n",
" tensors, and an error will be raised if an index is encountered more than\n",
" twice.\n",
"- On the other hand, it is a perfectly valid 'sum-of-products-of-tensor-entries'\n",
" expression to have indices appear an arbitrary number of times.\n",
" It can be very beneficial to do so in fact. `quimb`\n",
" supports this, with the main difference being that you will have to\n",
" explicitly name the `output_inds` for contractions, which can no\n",
" longer be automatically inferred."
]
},
{
"cell_type": "code",
"execution_count": 46,
"metadata": {},
"outputs": [
{
"data": {
"text/html": [
"Tensor (shape=(2 , 2 , 2 ), inds=[a , b , c ], tags={}), backend=numpy , dtype=float64 , data=array([[[ 0.20015531, -0.0434293 ],\n",
" [-0.18780095, 0.06465794]],\n",
"\n",
" [[ 0.17666851, -0.26065482],\n",
" [-0.12846997, 0.18518863]]]) "
],
"text/plain": [
"Tensor(shape=(2, 2, 2), inds=('a', 'b', 'c'), tags=oset([]))"
]
},
"execution_count": 46,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"T1 = qtn.rand_tensor((2, 2), ('a', 'x'))\n",
"T2 = qtn.rand_tensor((2, 2), ('b', 'x'))\n",
"T3 = qtn.rand_tensor((2, 2), ('c', 'x'))\n",
"\n",
"# we'll get an error due to 'x' unless we specify output_inds\n",
"qtn.tensor_contract(T1, T2, T3, output_inds=['a', 'b', 'c'])"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
"The above has performed the summation:\n",
"\n",
"$$\n",
"T[a, b, c] = \\sum_{x} T_1[a, x] T_2[b, x] T_2[c, x]\n",
"$$\n",
"\n",
"I.e. everything that doesn't appear in the `output_inds` is\n",
"summed over and removed. Indices such as 'x' can be thought\n",
"of as hyperedges:"
]
},
{
"cell_type": "code",
"execution_count": 47,
"metadata": {},
"outputs": [
{
"data": {
"image/svg+xml": [
"\n",
"\n",
"\n",
" \n",
" \n",
" \n",
" \n",
" 2023-11-29T16:12:23.123735 \n",
" image/svg+xml \n",
" \n",
" \n",
" Matplotlib v3.7.2, https://matplotlib.org/ \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n",
" \n"
],
"text/plain": [
""
]
},
"metadata": {
"needs_background": "light"
},
"output_type": "display_data"
}
],
"source": [
"(T1 | T2 | T3).draw(highlight_inds=['x'])"
]
},
{
"cell_type": "markdown",
"metadata": {
"raw_mimetype": "text/restructuredtext"
},
"source": [
":::{note}\n",
"Contractions of 'standard' tensor networks can be performed just using\n",
"pairwise 'matrix-multiplication equivalent' contractions,\n",
"i.e. using {func}`~numpy.tensordot`. While 'hyper' tensor networks\n",
"are still contracted using pairwise operations, these might require\n",
"routines equivalent to 'batched-matrix-multiplication' and so\n",
"{func}`~numpy.einsum`.\n",
":::"
]
}
],
"metadata": {
"kernelspec": {
"display_name": "Python [conda env:cupy]",
"language": "python",
"name": "conda-env-cupy-py"
},
"language_info": {
"codemirror_mode": {
"name": "ipython",
"version": 3
},
"file_extension": ".py",
"mimetype": "text/x-python",
"name": "python",
"nbconvert_exporter": "python",
"pygments_lexer": "ipython3",
"version": "3.10.10"
},
"vscode": {
"interpreter": {
"hash": "39c10650315d977fb13868ea1402e99f3e10e9885c2c202e692ae90b8995050d"
}
}
},
"nbformat": 4,
"nbformat_minor": 4
}