10. Exploring circuit sampling with various methods¶
This example shows some typical explorations and methods for simulating a quantum circuit on an arbitrary set of edges.
%config InlineBackend.figure_formats = ['svg']
import tqdm
import numpy as np
import cotengra as ctg
import quimb as qu
import quimb.tensor as qtn
First we define some edges and the gates that go into each layer:
# define some edges
edges = qtn.edges_2d_hexagonal(5, 5)
# get all unique sites, sorted
sites = sorted({site for edge in edges for site in edge})
N = len(sites)
# map each site coordinate to a linear register
sitemap = {
site: i for i, site in enumerate(sites)
}
def gen_layer_gates():
for site in sites:
yield qtn.Gate("RX", params=(0.2,), qubits=(sitemap[site],))
for sitea, siteb in edges:
yield qtn.Gate("RZZ", params=(-0.3,), qubits=(sitemap[sitea], sitemap[siteb]))
for site in sites:
yield qtn.Gate("RZ", params=(0.4,), qubits=(sitemap[site],))
Then create an actual sequence of gates to test with:
depth = 6
gates = []
for _ in range(depth):
gates.extend(gen_layer_gates())
10.1. check exact contraction¶
The first thing we might check is whether the circuit is exactly contractible.
Create our circuit:
circ = qtn.Circuit(N)
circ.apply_gates(gates)
circ
<Circuit(n=50, num_gates=990, gate_opts={'contract': 'auto-split-gate', 'propagate_tags': 'register'})>
Crucially, we create a high qualty contraction optimizer that we can reuse:
opt = ctg.ReusableHyperOptimizer(
minimize="combo", # target both flops and mops
reconf_opts={}, # perform subtree reconfiguration refinement
parallel=True, # optimize in parallel
optlib="cmaes", # an efficient parallel meta-optimizer
max_time="rate:1e8", # don't spend time optimizing cheap contractions
hash_method="b", # most generous cache hits
directory=True, # cache paths to disk
progbar=True, # show live progress
)
# # the random greedy optimizer might be better for targeting higher
# # throughput (i.e. cheap but many) contractions, especially if you have
# # `cotengrust` installed (https://github.com/jcmgray/cotengrust)
# opt = "random-greedy"
The most important cost to check is that of a single amplitude:
rehs = circ.amplitude_rehearse(optimize=opt)
F=9.48 C=9.89 S=24.00 P=25.00: 98%|█████████▊| 125/128 [00:31<00:00, 4.02it/s]
To actually sample unbiasedly, but bounded by this cost, one can use ‘gate by gate’ sampling.
Here there is a number of contractions proportional to circuit depth, divided by group_size
.
group_size
being the maximum number of qubits we want to compute marginals for simultaneously.
First ‘rehearsing’ the costs:
rehs = circ.sample_gate_by_gate_rehearse(group_size=17, optimize=opt)
F=5.15 C=6.94 S=17.00 P=17.03: 0%| | 0/128 [00:00<?, ?it/s]
F=6.00 C=7.22 S=17.00 P=17.58: 0%| | 0/128 [00:00<?, ?it/s]
F=5.13 C=6.94 S=17.00 P=17.01: 0%| | 0/128 [00:00<?, ?it/s]
F=6.93 C=7.58 S=17.00 P=18.04: 0%| | 0/128 [00:00<?, ?it/s]
F=6.92 C=7.54 S=17.00 P=17.59: 0%| | 0/128 [00:00<?, ?it/s]
F=5.44 C=7.24 S=17.00 P=18.00: 0%| | 0/128 [00:00<?, ?it/s]
F=5.16 C=6.95 S=17.00 P=17.02: 0%| | 0/128 [00:00<?, ?it/s]
F=6.30 C=7.29 S=17.00 P=17.59: 0%| | 0/128 [00:00<?, ?it/s]
F=6.76 C=7.40 S=17.00 P=17.33: 0%| | 0/128 [00:00<?, ?it/s]
F=5.81 C=7.04 S=17.00 P=17.18: 0%| | 0/128 [00:00<?, ?it/s]
F=6.01 C=7.11 S=17.00 P=17.11: 0%| | 0/128 [00:00<?, ?it/s]
F=7.67 C=8.13 S=19.00 P=19.62: 2%|▏ | 2/128 [00:00<00:49, 2.55it/s]
F=8.07 C=8.32 S=18.00 P=19.32: 5%|▍ | 6/128 [00:01<00:23, 5.15it/s]
F=7.04 C=7.71 S=17.00 P=18.01: 0%| | 0/128 [00:01<?, ?it/s]
F=7.84 C=8.24 S=19.00 P=20.04: 2%|▏ | 2/128 [00:01<01:03, 1.98it/s]
F=8.53 C=8.77 S=20.00 P=21.00: 15%|█▍ | 19/128 [00:03<00:22, 4.89it/s]
F=8.36 C=8.79 S=21.00 P=21.64: 9%|▊ | 11/128 [00:02<00:24, 4.78it/s]
F=7.54 C=8.07 S=18.00 P=18.81: 0%| | 0/128 [00:02<?, ?it/s]
F=8.66 C=9.09 S=22.00 P=22.70: 20%|██ | 26/128 [00:08<00:31, 3.19it/s]
F=9.09 C=9.44 S=23.00 P=23.59: 33%|███▎ | 42/128 [00:12<00:25, 3.35it/s]
F=8.47 C=9.07 S=20.00 P=21.00: 8%|▊ | 10/128 [00:03<00:44, 2.68it/s]
F=9.19 C=9.49 S=22.00 P=23.17: 51%|█████ | 65/128 [00:23<00:23, 2.73it/s]
F=9.78 C=10.09 S=25.00 P=25.59: 100%|██████████| 128/128 [00:53<00:00, 2.41it/s]
F=9.44 C=9.76 S=24.00 P=24.36: 71%|███████ | 91/128 [00:27<00:11, 3.25it/s]
The we can actually run the sampling algorithm. Simplifications are data dependent, i.e. might produce slightly different TN geometries, so we might see the optimizer running a few more times here:
%%time
rng = np.random.default_rng(42)
for b in circ.sample_gate_by_gate(10, group_size=17, optimize=opt, seed=rng):
print(b)
F=7.69 C=8.22 S=19.00 P=20.09: 1%| | 1/128 [00:00<01:24, 1.50it/s]
F=7.86 C=8.26 S=19.00 P=19.81: 4%|▍ | 5/128 [00:01<00:32, 3.84it/s]
00000000000000000111000111010010111001000101001100
F=7.75 C=8.15 S=18.00 P=19.32: 1%| | 1/128 [00:00<01:23, 1.52it/s]
00000000110000000000010001010100001110000010001001
F=7.66 C=8.10 S=18.00 P=19.59: 2%|▏ | 3/128 [00:00<00:29, 4.24it/s]
F=7.91 C=8.39 S=20.00 P=20.59: 2%|▏ | 3/128 [00:01<00:44, 2.82it/s]
00100000000001110100101011001010110010010000010010
F=7.74 C=8.25 S=19.00 P=20.17: 2%|▏ | 2/128 [00:00<00:42, 2.96it/s]
00000000000000001100000000001000000000010000000100
01100010111100000000000110001000001100010101000101
F=7.72 C=8.25 S=18.00 P=19.05: 1%| | 1/128 [00:00<01:23, 1.51it/s]
00000000010000000011100000000000000000011000000000
00000000000000001001100001100110000000110000000000
10001000000000000010100000100000000000000001000001
F=7.46 C=8.04 S=17.00 P=18.09: 0%| | 0/128 [00:00<?, ?it/s]
10000000110100000001000000000000111100010001011010
00000011000001100000100100000000000000000000000000
CPU times: user 3min 50s, sys: 18.4 s, total: 4min 8s
Wall time: 42.7 s
We can also check ‘qubit-by-qubit’ sampling:
rehs = circ.sample_rehearse(group_size=17, optimize=opt)
F=13.68 C=13.81 S=36.00 P=36.58: 100%|██████████| 128/128 [00:55<00:00, 2.30it/s]
F=12.52 C=12.75 S=32.00 P=33.00: 100%|██████████| 128/128 [01:49<00:00, 1.17it/s]
F=9.60 C=9.88 S=24.00 P=24.64: 88%|████████▊ | 112/128 [00:40<00:05, 2.77it/s]
For shallow circuits the total cost might be better, but here these costs are much higher.
10.1.1. check advanced simplifications¶
One thing worth checking, is if advanced TN simplification methods might reduce the contraction complexity more than the default setting (usually not the case).
# explicitly get the unsimplified TN
tn = circ.amplitude_tn(simplify_sequence="")
tn.draw(tn.site_tags, layout="neato")
# compare with default simplifications
tn_fs = circ.amplitude_tn()
tn_fs.draw(tn.site_tags, layout="neato")
Get a rough estimate of its contraction complexity, the 'random-greedy'
optimizer is very quick and medium quality for such estimates:
qu.log10(tn_fs.contraction_cost(optimize="random-greedy"))
9.886418822840847
Then we can try the more advanced method (this is also what you would use to prepare the TN for compressed contraction, thus the name):
tn_cs = tn.compress_simplify(output_inds=(), progbar=True)
S 636, 296: : 10it [00:00, 60.87it/s]
L 680, 1020: : 6it [00:00, 9.89it/s]
S 565, 293: : 10it [00:00, 33.57it/s]
L 672, 1008: : 6it [00:00, 9.59it/s]
S 450, 293: : 10it [00:00, 82.32it/s]
L 672, 1008: : 6it [00:00, 9.42it/s]
S 440, 293: : 10it [00:00, 73.06it/s]
L 672, 1008: : 6it [00:00, 9.61it/s]
S 427, 293: : 10it [00:00, 75.72it/s]
L 672, 1008: : 6it [00:00, 9.77it/s]
S 421, 293: : 10it [00:00, 78.35it/s]
L 672, 1008: : 6it [00:00, 9.60it/s]
S 420, 293: : 10it [00:00, 80.25it/s]
tn_cs.draw(tn.site_tags, layout="neato")
qu.log10(tn_cs.contraction_cost(optimize="random-greedy"))
10.599551819716632
So the advanced simplifications have made a smaller TN but not probably not improved the overall contraction cost.
10.2. check MPS simulation¶
If the geometry is close-ish to 1D, the simulating with an MPS might be more efficient.
circ = qtn.CircuitMPS(len(sites))
circ.apply_gates(gates, progbar=True)
max_bond=286, error~=1.39e-07: 100%|##########| 990/990 [00:20<00:00, 48.35it/s]
circ.psi.show()
2 4 8 15 25 43 72 117 190 273 270 267 255 262 252 261 256 268 271 286
... >─>─>─>──>──>──>──>━━━>━━━>━━━>━━━>━━━>━━━>━━━>━━━>━━━>━━━>━━━>━━━>━━━ ...
│ │ │ │ │ │ │ │ │ │ │ │ │ │ │ │ │ │ │ │
...
273 268 257 262 254 262 257 268 271 286 273 269 257 261 252 261 256 2
... >━━━>━━━>━━━>━━━>━━━>━━━>━━━>━━━>━━━>━━━>━━━>━━━>━━━>━━━>━━━>━━━>━━━>━ ...
│ │ │ │ │ │ │ │ │ │ │ │ │ │ │ │ │ │
...
68 271 284 199 122 73 43 25 15 8 4 2
━━>━━━>━━━>━━━>━━━>──>──>──>──>─>─>─●
│ │ │ │ │ │ │ │ │ │ │ │
%%time
rng = np.random.default_rng(42)
for b in circ.sample(10, seed=rng):
print(b)
10100110000111000000101100000011000000000100010000
00100000000000000000000000000100000000000000000001
10011000010000000000000000000001000010000011000000
11000000010100000001000100011011000000001001000100
11000000000000000010000101000000101100001100010000
01000010000100100000000000000100000000001010111100
10001000000000000001100000000010010000010110011000
00000011001000001001100000100001100111010000001000
01101101010001000001000000000111000000001000001100
01100010100010000100000000000000000000000000000000
CPU times: user 6.97 s, sys: 93 ms, total: 7.06 s
Wall time: 1.1 s
For this particular circuit, it is very efficiently handled by an MPS.
10.3. check PEPS simulation¶
If the geometry is not close to 1D we can try simulatig with a PEPS. In particular, if the amount of loop correlations is small, we can try using the very cheap “simple update” style evolution method.
First we setup the PEPS:
# vacuum state
peps = qtn.TN_from_sites_product_state({site: [1.0, 0.0] for site in sitemap.values()})
# create size 1 bonds
for sitea, siteb in edges:
i = sitemap[sitea]
j = sitemap[siteb]
peps[i].new_bond(peps[j])
Next we apply the gates to the PEPS, using the SU gauge:
gauges = {}
for gate in tqdm.tqdm(gates):
peps.gate_simple_(
gate.array,
gate.qubits,
gauges=gauges,
# this is a crucial paramter - the maximum bond dimension
max_bond=16,
cutoff=1e-6,
cutoff_mode="rel",
renorm=False,
)
100%|██████████| 990/990 [00:00<00:00, 2301.38it/s]
Now we have the PEPS. Note if we didn’t have commuting two qubit gates we might want to use
qtn.edge_coloring
to get an effective trotterization.
# the state normalization gives us an idea of fidelity
peps.normalize_simple(gauges)
0.9999970536698275
peps.draw(peps.site_tags, show_inds="bond-size")
# equilibrate gauges
peps.gauge_all_simple_(100, tol=1e-9, gauges=gauges, progbar=True)
max|dS|=1.49e-10, nfact=0.00: : 1it [00:00, 24.01it/s]
TensorNetworkGenVector(tensors=50, indices=115)
Tensor(shape=(2, 2), inds=[_efc85cAAAAA, k0], tags={I0}),
backend=numpy, dtype=complex128, data=array([[-0.395396 +0.77298075j, 0.27280399+0.41441615j], [-0.27434728-0.41339609j, 0.86610505-0.06082173j]])Tensor(shape=(6, 6, 2, 2), inds=[_efc85cAAAAC, _efc85cAAAAB, _efc85cAAAAA, k1], tags={I1}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(6, 6, 2), inds=[_efc85cAAAAD, _efc85cAAAAB, k2], tags={I2}),
backend=numpy, dtype=complex128, data=array([[[ 4.85363805e-01-7.66458608e-01j, -4.31872925e-01-7.14173083e-02j], [ 1.96617114e-01-4.47278227e-01j, 7.72013312e-01-2.77961572e-02j], [ 4.32727285e-02-7.10562170e-02j, 4.91558441e-01+1.15795030e-02j], [-7.00556342e-02-1.20746786e-01j, 4.63242750e-02+2.18135370e-01j], [-1.17768936e-02+3.39919715e-03j, -7.97033431e-02-2.62602305e-02j], [ 1.97933489e-02+1.84425500e-02j, -1.67954277e-02-3.98424961e-02j]], [[ 2.92546004e-02+4.87593029e-01j, -6.74441017e-01+3.77233463e-01j], [-1.65864303e-01+1.84472029e-01j, 1.51169601e+00-7.05191429e-01j], [ 1.18470655e+00-9.56476051e-01j, -2.41440356e+00+9.37418460e-01j], [-2.79081178e+00-1.19290272e+00j, -1.06748243e+00-8.23806130e-01j], [-2.00290135e-01+1.84914622e-01j, 4.85894357e-01-1.24562665e-01j], [ 5.60173414e-01+7.73402936e-02j, 2.18925254e-01+1.35717807e-01j]], [[ 4.40824048e-02-7.04038548e-02j, 4.92187270e-01+2.73841207e-02j], [-1.47203106e+00+2.78106239e-01j, 2.58215803e+00+3.50122012e-01j], [ 4.63956786e+00-1.06741072e+01j, -2.17065360e+01+1.67489601e+01j], [-2.35946256e+01+2.43686621e+00j, -1.36178912e+01-5.30646400e+00j], [ 2.97806679e+01-5.50862510e+01j, -1.25989157e+02+6.66666214e+01j], [-1.18671773e+02+4.90611014e+01j, -6.85356223e+01-1.95958767e+01j]], [[ 5.19876488e-02-1.29418865e-01j, -1.44787044e-01+1.65062232e-01j], [-6.60123011e-01+2.97647916e+00j, -6.16906260e-01+1.17618307e+00j], [-1.70637700e+01-1.72224309e+01j, -4.33131349e+00-1.37065574e+01j], [-1.62130085e+01+4.59910646e+01j, 5.74141375e+00-5.45284138e+00j], [-1.14743867e+02-7.24993293e+01j, -5.57668322e+01-7.82890404e+01j], [-2.27742533e+01+2.34189675e+02j, 4.49410924e+01-3.67023920e+01j]], [[ 1.16536309e-02-3.72320908e-03j, 7.95496678e-02+2.33140016e-02j], [-2.60202972e-01+8.09041381e-02j, 4.87484585e-01+9.37081187e-02j], [-2.72118799e+01+5.66209507e+01j, 1.23248812e+02-7.54474011e+01j], [ 1.21301337e+02-4.88015286e+01j, 9.69851344e+01-1.06542598e+00j], [-3.13024687e+02+2.68649711e+03j, 4.26221028e+03-3.67051779e+03j], [ 3.66954975e+03-3.66783484e+03j, 3.22280833e+03-5.04472613e+02j]], [[ 1.12718387e-03-2.73801965e-02j, -1.89146131e-02+3.78266185e-02j], [ 5.54566181e-02+5.71458106e-01j, -7.92493173e-02+2.39691571e-01j], [ 1.16517031e+02+5.59777355e+01j, 3.18416710e+01+6.38141928e+01j], [ 1.50381910e+01-2.34489575e+02j, -4.35889529e+01+4.11115404e+01j], [ 5.24538330e+03+5.79104723e+02j, 2.28507674e+03+2.15448419e+03j], [-3.45266586e+03-7.52794882e+03j, -9.43026933e+02+1.11684794e+03j]]])Tensor(shape=(7, 6, 6, 2), inds=[_efc85cAAAAF, _efc85cAAAAE, _efc85cAAAAD, k3], tags={I3}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(6, 6, 2), inds=[_efc85cAAAAG, _efc85cAAAAE, k4], tags={I4}),
backend=numpy, dtype=complex128, data=array([[[ 4.85386671e-01-7.66503867e-01j, -4.31765571e-01-7.16347286e-02j], [-5.88061317e-02-4.84939897e-01j, 6.50333733e-01-4.16991693e-01j], [ 6.24397543e-03+8.28887035e-02j, -4.12532608e-01+2.70145816e-01j], [-1.39603567e-01+1.15959844e-03j, 2.07805986e-01+7.11011186e-02j], [-5.99115218e-03-1.04199791e-02j, 7.37307531e-03-8.49198389e-02j], [-2.66800136e-02-5.81809532e-03j, 3.37253757e-02+2.49004466e-02j]], [[ 3.84475854e-01+3.01335667e-01j, -1.65024415e-01+7.54707211e-01j], [ 1.49971099e-01+1.97690244e-01j, -4.06945716e-01-1.61786747e+00j], [ 8.09538211e-01+1.25940495e+00j, -6.43677544e-01-2.52004159e+00j], [-2.39258550e-01+3.03788830e+00j, -4.46293394e-01+1.24874318e+00j], [-2.56048076e-01+6.57905195e-02j, 4.93227148e-01+1.21371581e-01j], [-1.84361780e-01+5.41057571e-01j, -1.68890652e-01+1.84379499e-01j]], [[ 7.86602051e-02-2.68718608e-02j, 3.65139086e-01+3.31414160e-01j], [-1.49595505e+00+5.97683447e-02j, 2.49872983e+00+7.21996304e-01j], [-5.14467295e+00+1.02944187e+01j, 2.29636405e+01-1.54646304e+01j], [-2.11308916e+01+1.11190297e+01j, -1.44146443e+01+2.08619597e-02j], [ 3.57599437e+01+4.90802224e+01j, -7.31747620e+00-1.43508640e+02j], [ 1.24616206e+02-3.31165124e+01j, 6.32230278e+01+2.71279978e+01j]], [[-4.07730581e-02-1.33521753e-01j, -9.13386239e-03+2.19443130e-01j], [ 2.55578740e+00+1.65950588e+00j, 8.84585554e-01+9.87950582e-01j], [ 2.07399347e+01-1.18322548e+01j, 1.44055693e+01-5.11788005e-01j], [ 4.84316676e+01+8.19633293e+00j, -6.27787775e+00-4.68849576e+00j], [-4.63657911e+01-1.21501067e+02j, 2.00292312e+00-9.60600528e+01j], [-2.16313962e+02-1.04638189e+02j, 1.92004475e+01+5.57270808e+01j]], [[ 1.03134469e-02+6.17233087e-03j, 3.49973198e-02+7.77237909e-02j], [-2.64283823e-01-6.38912534e-03j, 4.41569632e-01+2.51042616e-01j], [ 3.84402247e+01-4.70093660e+01j, -1.36992810e+02+4.33749660e+01j], [ 1.08335861e+02-7.19383579e+01j, 9.41514030e+01-1.91612616e+01j], [-1.77717014e+03-1.86828945e+03j, 1.69588120e+02+5.62276794e+03j], [-4.74569667e+03+1.95811049e+03j, -2.91104755e+03-6.50183912e+02j]], [[-1.84450120e-02-2.01360469e-02j, 1.31995913e-02+3.97887236e-02j], [ 5.64535506e-01+8.96363620e-02j, 2.10419148e-01+1.35060996e-01j], [-7.37345446e+01+1.05779127e+02j, -6.57845171e+01+2.01321436e+01j], [-2.34623746e+02+5.18967517e+01j, 4.95993265e+01+3.18422077e+01j], [ 3.81580739e+03+3.43445496e+03j, 6.70928815e+02+2.90654012e+03j], [ 8.48227957e+03-7.93149229e+02j, -6.37360117e+02-1.35546475e+03j]]])Tensor(shape=(7, 6, 6, 2), inds=[_efc85cAAAAI, _efc85cAAAAH, _efc85cAAAAG, k5], tags={I5}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(6, 6, 2), inds=[_efc85cAAAAJ, _efc85cAAAAH, k6], tags={I6}),
backend=numpy, dtype=complex128, data=array([[[ 4.85386484e-01-7.66503507e-01j, -4.31766877e-01-7.16328910e-02j], [ 5.10921752e-02-4.85812488e-01j, 7.27093844e-01-2.61055163e-01j], [ 4.23069977e-02+7.15549613e-02j, -2.50236501e-01+4.24907827e-01j], [-1.09252665e-01-8.69153834e-02j, 1.16814137e-01+1.85990722e-01j], [-4.99851230e-03-1.09287233e-02j, 1.52084365e-02-8.38546339e-02j], [-2.55038381e-02-9.74709080e-03j, 2.96107908e-02+2.96653677e-02j]], [[ 4.01913636e-01+2.77651839e-01j, -1.19304283e-01+7.63268180e-01j], [ 1.15417565e-01+2.19675175e-01j, -1.35288446e-01-1.66277952e+00j], [ 1.32763537e+00+6.92007710e-01j, -1.80893480e+00-1.86892548e+00j], [-1.96012931e+00+2.33326574e+00j, -1.08889508e+00+7.56899597e-01j], [-2.58026755e-01+5.73962240e-02j, 4.88923722e-01+1.37358800e-01j], [-2.32279094e-01+5.22197457e-01j, -1.84768632e-01+1.68406427e-01j]], [[ 8.02609518e-02-2.16431374e-02j, 3.42573598e-01+3.54690564e-01j], [-1.45017554e+00-3.72277935e-01j, 2.18620187e+00+1.40900200e+00j], [-8.10657779e-01+1.14801534e+01j, 1.52917793e+01-2.30791822e+01j], [-2.30653013e+01-6.17989458e+00j, -1.05989748e+01-9.77004375e+00j], [ 2.75623532e+01+5.41234932e+01j, 1.54501361e+01-1.42895840e+02j], [ 1.28840322e+02-5.62109295e+00j, 5.59468409e+01+4.00665654e+01j]], [[ 1.88290703e-02-1.38331894e-01j, -1.00104891e-01+1.95501666e-01j], [ 1.00988711e+00+2.87506084e+00j, 9.69167888e-02+1.32259979e+00j], [ 2.04127023e+01-1.23886720e+01j, 1.43863075e+01-9.00858995e-01j], [ 1.40962144e+01+4.70547315e+01j, 1.42459589e+00-7.70470580e+00j], [ 2.07491130e+01-1.28409746e+02j, 4.98447997e+01-8.21202821e+01j], [-1.23163537e+02-2.06387890e+02j, -1.46123924e+01+5.70987136e+01j]], [[ 1.00447230e-02+6.61020431e-03j, 3.16755155e-02+7.91263431e-02j], [-2.53274590e-01-7.62124641e-02j, 3.59372886e-01+3.58894175e-01j], [ 1.61626069e+01-5.85615486e+01j, -1.07733298e+02+9.50663987e+01j], [ 1.28884463e+02+1.77580358e+01j, 8.32882645e+01+4.79089295e+01j], [-1.50911600e+03-2.09215850e+03j, -5.90522094e+02+5.59330463e+03j], [-5.03519064e+03+1.01265876e+03j, -2.73356555e+03-1.19698934e+03j]], [[-1.47760548e-02-2.29618426e-02j, 6.30151559e-03+4.14639264e-02j], [ 4.86458268e-01+3.00065045e-01j, 1.42313206e-01+2.05723184e-01j], [-4.05369854e+01+1.22377120e+02j, -5.73403445e+01+3.80699832e+01j], [-1.93859229e+02-1.42071418e+02j, 8.85488521e+00+5.83003083e+01j], [ 2.79349118e+03+4.30731354e+03j, -1.06675740e+02+2.98020416e+03j], [ 8.30333324e+03+1.91641264e+03j, -1.79645645e+02-1.48948260e+03j]]])Tensor(shape=(7, 6, 6, 2), inds=[_efc85cAAAAL, _efc85cAAAAK, _efc85cAAAAJ, k7], tags={I7}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(6, 6, 2), inds=[_efc85cAAAAM, _efc85cAAAAK, k8], tags={I8}),
backend=numpy, dtype=complex128, data=array([[[ 8.75280245e-01-2.68410601e-01j, -2.49599652e-01-3.41992243e-01j], [ 4.09914572e-01-2.75444759e-01j, 6.27203889e-01+4.08395122e-01j], [-3.47382095e-02+6.98402729e-02j, -5.40897730e-01+4.59759400e-02j], [ 2.67609755e-02-1.49701340e-01j, -1.11226060e-01+2.46193860e-01j], [ 6.87010385e-03-3.61192582e-03j, 5.86191803e-02+3.35894590e-03j], [ 6.67549387e-03-1.72050356e-02j, -1.87048202e-02+2.28174123e-02j]], [[ 7.14275443e-02-4.69326897e-01j, 7.65032920e-01-1.83006517e-01j], [ 1.53025630e-01-1.24360947e-01j, -1.56688794e+00+5.09168458e-01j], [ 4.65861456e-01-1.38044658e+00j, -1.27022283e+00+1.84243411e+00j], [ 2.25452671e+00+1.53540856e+00j, 8.60276742e-01+8.82338663e-01j], [-1.77393330e-02+1.45635555e-01j, 1.42341048e-01-2.27585060e-01j], [ 2.16865717e-01+2.36033225e-01j, 5.63492025e-02+1.22101373e-01j]], [[ 3.67786835e-02+5.90920810e-02j, -2.27782744e-01+4.23137325e-01j], [ 2.37833823e-01-1.52106321e+00j, -1.07980951e+00+2.38101934e+00j], [-1.16543836e+01+1.07331456e+00j, 2.46758545e+01+9.48260571e+00j], [ 3.40593809e+00-2.22132877e+01j, 8.14381589e+00-1.24280037e+01j], [-1.02943850e+02+9.53807421e+00j, 2.00175643e+02+8.86948427e+01j], [-8.08412202e+01+1.98180349e+02j, -8.70568861e+01+6.26218796e+01j]], [[-1.18976251e-01+7.34606864e-02j, 2.16835080e-01-5.97575597e-02j], [ 1.95214197e+00-2.32715254e+00j, 1.06779242e+00-8.29018828e-01j], [-2.10011624e+01-8.61035740e+00j, -1.12376486e+01-1.08524847e+01j], [ 2.81369941e+01-3.30886790e+01j, -7.34887208e+00+4.85802205e+00j], [-1.49162025e+02+2.37865458e+01j, -2.51217197e+02-4.14446144e+01j], [-2.23340796e+02+1.87789521e+02j, 1.03597381e+02-5.00895518e+01j]], [[ 3.46009782e-03-8.96673845e-03j, 6.90557513e-02-5.22512390e-02j], [-5.48191573e-02+2.14085195e-01j, 3.25058597e-01-4.04389165e-01j], [-5.64368878e+01-1.48959026e+01j, 1.12793994e+02+8.64520987e+01j], [ 1.39964210e+01-1.08879891e+02j, 4.01810686e+01-9.11455469e+01j], [-3.75010930e+03+4.23672107e+02j, 7.81312701e+03+3.46334965e+03j], [-2.04825810e+03+6.92801620e+03j, -3.28771230e+03+3.28102159e+03j]], [[-1.54205561e-02+2.26011016e-02j, 2.94449894e-02-2.46754957e-02j], [ 1.82684265e-01-5.46878640e-01j, 1.25855350e-01-1.93232135e-01j], [ 1.25848252e+02-1.02294279e+00j, 5.69861860e+01+2.92788638e+01j], [-7.27117849e+01+2.29366916e+02j, 4.20097031e+01-3.61757496e+01j], [ 6.90915058e+03-3.08194553e+03j, 4.58433724e+03+3.86305595e+02j], [ 5.88208993e+03-1.30306783e+04j, -2.01884049e+03+7.27676682e+02j]]])Tensor(shape=(6, 6, 2), inds=[_efc85cAAAAN, _efc85cAAAAM, k9], tags={I9}),
backend=numpy, dtype=complex128, data=array([[[ 4.88011981e-01-7.74599293e-01j, -4.12309598e-01-9.62261695e-02j], [-2.06645918e-01+4.27395712e-01j, -7.84980232e-01-5.07226826e-02j], [-5.22580407e-03+6.94063151e-02j, -4.33732721e-01+2.06896352e-01j], [-8.92259185e-02-1.07659445e-01j, 8.92037193e-02+2.06473035e-01j], [ 5.36877010e-03-7.97310288e-03j, 7.89503870e-02-3.55741366e-02j], [-4.49905924e-03-2.69878035e-02j, -4.19596592e-03+3.81903423e-02j]], [[-2.85939587e-01+4.02664034e-01j, -7.31131124e-01-1.60058034e-01j], [ 1.97074414e-01-6.64558875e-03j, -1.55621130e+00-5.40922088e-01j], [-1.39653814e+00+6.47981061e-01j, 2.57944949e+00-4.26241719e-01j], [-1.86791462e+00-2.39528541e+00j, -5.87619849e-01-1.21743820e+00j], [ 1.62253742e-01-1.50066566e-01j, -4.92504216e-01+1.63175337e-01j], [-5.30565083e-02-5.74130423e-01j, 3.73850553e-02-2.27570233e-01j]], [[ 6.38464991e-02-4.48110643e-02j, 4.97321677e-01+2.17612709e-01j], [ 1.40709841e+00+3.77801172e-01j, -1.91893465e+00-1.15139885e+00j], [-2.06710223e+00+1.15197111e+01j, 1.57078586e+01-2.12621915e+01j], [-2.25547676e+01+2.54352247e+00j, -1.50706528e+01-4.11540996e+00j], [-6.98142711e+00+5.79536995e+01j, 7.01026962e+01-1.23618625e+02j], [ 1.14064590e+02+5.31765674e+01j, 3.88638476e+01+5.09394912e+01j]], [[ 1.34876465e-01-7.02494137e-02j, -2.62869107e-01+6.23095975e-02j], [ 2.22062569e+00-1.58404283e+00j, 1.12549195e+00-5.01861997e-01j], [ 1.50334698e+01+1.67040517e+01j, 4.46016243e+00+1.41733525e+01j], [-2.52286268e+01+3.53562917e+01j, 6.91027253e+00-5.46392017e+00j], [ 1.00852541e+02+4.33752446e+01j, 7.65946256e+01+6.36779363e+01j], [ 2.21337353e+02-9.43587620e+01j, -5.48200786e+01-8.28379927e+00j]], [[ 7.30293708e-03+2.62882930e-03j, 3.71021489e-02+4.55074722e-02j], [ 1.35675275e-01+5.58246695e-02j, -1.81544602e-01-1.97731304e-01j], [ 5.04062007e+01-9.02636369e+01j, -1.86424196e+02+1.14818127e+02j], [ 1.43030915e+02-4.85495953e+01j, 2.54609159e+02-1.42366094e+00j], [ 3.07653362e+02-3.76164807e+03j, -4.90607907e+03+6.99767530e+03j], [-7.08818392e+03-2.64450492e+03j, -2.99866670e+03-3.48866984e+03j]], [[ 1.59256706e-02-9.32466260e-03j, -2.88364991e-02+6.24100358e-03j], [ 2.58519972e-01-1.89499654e-01j, 1.27572255e-01-4.25332199e-02j], [-1.44380618e+02-1.58003271e+02j, -1.92485832e+01-1.05497613e+02j], [ 1.51020945e+02-2.49678627e+02j, -8.19631026e+01+8.07678066e+01j], [-6.85121028e+03-2.29359934e+03j, -3.16181266e+03-3.40254405e+03j], [-1.23446103e+04+7.21115603e+03j, 2.06706982e+03+5.77319728e+02j]]])Tensor(shape=(6, 6, 2), inds=[_efc85cAAAAO, _efc85cAAAAC, k10], tags={I10}),
backend=numpy, dtype=complex128, data=array([[[ 4.85363805e-01-7.66458608e-01j, -4.31872925e-01-7.14173083e-02j], [ 9.11875151e-02-4.80000978e-01j, 7.46068065e-01-2.00398659e-01j], [ 7.98599076e-02+2.33218814e-02j, 1.16552281e-01+4.77681223e-01j], [-1.38369080e-01-1.84817770e-02j, 2.00872356e-01+9.68466544e-02j], [ 9.26891328e-03-8.02104503e-03j, 8.33724702e-02-9.55261353e-03j], [-2.56886492e-02-8.48513134e-03j, 3.19066102e-02+2.91801397e-02j]], [[-4.88414971e-01-7.32185715e-03j, -3.25696130e-01-7.00783705e-01j], [-2.07392453e-01-1.36122250e-01j, 9.25270838e-01+1.38794584e+00j], [-9.85122818e-01+1.16099487e+00j, 2.19672358e+00-1.37204354e+00j], [-1.27471262e+00-2.75440647e+00j, -2.41551051e-01-1.32658523e+00j], [ 2.43304037e-01+1.22933978e-01j, -2.86692133e-01-4.11602669e-01j], [-1.22824954e-01-5.51987083e-01j, 5.08701639e-02-2.52507157e-01j]], [[ 7.77644500e-02+2.92001969e-02j, 7.10471844e-02+4.87801698e-01j], [-8.62317675e-01-1.22500072e+00j, 7.49782139e-01+2.49558655e+00j], [ 6.17290706e-01+1.16224397e+01j, 1.19384852e+01-2.46814493e+01j], [-2.22379698e+01-8.25332498e+00j, -9.86495158e+00-1.07837030e+01j], [-6.20346206e+01+8.54935207e+00j, 1.28254250e+02+6.21976679e+01j], [ 1.09642670e+02+6.68465700e+01j, 2.46284586e+01+6.68922461e+01j]], [[ 1.14781031e-01-7.92292441e-02j, -2.11838509e-01+5.77345169e-02j], [-1.66184645e+00+2.55606279e+00j, -9.90142675e-01+8.85210383e-01j], [ 2.17350993e+01-1.07410347e+01j, 1.43728547e+01+2.25923261e-01j], [-5.50790994e-01+4.87620374e+01j, 3.68130673e+00-7.01037027e+00j], [ 1.02415958e+02+8.90690533e+01j, 4.32431715e+01+8.58436836e+01j], [ 5.83368187e+01-2.27947948e+02j, -5.00274304e+01+2.93943347e+01j]], [[ 1.21919065e-02-1.01331595e-03j, 7.22897075e-02+4.05695769e-02j], [-2.60194313e-01+8.09319707e-02j, 4.87494608e-01+9.36559684e-02j], [-5.75776561e+01-2.51243295e+01j, 7.99014683e+01+1.20409034e+02j], [ 6.29510565e+01-1.14598287e+02j, 7.43724596e+01-6.22574337e+01j], [-2.43209793e+02-2.69371627e+03j, -3.42113620e+03+4.46485785e+03j], [-2.84717060e+03+4.33730630e+03j, -3.05292614e+03+1.14918588e+03j]], [[ 1.29289957e-02+2.41616796e-02j, -2.90755657e-03-4.21919635e-02j], [-4.33473662e-01-3.76484241e-01j, -1.07435181e-01-2.28451714e-01j], [-4.40195974e+01-1.21540117e+02j, 2.19336666e+01-6.78605742e+01j], [ 2.30739654e+02+4.43927757e+01j, -3.53066332e+01-4.84107146e+01j], [ 3.46352808e+03-3.98163346e+03j, 3.07236106e+03-6.51131252e+02j], [-8.15222320e+03-1.46018652e+03j, 3.79065823e+02+1.41172286e+03j]]])Tensor(shape=(6, 7, 6, 2), inds=[_efc85cAAAAQ, _efc85cAAAAP, _efc85cAAAAO, k11], tags={I11}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(7, 7, 7, 2), inds=[_efc85cAAAAR, _efc85cAAAAP, _efc85cAAAAF, k12], tags={I12}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(7, 7, 7, 2), inds=[_efc85cAAAAT, _efc85cAAAAS, _efc85cAAAAR, k13], tags={I13}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(7, 7, 7, 2), inds=[_efc85cAAAAU, _efc85cAAAAS, _efc85cAAAAI, k14], tags={I14}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(7, 7, 7, 2), inds=[_efc85cAAAAW, _efc85cAAAAV, _efc85cAAAAU, k15], tags={I15}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(7, 7, 7, 2), inds=[_efc85cAAAAX, _efc85cAAAAV, _efc85cAAAAL, k16], tags={I16}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(7, 7, 7, 2), inds=[_efc85cAAAAZ, _efc85cAAAAY, _efc85cAAAAX, k17], tags={I17}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(6, 7, 6, 2), inds=[_efc85cAAAAa, _efc85cAAAAY, _efc85cAAAAN, k18], tags={I18}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(6, 6, 2), inds=[_efc85cAAAAb, _efc85cAAAAa, k19], tags={I19}),
backend=numpy, dtype=complex128, data=array([[[ 4.85386484e-01-7.66503507e-01j, -4.31766877e-01-7.16328910e-02j], [ 4.75603210e-01-1.11475115e-01j, 4.85596804e-01+6.00839055e-01j], [ 1.35891676e-02+8.20096354e-02j, -3.86869666e-01+3.05767768e-01j], [-1.18813498e-01-7.33048404e-02j, 1.38115004e-01+1.70780956e-01j], [ 6.60244682e-03-1.00500311e-02j, 7.91013866e-02-3.17370924e-02j], [-1.67465605e-02-2.15667824e-02j, 9.93958402e-03+4.07458475e-02j]], [[-1.98973787e-01+4.46131837e-01j, -7.72191330e-01+2.31474872e-02j], [-2.09318455e-01-1.33282337e-01j, 9.43655989e-01+1.37573692e+00j], [-9.24104304e-01+1.17797758e+00j, 2.16930958e+00-1.43487324e+00j], [-2.81324465e+00-1.17111118e+00j, -1.04836514e+00-8.12153658e-01j], [ 1.50705848e-01-2.17361389e-01j, -4.52333124e-01+2.30969914e-01j], [-4.06732377e-01-4.01556189e-01j, -9.25291704e-02-2.32411913e-01j]], [[ 5.63353000e-02-6.11255145e-02j, 4.78979465e-01+1.17233269e-01j], [-9.01857312e-01-1.19504968e+00j, 8.31155976e-01+2.46461098e+00j], [ 1.40645091e+00+1.14224771e+01j, 1.05809168e+01-2.55838106e+01j], [-2.18011885e+01+9.73990413e+00j, -1.43860830e+01-9.04504229e-01j], [ 1.43916184e+00+6.07347133e+01j, 7.56872153e+01-1.22128657e+02j], [ 1.28677929e+02-7.82988979e+00j, 5.66489091e+01+3.90933147e+01j]], [[ 9.63471897e-02-1.01032904e-01j, -1.95701642e-01+9.96943304e-02j], [-3.04404982e+00+1.41409275e-01j, -1.29466454e+00-2.87111906e-01j], [ 2.14952871e+01+1.03996110e+01j, 8.57664738e+00+1.15858976e+01j], [-1.05328173e+01+4.79782324e+01j, 4.98747596e+00-6.04293377e+00j], [ 1.26801498e+02+2.91264934e+01j, 7.87075648e+01+5.51110338e+01j], [ 1.42301959e+02-1.93688547e+02j, -5.83121564e+01+8.78360454e+00j]], [[ 1.13220347e-02+4.02857835e-03j, 4.95566759e-02+6.93330846e-02j], [-7.37682342e-02-2.53827693e-01j, -1.05829295e-01+4.96704708e-01j], [ 2.35952743e+01-5.59661181e+01j, -1.19182368e+02+8.03343789e+01j], [ 1.26874205e+02-2.86723263e+01j, 9.47897129e+01+1.55964429e+01j], [ 3.73105821e+02-2.55251752e+03j, -4.31776520e+03+3.60425665e+03j], [-5.10450057e+03-5.48385421e+02j, -2.24548372e+03-1.96205407e+03j]], [[-9.33116524e-03+2.56590244e-02j, 2.91794521e-02-3.00888554e-02j], [ 5.18903112e-01-2.39552512e-01j, 2.49949834e-01-4.83520542e-03j], [ 1.28863548e+02+5.07013115e+00j, 5.24352299e+01+4.45620594e+01j], [ 3.27688353e+01+2.38100739e+02j, 3.57514216e+01-4.68561948e+01j], [ 5.01718748e+03-1.09834599e+03j, 2.75372227e+03+1.15033606e+03j], [ 1.03806740e+03-8.45815211e+03j, -1.46255821e+03+3.34337095e+02j]]])Tensor(shape=(6, 6, 2), inds=[_efc85cAAAAc, _efc85cAAAAQ, k20], tags={I20}),
backend=numpy, dtype=complex128, data=array([[[ 4.85386671e-01-7.66503867e-01j, -4.31765571e-01-7.16347286e-02j], [ 1.02095427e-01-4.77704290e-01j, 7.50591749e-01-1.82833427e-01j], [ 7.48026610e-02+3.62503283e-02j, 2.62830172e-02+4.92413564e-01j], [-1.37190392e-01+2.58707741e-02j, 2.17114557e-01+3.31657309e-02j], [-9.59929939e-03+7.23348726e-03j, -8.51798753e-02+3.18278112e-03j], [-2.54251313e-02-9.96172730e-03j, 2.93659936e-02+2.99177529e-02j]], [[-4.12173672e-01+2.62178813e-01j, -6.56671654e-01-4.06937799e-01j], [-2.48079954e-01-5.39439737e-03j, 1.52032337e+00+6.86816449e-01j], [-4.34796703e-01+1.43262164e+00j, 1.50250510e+00-2.12306590e+00j], [-2.95993906e+00-7.24410628e-01j, -1.16077847e+00-6.41194709e-01j], [-2.64194448e-01+9.50318830e-03j, 4.55853705e-01+2.24056867e-01j], [-3.87184002e-01-4.20500998e-01j, -8.12587803e-02-2.36467489e-01j]], [[ 6.17377861e-02+5.56594155e-02j, -1.11977483e-01+4.80232175e-01j], [-3.21813740e-01-1.46215241e+00j, -2.71497565e-01+2.58673897e+00j], [-5.33434361e+00+1.01974259e+01j, 2.32459601e+01-1.50369185e+01j], [-2.19892316e+01-9.30704541e+00j, -9.14727294e+00-1.11404580e+01j], [ 5.78454717e+01+1.84808895e+01j, -9.06713160e+01-1.11476400e+02j], [ 7.37573063e+01+1.05762780e+02j, -3.85757803e+00+6.86891485e+01j]], [[ 1.12930808e-01-8.20800421e-02j, -2.10429012e-01+6.29153858e-02j], [-1.65721390e+00+2.55727415e+00j, -9.87157139e-01+8.85470913e-01j], [ 2.26344277e+01-7.60459255e+00j, 1.42332608e+01+2.27961446e+00j], [ 1.64994717e+01+4.62663348e+01j, 1.02554052e+00-7.76801185e+00j], [-9.51395609e+01-8.86609470e+01j, -4.05097500e+01-8.71235007e+01j], [ 6.57454184e+01-2.31124137e+02j, -5.15784762e+01+2.85277877e+01j]], [[ 1.20188095e-02+1.14344270e-04j, 6.94772021e-02+4.93833831e-02j], [-2.60351645e-01+4.58670010e-02j, 4.82412746e-01+1.59009685e-01j], [-4.71372409e+01-3.82833113e+01j, 4.38313746e+01+1.36847464e+02j], [ 3.55332218e+01-1.25096669e+02j, 5.90728682e+01-7.57762277e+01j], [-4.02277514e+02+2.54696095e+03j, 4.36030595e+03-3.55415219e+03j], [-3.71267337e+03+3.54569020e+03j, -2.94867981e+03+4.49639529e+02j]], [[-9.24796651e-03-2.56934616e-02j, -3.17373916e-03+4.18007074e-02j], [ 3.62085340e-01+4.42299922e-01j, 6.73748275e-02+2.40786885e-01j], [ 5.09157482e+00+1.28841308e+02j, -4.02840199e+01+5.57683037e+01j], [-2.40294753e+02-7.94528111e-02j, 4.15379444e+01+4.18164872e+01j], [ 4.33373178e+03-2.75220639e+03j, 2.97991333e+03+1.35031343e+02j], [ 7.62774554e+03+3.79415171e+03j, 1.72286886e+02-1.48790515e+03j]]])Tensor(shape=(6, 7, 6, 2), inds=[_efc85cAAAAe, _efc85cAAAAd, _efc85cAAAAc, k21], tags={I21}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(7, 7, 7, 2), inds=[_efc85cAAAAf, _efc85cAAAAd, _efc85cAAAAT, k22], tags={I22}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(7, 7, 7, 2), inds=[_efc85cAAAAh, _efc85cAAAAg, _efc85cAAAAf, k23], tags={I23}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(7, 7, 7, 2), inds=[_efc85cAAAAi, _efc85cAAAAg, _efc85cAAAAW, k24], tags={I24}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(7, 7, 7, 2), inds=[_efc85cAAAAk, _efc85cAAAAj, _efc85cAAAAi, k25], tags={I25}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(7, 7, 7, 2), inds=[_efc85cAAAAl, _efc85cAAAAj, _efc85cAAAAZ, k26], tags={I26}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(7, 7, 7, 2), inds=[_efc85cAAAAn, _efc85cAAAAm, _efc85cAAAAl, k27], tags={I27}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(6, 7, 6, 2), inds=[_efc85cAAAAo, _efc85cAAAAm, _efc85cAAAAb, k28], tags={I28}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(6, 6, 2), inds=[_efc85cAAAAp, _efc85cAAAAo, k29], tags={I29}),
backend=numpy, dtype=complex128, data=array([[[ 4.85386671e-01-7.66503867e-01j, -4.31765571e-01-7.16347286e-02j], [ 4.76944330e-01-1.05588692e-01j, 4.78134023e-01+6.06798062e-01j], [ 1.45768338e-02+8.18354488e-02j, -3.83164733e-01+3.10397625e-01j], [-1.16634833e-01-7.67255922e-02j, 1.33087911e-01+1.74717835e-01j], [ 6.51976011e-03-1.00976138e-02j, 7.88693805e-02-3.23322806e-02j], [-1.64768936e-02-2.17758158e-02j, 9.41299367e-03+4.08512043e-02j]], [[-4.72058892e-01+1.25639434e-01j, -5.03337012e-01-5.86061335e-01j], [-2.48398261e-02-2.46892175e-01j, -4.96669904e-01+1.59261388e+00j], [-1.49662432e+00+3.96156377e-02j, 2.48877948e+00+7.55583316e-01j], [-7.71967180e-01-2.94789344e+00j, 1.13129043e-02-1.32605050e+00j], [ 2.63492549e-01-2.14552050e-02j, -4.65531334e-01-2.03187484e-01j], [ 6.38083256e-02-5.68032932e-01j, 1.25331445e-01-2.16359892e-01j]], [[ 5.66509163e-02-6.08292535e-02j, 4.78357759e-01+1.19732070e-01j], [-8.80691478e-01-1.21071692e+00j, 7.87672327e-01+2.47881074e+00j], [ 1.48499544e+00+1.14121648e+01j, 1.04044967e+01-2.56559939e+01j], [-2.21210544e+01+8.98921803e+00j, -1.43469612e+01-1.39541324e+00j], [ 1.59607975e+00+6.07048295e+01j, 7.54053202e+01-1.22320620e+02j], [ 1.28821411e+02-5.56415810e+00j, 5.59164123e+01+4.00800357e+01j]], [[ 9.62648951e-02-1.01111674e-01j, -1.95620843e-01+9.98558136e-02j], [-3.04544197e+00+1.06269493e-01j, -1.29124699e+00-3.02024096e-01j], [ 2.16269582e+01+1.01203823e+01j, 8.72557102e+00+1.14737395e+01j], [-1.18779603e+01+4.76625677e+01j, 5.15509135e+00-5.90074357e+00j], [ 1.26996792e+02+2.80004820e+01j, 7.91680567e+01+5.44423027e+01j], [ 1.44523759e+02-1.91973604e+02j, -5.83892551e+01+8.05239799e+00j]], [[ 1.13127906e-02+4.06034906e-03j, 4.93734418e-02+6.94841642e-02j], [-6.99341800e-02-2.54943795e-01j, -1.13353319e-01+4.95132977e-01j], [ 2.30679416e+01-5.61731212e+01j, -1.18400018e+02+8.14237833e+01j], [ 1.27687875e+02-2.46511720e+01j, 9.42662242e+01+1.85881185e+01j], [ 3.59890733e+02-2.55329511e+03j, -4.30063264e+03+3.62613080e+03j], [-5.09545538e+03-6.26350127e+02j, -2.21595087e+03-1.99663145e+03j]], [[-9.32651126e-03+2.56650380e-02j, 2.91777998e-02-3.01004308e-02j], [ 5.21840419e-01-2.33274985e-01j, 2.50029471e-01-1.80111774e-03j], [ 1.28894773e+02+3.48345842e+00j, 5.29674932e+01+4.39018937e+01j], [ 2.59409712e+01+2.38890112e+02j, 3.70728751e+01-4.58219245e+01j], [ 5.00614769e+03-1.13764660e+03j, 2.76163853e+03+1.12755594e+03j], [ 1.14204006e+03-8.44238316e+03j, -1.46436299e+03+3.14936463e+02j]]])Tensor(shape=(6, 6, 2), inds=[_efc85cAAAAq, _efc85cAAAAe, k30], tags={I30}),
backend=numpy, dtype=complex128, data=array([[[ 4.85386484e-01-7.66503507e-01j, -4.31766877e-01-7.16328910e-02j], [-3.33126412e-02-4.87354544e-01j, 6.71290834e-01-3.82339997e-01j], [ 7.01640426e-02+4.45757968e-02j, -3.02464573e-02+4.92189109e-01j], [-1.36973617e-01+2.69936487e-02j, 2.17377392e-01+3.13872676e-02j], [-9.66965918e-03+7.13581308e-03j, -8.51907203e-02+2.33180285e-03j], [-2.57550073e-02-9.06262073e-03j, 3.03923513e-02+2.88641346e-02j]], [[-3.25831616e-01+3.63948997e-01j, -7.42299138e-01-2.14018262e-01j], [-2.13962407e-01+1.25692025e-01j, 1.65450628e+00-2.13887071e-01j], [-2.01927342e-01+1.48348094e+00j, 1.14654849e+00-2.33464218e+00j], [-3.04461476e+00+1.28673670e-01j, -1.29345211e+00-2.92525439e-01j], [-2.52526727e-01+7.81178090e-02j, 4.98447080e-01+9.72849126e-02j], [-4.96697143e-01-2.82729004e-01j, -1.49330061e-01-2.00501161e-01j]], [[ 6.53644821e-03-8.28705052e-02j, 4.49060325e-01-2.03732400e-01j], [-6.10992651e-01+1.36685330e+00j, 1.75678597e+00-1.91793334e+00j], [ 1.14892184e+01-6.70035209e-01j, -2.48515346e+01-1.22019832e+01j], [ 6.20125542e+00+2.30595675e+01j, -3.14908520e+00+1.40668150e+01j], [-2.00083927e+01-5.73471891e+01j, -3.45874219e+01+1.39504975e+02j], [ 3.46223094e+01-1.24228518e+02j, 5.54447955e+01-4.07584822e+01j]], [[ 1.14163386e-01-8.03552622e-02j, -2.11365764e-01+5.97193813e-02j], [-9.33580423e-01+2.90073687e+00j, -7.22340674e-01+1.11215422e+00j], [ 2.34293913e+01-4.60665078e+00j, 1.38175068e+01+4.10534994e+00j], [ 1.61707512e+01+4.63827320e+01j, 1.08123461e+00-7.76034126e+00j], [-9.28722218e+01-9.10732823e+01j, -3.82784276e+01-8.81079202e+01j], [ 6.11982026e+01-2.32421908e+02j, -5.10082240e+01+2.95288856e+01j]], [[ 1.20098400e-02+5.95823250e-04j, 6.74593138e-02+5.20918180e-02j], [-2.42232768e-01+1.06205669e-01j, 5.06297996e-01+4.02025089e-02j], [-4.06950561e+01-4.51065105e+01j, 2.22582239e+01+1.41945959e+02j], [ 3.95224563e+01-1.23953734e+02j, 6.14376670e+01-7.38756619e+01j], [-5.29325478e+02+2.52475165e+03j, 4.53118400e+03-3.33198604e+03j], [-3.73112280e+03+3.52950066e+03j, -2.95220809e+03+4.35410747e+02j]], [[-1.03090917e-02-2.52843952e-02j, -1.41836553e-03+4.19160417e-02j], [ 4.82939993e-01+3.05695338e-01j, 1.39913996e-01+2.07362368e-01j], [-4.33873444e+00+1.28843246e+02j, -4.42599567e+01+5.27095344e+01j], [-2.40050428e+02+1.18945403e+01j, 4.36215940e+01+3.96798629e+01j], [ 4.24443267e+03-2.88813757e+03j, 2.98183448e+03+4.06922466e+01j], [ 7.89855274e+03+3.19855931e+03j, 5.68703850e+01-1.49920912e+03j]]])Tensor(shape=(6, 7, 6, 2), inds=[_efc85cAAAAs, _efc85cAAAAr, _efc85cAAAAq, k31], tags={I31}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(7, 7, 7, 2), inds=[_efc85cAAAAt, _efc85cAAAAr, _efc85cAAAAh, k32], tags={I32}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(7, 7, 7, 2), inds=[_efc85cAAAAv, _efc85cAAAAu, _efc85cAAAAt, k33], tags={I33}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(7, 7, 7, 2), inds=[_efc85cAAAAw, _efc85cAAAAu, _efc85cAAAAk, k34], tags={I34}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(7, 7, 7, 2), inds=[_efc85cAAAAy, _efc85cAAAAx, _efc85cAAAAw, k35], tags={I35}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(7, 7, 7, 2), inds=[_efc85cAAAAz, _efc85cAAAAx, _efc85cAAAAn, k36], tags={I36}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(7, 7, 7, 2), inds=[_efc85cAAABB, _efc85cAAABA, _efc85cAAAAz, k37], tags={I37}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(6, 7, 6, 2), inds=[_efc85cAAABC, _efc85cAAABA, _efc85cAAAAp, k38], tags={I38}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(6, 6, 2), inds=[_efc85cAAABD, _efc85cAAABC, k39], tags={I39}),
backend=numpy, dtype=complex128, data=array([[[ 4.85363805e-01-7.66458608e-01j, -4.31872925e-01-7.14173083e-02j], [ 2.57620377e-02-4.87790028e-01j, 7.12583967e-01-2.98998096e-01j], [ 2.10258018e-02-8.03609161e-02j, 4.77875321e-01-1.20968495e-01j], [-1.37921361e-01+2.07281579e-02j, 2.15658884e-01+4.12319450e-02j], [-1.09116584e-02+5.53251673e-03j, -8.22478554e-02-1.03414545e-02j], [-2.30155231e-02-1.48737150e-02j, 2.29834719e-02+3.55023887e-02j]], [[-4.40918323e-01-2.10492600e-01j, -3.66755737e-03-7.72504841e-01j], [-1.58773582e-01-1.90609145e-01j, 4.81067676e-01+1.59721435e+00j], [ 6.63716315e-01+1.34301869e+00j, -3.57162271e-01-2.58119364e+00j], [-8.47404330e-01-2.92866805e+00j, -2.23775463e-02-1.32796053e+00j], [-1.13650519e-01-2.47661607e-01j, -3.03317849e-02+4.95480691e-01j], [ 2.58528876e-01-5.12641013e-01j, 1.95083306e-01-1.60239061e-01j]], [[ 7.58478452e-02-3.41850757e-02j, 3.99061949e-01+2.87251365e-01j], [-1.52175010e+00-5.15042183e-02j, 2.43939754e+00+8.70308910e-01j], [ 7.55882904e+00-8.85021224e+00j, -2.56562761e+01+9.66730897e+00j], [-1.83611592e+01+1.58319979e+01j, -1.39732604e+01+3.37314089e+00j], [ 4.87507231e+01-3.96215349e+01j, -1.43650541e+02+1.57160962e+01j], [ 1.29259333e+02+1.29346403e+00j, 5.58956291e+01+4.42939151e+01j]], [[ 5.50603451e-02-1.28280695e-01j, -1.47361832e-01+1.67372247e-01j], [ 3.72401850e-01+3.01213614e+00j, -1.75223627e-01+1.33696372e+00j], [-2.05287508e+01-1.18833960e+01j, -7.90625707e+00-1.22921331e+01j], [ 3.60464080e+01+3.28435095e+01j, -2.78646951e+00-7.41167316e+00j], [-1.21000462e+02-4.95463780e+01j, -7.02915410e+01-6.68301801e+01j], [-1.01017264e+01-2.34752985e+02j, -3.89486119e+01+4.55330865e+01j]], [[ 1.19450513e-02+2.74988421e-03j, 5.69101918e-02+6.16724555e-02j], [-2.70449257e-01+3.41336080e-02j, 4.67843921e-01+1.80921590e-01j], [-4.03109386e+01+4.79200352e+01j, 1.36873017e+02-3.97943690e+01j], [ 6.97380072e+01-1.16441176e+02j, 7.69354808e+01-5.76204821e+01j], [-1.20859034e+03+2.41961954e+03j, 5.25693248e+03-2.00094585e+03j], [-4.81649608e+03+2.15657642e+03j, -3.06126380e+03-7.00751472e+02j]], [[ 8.90325796e-03-2.55468274e-02j, -2.96964445e-02+3.14258502e-02j], [ 9.28000013e-02+5.57822357e-01j, -6.44465671e-02+2.49384322e-01j], [ 1.15035385e+02+5.70694756e+01j, 3.03460628e+01+6.44989080e+01j], [-1.76608345e+02-1.55478188e+02j, 3.47635008e+00+5.79189913e+01j], [ 5.07337432e+03+1.08595713e+03j, 2.17550726e+03+2.43092251e+03j], [ 2.54361275e+03+7.88167974e+03j, 1.06780138e+03-9.98242208e+02j]]])Tensor(shape=(6, 6, 2), inds=[_efc85cAAABE, _efc85cAAAAs, k40], tags={I40}),
backend=numpy, dtype=complex128, data=array([[[ 8.75280245e-01-2.68410601e-01j, -2.49599652e-01-3.41992243e-01j], [ 4.93788112e-01+8.54826660e-03j, 2.80847824e-01+6.93754852e-01j], [-4.84615478e-02+6.11218889e-02j, -5.38678790e-01-6.71513446e-02j], [-9.60357150e-02-1.17914301e-01j, 1.14102685e-01+2.44873907e-01j], [-7.30578421e-03-2.62104241e-03j, -3.71506444e-02-4.54677957e-02j], [-5.08177953e-03-1.77412216e-02j, -1.10877008e-03+2.94834757e-02j]], [[ 1.64455233e-01-4.45335949e-01j, 7.86212443e-01-2.52338793e-02j], [ 1.95541475e-01+2.54169209e-02j, -1.44782990e+00-7.86243540e-01j], [ 9.79164227e-01-1.07884062e+00j, -1.90117017e+00+1.18050052e+00j], [ 2.72283076e+00-1.62980984e-01j, 1.22019001e+00+1.72435265e-01j], [ 1.34057871e-01-5.96061174e-02j, -2.68382225e-01-5.18799843e-03j], [ 2.97511172e-01+1.19287521e-01j, 1.03538447e-01+8.58124123e-02j]], [[-5.33290588e-02+4.47275875e-02j, -4.50978699e-01-1.65977069e-01j], [ 1.25041007e+00+8.98149700e-01j, -1.63832181e+00-2.03743547e+00j], [-2.83193230e-01-1.17002762e+01j, -1.11286532e+01+2.39785463e+01j], [ 1.48415129e+01-1.68748351e+01j, 1.35581724e+01-6.07889574e+00j], [-5.79415215e+01+8.56223753e+01j, 1.95560946e+02-9.84529309e+01j], [-1.96892032e+02+8.39302872e+01j, -1.05921681e+02-1.67629582e+01j]], [[-1.39129795e-01-1.39542416e-02j, 2.08522467e-01+8.43018301e-02j], [ 2.81027455e+00+1.15275937e+00j, 1.11426595e+00+7.65417704e-01j], [-7.13525578e+00-2.15470462e+01j, 9.22583092e-01-1.55951919e+01j], [ 2.07251108e+01-3.81709430e+01j, -6.19179247e+00+6.26641885e+00j], [ 3.69793994e+01+1.46450143e+02j, -1.89013528e+01+2.53910347e+02j], [-2.24094042e+02+1.86889682e+02j, 1.03798237e+02-4.96724154e+01j]], [[ 5.18360511e-03+8.09351062e-03j, 1.28745733e-03+8.65865538e-02j], [-1.72809262e-02-2.20315639e-01j, -1.76990401e-01+4.87717222e-01j], [ 5.20127118e+01-2.64893958e+01j, -1.41713097e+02+1.06699941e+01j], [ 1.09231112e+02-1.09222153e+01j, 9.78597299e+01+1.85876874e+01j], [-3.68140926e+03+8.30685758e+02j, 8.14462374e+03+2.58936623e+03j], [-6.96030944e+03-1.93569992e+03j, -3.33386203e+03-3.23410081e+03j]], [[-1.26235897e-02+2.42744377e-02j, 2.63018442e-02-2.80017949e-02j], [ 4.18472830e-01-3.96648808e-01j, 2.01960892e-01-1.11310495e-01j], [ 1.25434747e+02+1.02446570e+01j, 5.41370129e+01+3.42616580e+01j], [ 1.49907597e+02+1.88212375e+02j, -6.62042971e+00-5.50424376e+01j], [-7.25007595e+03-2.16128793e+03j, -3.22708475e+03-3.27891633e+03j], [-4.82515576e+03-1.34579060e+04j, -9.46118993e+02+1.92617336e+03j]]])Tensor(shape=(6, 6, 2), inds=[_efc85cAAABF, _efc85cAAABE, k41], tags={I41}),
backend=numpy, dtype=complex128, data=array([[[ 4.88011981e-01-7.74599293e-01j, -4.12309598e-01-9.62261695e-02j], [-1.16366233e-01+4.60248336e-01j, -7.79118391e-01+1.08357223e-01j], [ 6.79804661e-02+1.49399410e-02j, 1.43807605e-01+4.58529767e-01j], [-1.36283188e-01-3.12843944e-02j, 1.96310291e-01+1.09775623e-01j], [-9.60439606e-03+3.86835907e-04j, -7.53810366e-02-4.26190364e-02j], [-1.25794981e-03-2.73313124e-02j, -8.70750337e-03+3.74204201e-02j]], [[-4.35255382e-01+2.33350648e-01j, -5.83519937e-01-4.68695826e-01j], [ 1.92226796e-01+4.39470968e-02j, -1.36625972e+00-9.20719668e-01j], [-2.57328719e-01+1.51788674e+00j, 1.11026791e+00-2.36696997e+00j], [-2.29051120e+00-1.99500769e+00j, -8.09053066e-01-1.08299861e+00j], [-1.75581377e-01-1.34229031e-01j, 2.44921135e-01+4.57384122e-01j], [ 2.70450147e-01-5.09212582e-01j, 1.56049410e-01-1.69807061e-01j]], [[ 2.50863113e-02-7.38585399e-02j, 5.29890706e-01-1.17898184e-01j], [ 1.22517409e+00-7.88421285e-01j, -2.14204830e+00+6.47810593e-01j], [ 1.11362606e+01-3.60004929e+00j, -2.61537012e+01-3.84727488e+00j], [-3.92940384e+00+2.23550181e+01j, -8.20381027e+00+1.32950620e+01j], [ 2.41431798e+01-5.31458242e+01j, -1.04136286e+02+9.67036044e+01j], [ 1.25488098e+02-9.55068150e+00j, 5.88417186e+01+2.53551350e+01j]], [[ 1.25902478e-01-8.52948242e-02j, -2.53958993e-01+9.21275376e-02j], [ 1.61406012e+00-2.19890429e+00j, 9.12165456e-01-8.28584199e-01j], [ 1.70880290e+01-1.45955411e+01j, 1.42841761e+01-4.09136689e+00j], [ 6.43072675e+00+4.29557531e+01j, 1.17158966e+00-8.73118937e+00j], [-3.68152956e+01-1.03427663e+02j, -5.43224143e+00-9.94590723e+01j], [ 2.21706856e+02-9.34872725e+01j, -5.47870665e+01-8.49940180e+00j]], [[ 1.41806251e-03+7.63103815e-03j, -2.05484834e-02+5.50023574e-02j], [ 4.91520641e-02+1.38232600e-01j, 2.57494115e-02-2.67194834e-01j], [ 1.25007787e+01-1.02625726e+02j, -1.29112233e+02+1.76825379e+02j], [ 1.50763314e+02+9.23795452e+00j, 2.36181225e+02+9.51118752e+01j], [-3.32963889e+03-1.77711766e+03j, 4.47187780e+03+7.28281199e+03j], [-4.00582998e+02-7.55482020e+03j, 2.00396417e+03-4.14088593e+03j]], [[ 8.45343415e-03-1.64047471e-02j, -2.10075912e-02+2.07165399e-02j], [ 5.43089078e-02-3.15900657e-01j, 6.22151964e-02-1.19218400e-01j], [-8.45920574e+01+1.96608728e+02j, -8.82304110e+01+6.09561055e+01j], [-1.82208172e+02-2.27918412e+02j, 4.75137954e+01+1.04803962e+02j], [ 5.54598908e+03+4.63052470e+03j, 1.70603716e+03+4.32014955e+03j], [-8.02540904e+03+1.18314323e+04j, 2.11388393e+03-3.70979037e+02j]]])Tensor(shape=(6, 6, 7, 2), inds=[_efc85cAAABG, _efc85cAAABF, _efc85cAAAAv, k42], tags={I42}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(6, 6, 2), inds=[_efc85cAAABH, _efc85cAAABG, k43], tags={I43}),
backend=numpy, dtype=complex128, data=array([[[ 4.85386484e-01-7.66503507e-01j, -4.31766877e-01-7.16328910e-02j], [ 1.02821650e-01-4.77548765e-01j, 7.50865678e-01-1.81693587e-01j], [ 5.76037369e-02+5.99337573e-02j, -1.46004069e-01+4.71004117e-01j], [-1.27003976e-01+5.79675499e-02j, 2.18746829e-01-1.97917686e-02j], [-9.93436389e-03+6.77523775e-03j, -8.52269009e-02-8.04925501e-04j], [-2.72133962e-02-2.23706639e-03j, 3.67459633e-02+2.02176571e-02j]], [[-4.62494944e-01+1.57234254e-01j, -5.41769951e-01-5.50727316e-01j], [-2.39676886e-01-6.42914224e-02j, 1.31315321e+00+1.02896421e+00j], [-1.16098065e+00+9.45369712e-01j, 2.43218625e+00-9.21539648e-01j], [-2.96205852e+00-7.15580976e-01j, -1.16272729e+00-6.37752531e-01j], [-2.56272716e-01-6.54409894e-02j, 3.74217323e-01+3.43385727e-01j], [-4.08721724e-01-3.99531166e-01j, -9.36823265e-02-2.31949496e-01j]], [[ 6.27917685e-02+5.44719043e-02j, -1.02778292e-01+4.82287871e-01j], [-3.47504778e-01-1.45627291e+00j, -2.25930241e-01+2.59115478e+00j], [-8.38973948e+00+7.87803016e+00j, 2.68966957e+01-6.56162522e+00j], [-2.36457547e+01-3.32202294e+00j, -1.17104563e+01-8.40492269e+00j], [ 5.73297425e+01+2.01016751e+01j, -8.75312779e+01-1.13939631e+02j], [ 1.02517071e+02+7.81637179e+01j, 1.73320091e+01+6.66107188e+01j]], [[ 1.14042473e-01-8.05279024e-02j, -2.11268715e-01+6.00301770e-02j], [-1.69590795e+00+2.53182387e+00j, -1.00051837e+00+8.70374931e-01j], [ 2.38562006e+01+1.03969687e+00j, 1.24628029e+01+7.24364866e+00j], [ 2.65308850e+01+4.13396077e+01j, -7.56312186e-01-7.79871508e+00j], [-8.96028283e+01-9.43308345e+01j, -3.51592186e+01-8.94199988e+01j], [-1.41088441e-01-2.40343673e+02j, -4.17676887e+01+4.16283510e+01j]], [[ 1.20025253e-02-5.97761519e-04j, 7.22678626e-02+4.51696441e-02j], [-2.57239661e-01+6.08112950e-02j, 4.90703910e-01+1.30862822e-01j], [-3.39926357e+01-5.03333301e+01j, 2.19530137e+00+1.43712317e+02j], [-3.12923777e+00-1.30036063e+02j, 3.39553973e+01-8.98630543e+01j], [-3.71554569e+02+2.55274345e+03j, 4.31557371e+03-3.60687857e+03j], [-2.26848017e+03+4.60550271e+03j, -2.61401273e+03+1.43484748e+03j]], [[-1.26562307e-02-2.41924773e-02j, 2.54387502e-03+4.18366903e-02j], [ 4.18254109e-01+3.89498814e-01j, 9.91517175e-02+2.29493456e-01j], [-2.28304073e+01+1.26926342e+02j, -5.13809992e+01+4.57736164e+01j], [-2.23400298e+02+8.86454221e+01j, 5.40345243e+01+2.35360620e+01j], [ 4.07215838e+03-3.12987070e+03j, 2.98142895e+03-1.31498967e+02j], [ 8.51680177e+03+2.86262636e+02j, -4.62324275e+02-1.42728602e+03j]]])Tensor(shape=(6, 6, 7, 2), inds=[_efc85cAAABI, _efc85cAAABH, _efc85cAAAAy, k44], tags={I44}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(6, 6, 2), inds=[_efc85cAAABJ, _efc85cAAABI, k45], tags={I45}),
backend=numpy, dtype=complex128, data=array([[[ 4.85386671e-01-7.66503867e-01j, -4.31765571e-01-7.16347286e-02j], [ 1.22165471e-01-4.72969834e-01j, 7.57639275e-01-1.50992587e-01j], [ 6.71850320e-02+4.89458500e-02j, -6.14915724e-02+4.89265454e-01j], [-1.34415879e-01+3.77209792e-02j, 2.19178127e-01+1.41301988e-02j], [-1.14465779e-02+3.66673359e-03j, -8.14722514e-02-2.50604860e-02j], [-2.56805070e-02+9.28362725e-03j, 4.18057035e-02+3.11585588e-03j]], [[-4.87110397e-01+3.67195651e-02j, -3.87023318e-01-6.68602264e-01j], [-2.10860318e-01-1.30807834e-01j, 9.59754174e-01+1.36454091e+00j], [-1.23013783e+00+8.53354959e-01j, 2.49582978e+00-7.31958438e-01j], [-2.43672520e+00-1.82985773e+00j, -8.15185904e-01-1.04594926e+00j], [-1.85865653e-01-1.87996303e-01j, 1.43743286e-01+4.87178048e-01j], [-4.72537162e-01-3.21623297e-01j, -1.32819897e-01-2.11845346e-01j]], [[ 2.74258494e-02+7.84687657e-02j, -3.28831300e-01+3.67466856e-01j], [ 4.80126366e-01-1.41807320e+00j, -1.56924729e+00+2.07422113e+00j], [-1.05564238e+01+4.58308112e+00j, 2.75327526e+01+2.90363736e+00j], [-1.63930824e+01-1.73612770e+01j, -3.89738899e+00-1.38777820e+01j], [ 2.50336288e+01+5.53257748e+01j, 2.20142324e+01-1.41998840e+02j], [ 9.50098329e+01+8.71725164e+01j, 1.11862397e+01+6.78816817e+01j]], [[ 1.14841471e-01-7.93847413e-02j, -2.11858850e-01+5.79182638e-02j], [-1.82197539e+00+2.44262475e+00j, -1.04330992e+00+8.18560788e-01j], [ 2.36999659e+01-2.90844811e+00j, 1.34864145e+01+5.08909727e+00j], [ 1.94039316e+01+4.51253115e+01j, 5.30228696e-01-7.81745476e+00j], [-5.78747817e+01-1.16459103e+02j, -7.27670414e+00-9.58049850e+01j], [-9.78738265e+01-2.19457732e+02j, -2.12557602e+01+5.49758039e+01j]], [[ 1.17988587e-02+2.29185897e-03j, 5.93704943e-02+6.11631442e-02j], [-2.64029068e-01-1.32584451e-02j, 4.34892568e-01+2.62438743e-01j], [-3.05931637e+01-5.24557707e+01j, -7.28565249e+00+1.43510683e+02j], [ 4.72560141e+01-1.21155877e+02j, 6.60039262e+01-6.98220573e+01j], [-1.61081567e+03+2.01348166e+03j, 5.54818995e+03-9.28353985e+02j], [-1.37148282e+03+4.94722342e+03j, -2.30180031e+03+1.89702059e+03j]], [[ 6.82461234e-03+2.64405486e-02j, 7.03776639e-03-4.13261429e-02j], [-2.99205038e-01-4.87043084e-01j, -3.43254678e-02-2.47668636e-01j], [ 2.96158244e+01-1.25494636e+02j, 5.37536338e+01-4.29357541e+01j], [ 2.40289968e+02+1.44794931e+00j, -4.12993160e+01-4.20526157e+01j], [-5.08558178e+03+7.01884630e+02j, -2.65415392e+03-1.36141408e+03j], [-8.40239857e+03+1.40631617e+03j, 7.34002296e+02+1.30568622e+03j]]])Tensor(shape=(6, 6, 7, 2), inds=[_efc85cAAABK, _efc85cAAABJ, _efc85cAAABB, k46], tags={I46}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(6, 6, 2), inds=[_efc85cAAABL, _efc85cAAABK, k47], tags={I47}),
backend=numpy, dtype=complex128, data=array([[[ 4.85363805e-01-7.66458608e-01j, -4.31872925e-01-7.14173083e-02j], [-3.55678754e-01-3.34806538e-01j, 2.32153015e-01-7.37075809e-01j], [ 6.59761321e-02+5.04689129e-02j, -7.32604698e-02+4.87474207e-01j], [-1.29974542e-01+5.05823754e-02j, 2.19445163e-01-7.25588434e-03j], [ 6.07912591e-03-1.06168382e-02j, 7.44092205e-02-3.65366031e-02j], [-7.02426432e-03+2.64877595e-02j, 2.66509766e-02-3.28387735e-02j]], [[-4.75717574e-01-1.11395206e-01j, -1.68683962e-01-7.53871940e-01j], [-2.42753007e-01+5.11064995e-02j, 1.63664193e+00+3.22370903e-01j], [-1.44085126e+00+4.10080698e-01j, 2.60323983e+00+1.15186298e-01j], [-2.00804580e+00-2.29411043e+00j, -5.82712230e-01-1.19349335e+00j], [ 2.56559424e-01+9.18148802e-02j, -3.33165269e-01-3.67997334e-01j], [-5.73009101e-01+3.60354424e-02j, -2.36867194e-01-8.73383422e-02j]], [[ 9.27606076e-03+8.26768995e-02j, -3.95796317e-01+2.91734572e-01j], [-7.60307326e-01-1.31920765e+00j, 5.49976892e-01+2.53093330e+00j], [-1.14990652e+01+1.79823568e+00j, 2.57693682e+01+9.36167882e+00j], [-1.41226106e+01-1.97062475e+01j, -2.10717368e+00-1.42193496e+01j], [-5.11824082e+01-3.64260906e+01j, 5.34369409e+01+1.34264542e+02j], [ 1.29262446e+02+9.31729066e-01j, 5.60193666e+01+4.41373187e+01j]], [[ 1.14220310e-01-8.02577036e-02j, -2.13299301e-01+6.50567974e-02j], [ 1.28874525e+00+2.74786880e+00j, 2.48314323e-01+1.32533592e+00j], [ 2.36106638e+01-2.27623964e+00j, 1.35814612e+01+5.39901351e+00j], [ 2.34207333e+01+4.27727587e+01j, -2.27386587e-01-7.91489775e+00j], [ 1.22370898e+02+4.60577694e+01j, 7.21782178e+01+6.47879465e+01j], [-2.30501642e+02-4.56070152e+01j, 3.50504688e+01+4.85975417e+01j]], [[ 1.22571494e-02+9.24864312e-05j, 6.89362095e-02+4.78545641e-02j], [-9.55544056e-02+2.55298493e-01j, 3.77737888e-01-3.30037266e-01j], [-4.00072790e+01-4.81738677e+01j, 1.48772899e+01+1.41762044e+02j], [ 1.34059200e+01-1.35063819e+02j, 4.50063572e+01-8.49329880e+01j], [-8.94315772e+02-2.55253910e+03j, -2.22590051e+03+5.16570243e+03j], [ 4.30550031e+03+3.05158106e+03j, 9.52637371e+02+2.99246844e+03j]], [[ 1.35762574e-02+2.34007151e-02j, -4.48595871e-03-4.30039590e-02j], [-5.36158191e-01+1.79755253e-01j, -2.56433793e-01-2.42416698e-02j], [-4.39938045e-01-1.28412885e+02j, 4.41014313e+01-5.60005364e+01j], [ 2.20116049e+02-8.31437389e+01j, -5.27934343e+01-2.40738692e+01j], [ 1.72580131e+03-4.89285585e+03j, 2.68915995e+03-1.84679049e+03j], [-6.87129533e+02+8.25340339e+03j, 1.36983034e+03-5.10152231e+02j]]])Tensor(shape=(2, 6, 6, 2), inds=[_efc85cAAABM, _efc85cAAABL, _efc85cAAABD, k48], tags={I48}),
backend=numpy, dtype=complex128, data=...Tensor(shape=(2, 2), inds=[_efc85cAAABM, k49], tags={I49}),
backend=numpy, dtype=complex128, data=array([[-0.27548952+0.8233728j , 0.33164957+0.36901399j], [-0.25875775-0.4233287j , 0.86776957-0.02851685j]])To keep this cheap we can evaluate local observables or sample configurations using the ‘cluster’ method, where beyond a local region the simple update style gauges are used to approximate the environment (exact on trees).
peps.local_expectation_cluster(
G=qu.pauli("Z") & qu.pauli("Z"),
where=tuple(sitemap[site] for site in edges[10]),
gauges=gauges,
# increase this to include loops and check for convergence
max_distance=0,
)
array(0.34701078-1.83418952e-18j)
We can also sample via the simple update style approximation:
%%time
rng = np.random.default_rng(42)
for _ in range(10):
config, omega = peps.sample_configuration_cluster(
gauges=gauges,
seed=rng,
# single site clusters
max_distance=0,
)
x = "".join(map(str, (config[i] for i in range(N))))
print(x)
00000000000000000000000100110011001100000011000101
00001100000000000000100000000000000000000000000000
00000000000000110000000000010100100011010100000000
01011001000010001000001000000100000011001100010100
10011100011000000000000000000100001110001100000100
10000000010000000000100001110000000110011000000001
10000101111110000101000001100000000000000000000000
00000010000000011000100000000101000000000001000011
11011000100110011001001000010000000010001101000110
10000011000001110000000000000000000010000000000000
CPU times: user 23.3 s, sys: 32.4 ms, total: 23.4 s
Wall time: 23.4 s
If one contracts the amplitude associated with config
, then its ratio to the ‘local’ probability omega
allows one to do unbiased sampling.
# clean up the contraction cache for docs
!rm -rf ctg_cache/