Sequential Quantum Gate Decomposer  v1.9.3
Powerful decomposition of general unitarias into one- and two-qubit gates gates
test_VQE.py
Go to the documentation of this file.
1 from scipy.stats import unitary_group
2 import numpy as np
3 from squander.variational_quantum_eigensolver.qgd_Variational_Quantum_Eigensolver_Base import qgd_Variational_Quantum_Eigensolver_Base as Variational_Quantum_Eigensolver_Base
4 from squander import utils as utils
5 import time
6 import sys
7 from squander.gates.qgd_Circuit_Wrapper import qgd_Circuit_Wrapper as Circuit
8 import scipy as sp
9 import pytest
10 
11 
12 sigmax = sp.sparse.csr_matrix(np.array([[0,1],
13  [1,0]])+0.j)
14 sigmay = sp.sparse.csr_matrix(np.array([[0,0+-1j],
15  [0+1j,0]])+0.j)
16 sigmaz = sp.sparse.csr_matrix(np.array([[1,0],
17  [0,-1]])+0.j)
18 
19 
20 def generate_hamiltonian(topology, n):
21 
22  Hamiltonian = sp.sparse.coo_matrix((2**n, 2**n), dtype=np.complex128)
23  for i in topology:
24  if i[0]==0:
25  lhs_1 = sp.sparse.eye(1,format='coo')
26  else:
27  lhs_1 = sp.sparse.eye(2,format='coo')
28  for k in range(i[0]-1):
29  lhs_1 = sp.sparse.kron(lhs_1,sp.sparse.eye(2,format='coo'),format='coo')
30  if i[0]==n-1:
31  rhs_1 = sp.sparse.eye(1,format='coo')
32  else:
33  rhs_1 = sp.sparse.eye(2,format='coo')
34  for k in range(n-i[0]-2):
35  rhs_1 = sp.sparse.kron(rhs_1,sp.sparse.eye(2,format='coo'),format='coo')
36  if i[1]==0:
37  lhs_2 = sp.sparse.eye(1,format='coo')
38  else:
39  lhs_2 = sp.sparse.eye(2,format='coo')
40  for k in range(i[1]-1):
41  lhs_2 = sp.sparse.kron(lhs_2,sp.sparse.eye(2,format='coo'),format='coo')
42  if i[1]==n-1:
43  rhs_2 = sp.sparse.eye(1,format='coo')
44  else:
45  rhs_2 = sp.sparse.eye(2,format='coo')
46  for k in range(n-i[1]-2):
47  rhs_2 = sp.sparse.kron(rhs_2,sp.sparse.eye(2,format='coo'),format='coo')
48  Hamiltonian += -0.5*sp.sparse.kron(sp.sparse.kron(lhs_1,sigmax,format='coo'),rhs_1,format='coo')@sp.sparse.kron(sp.sparse.kron(lhs_2 ,sigmax,format='coo'),rhs_2 ,format='coo')
49  Hamiltonian += -0.5*sp.sparse.kron(sp.sparse.kron(lhs_1,sigmay,format='coo'),rhs_1,format='coo')@sp.sparse.kron(sp.sparse.kron(lhs_2 ,sigmay,format='coo'),rhs_2 ,format='coo')
50  Hamiltonian += -0.5*sp.sparse.kron(sp.sparse.kron(lhs_1,sigmaz,format='coo'),rhs_1,format='coo')@sp.sparse.kron(sp.sparse.kron(lhs_2 ,sigmaz,format='coo'),rhs_2 ,format='coo')
51  for i in range(n):
52 
53  if i==0:
54  lhs_1 = sp.sparse.eye(1,format='coo')
55  else:
56  lhs_1 = sp.sparse.eye(2,format='coo')
57  for k in range(i-1):
58  lhs_1 = sp.sparse.kron(lhs_1,sp.sparse.eye(2,format='coo'),format='coo')
59  if i==n-1:
60  rhs_1 = sp.sparse.eye(1,format='coo')
61  else:
62  rhs_1 = sp.sparse.eye(2,format='coo')
63  for k in range(n-i-2):
64  rhs_1 = sp.sparse.kron(rhs_1,sp.sparse.eye(2,format='coo'),format='coo')
65  Hamiltonian += -0.5*sp.sparse.kron(sp.sparse.kron(lhs_1,sigmaz,format='coo'),rhs_1,format='coo')
66 
67 
68  return Hamiltonian.tocsr()
69 
70 
71 
72 
73 
74 
75 class Test_VQE:
76 
77  def test_VQE_Identity(self):
78  layers = 1
79  blocks = 1
80 
81  qbit_num = 7
82  Hamiltonian = sp.sparse.eye(2**qbit_num,format="csr")
83 
84  config = {"agent_lifetime":500,
85  "optimization_tolerance": -7.1,
86  "max_inner_iterations":10,
87  "max_iterations":50,
88  "learning_rate": 2e-1,
89  "agent_num":64,
90  "agent_exploration_rate":0.3,
91  "max_inner_iterations_adam":50000}
92 
93  # initiate the VQE object with the Hamiltonian
94  VQE_eye = Variational_Quantum_Eigensolver_Base(Hamiltonian, qbit_num, config)
95 
96  # set the optimization engine
97  VQE_eye.set_Optimizer("GRAD_DESCEND")
98 
99  # set the ansatz variant (U3 rotations and CNOT gates)
100  VQE_eye.set_Ansatz("HEA")
101 
102  # generate the circuit ansatz for the optimization
103  VQE_eye.Generate_Circuit(layers, blocks)
104 
105  # create initial parameters
106  param_num = VQE_eye.get_Parameter_Num()
107  parameters = np.random.random( (param_num,) )
108 
109  VQE_eye.set_Optimized_Parameters(parameters)
110 
111  # start an etap of the optimization (max_inner_iterations iterations)
112  VQE_eye.Start_Optimization()
113 
114  # retrieve QISKIT format of the optimized circuit
115  quantum_circuit = VQE_eye.get_Qiskit_Circuit()
116 
117  # retrieve the optimized parameter
118  parameters = VQE_eye.get_Optimized_Parameters()
119 
120  # evaluate the VQE energy at the optimized parameters
121  Energy = VQE_eye.Optimization_Problem( parameters )
122 
123 
124  assert (abs(Energy-1)<1e-4)
125 
127 
128  layers = 1
129  blocks = 1
130 
131  qbit_num = 9
132  topology = []
133  for idx in range(qbit_num-1):
134  topology.append( (idx, idx+1) )
135 
136  Hamiltonian = generate_hamiltonian(topology, qbit_num)
137 
138  config = {"agent_lifetime":10,
139  "max_inner_iterations":1000,
140  "max_iterations":1000,
141  "agent_num":3,
142  "agent_exploration_rate":0.5,
143  "max_inner_iterations_adam":500,
144  "convergence_length": 300}
145 
146  # initiate the VQE object with the Hamiltonian
147  VQE_Heisenberg = Variational_Quantum_Eigensolver_Base(Hamiltonian, qbit_num, config)
148 
149  # set the optimization engine
150  VQE_Heisenberg.set_Optimizer("AGENTS")
151 
152  # set the ansatz variant (U3 rotations and CNOT gates)
153  VQE_Heisenberg.set_Ansatz("HEA_ZYZ")
154 
155  # generate the circuit ansatz for the optimization
156  VQE_Heisenberg.Generate_Circuit(layers, blocks)
157 
158  # create initial parameters
159  param_num = VQE_Heisenberg.get_Parameter_Num()
160  parameters = np.random.random( (param_num,) )
161 
162  VQE_Heisenberg.set_Optimized_Parameters(parameters)
163 
164  # start an etap of the optimization (max_inner_iterations iterations)
165  VQE_Heisenberg.Start_Optimization()
166 
167  # retrieve QISKIT format of the optimized circuit
168  quantum_circuit = VQE_Heisenberg.get_Qiskit_Circuit()
169 
170  # print the quantum circuit
171  lambdas, vecs = sp.sparse.linalg.eigs(Hamiltonian)
172 
173  # retrieve the optimized parameter
174  parameters = VQE_Heisenberg.get_Optimized_Parameters()
175 
176  # evaluate the VQE energy at the optimized parameters
177  Energy = VQE_Heisenberg.Optimization_Problem( parameters )
178 
179  print('Expected energy: ', np.real(lambdas[0]))
180  print('Obtained energy: ', Energy)
181  assert ((Energy - np.real(lambdas[0]))<1e-2)
182 
183 
184 
185 
186 
187 
def generate_hamiltonian(topology, n)
Definition: test_VQE.py:20
A base class to solve VQE problems This class can be used to approximate the ground state of the inpu...