Sequential Quantum Gate Decomposer  v1.9.3
Powerful decomposition of general unitarias into one- and two-qubit gates gates
test_parametric_circuit.py
Go to the documentation of this file.
1 '''
2 Copyright 2020 Peter Rakyta, Ph.D.
3 
4 Licensed under the Apache License, Version 2.0 (the "License");
5 you may not use this file except in compliance with the License.
6 You may obtain a copy of the License at
7 
8  http://www.apache.org/licenses/LICENSE-2.0
9 
10 Unless required by applicable law or agreed to in writing, software
11 distributed under the License is distributed on an "AS IS" BASIS,
12 WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
13 See the License for the specific language governing permissions and
14 limitations under the License.
15 
16 You should have received a copy of the GNU General Public License
17 along with this program. If not, see http://www.gnu.org/licenses/.
18 
19 '''
20 
21 import numpy as np
22 import time
23 import random
24 
25 
26 
27 
28 
29 
30 from qiskit import QuantumCircuit, transpile
31 import numpy as np
32 
33 
34 
35 
36 from squander import utils
37 from squander import N_Qubit_Decomposition_custom
38 from squander import Qiskit_IO
39 
40 
41 
42 try:
43  from mpi4py import MPI
44  MPI_imported = True
45 except ModuleNotFoundError:
46  MPI_imported = False
47 
48 
49 
52 def pauli_exponent( alpha=0.6217*np.pi ):
53  # creating Qiskit quantum circuit
54  qc_orig = QuantumCircuit(5)
55 
56  qc_orig.h(1)
57  qc_orig.cx(1,2)
58 
59  qc_orig.rx(np.pi/2,0)
60  qc_orig.rx(np.pi/2,1)
61  qc_orig.cx(2,4)
62  qc_orig.cx(0,1)
63 
64  qc_orig.rx(np.pi/2,0)
65  qc_orig.h(2)
66  qc_orig.cx(0,2)
67 
68 
69  qc_orig.rx(np.pi/2,0)
70  qc_orig.h(3)
71  qc_orig.rz(alpha,4)
72  qc_orig.cx(0,3)
73 
74  qc_orig.h(0)
75  qc_orig.rz(-alpha,1)
76  qc_orig.cx(2,4)
77 
78  qc_orig.cx(2,1)
79  qc_orig.rz(-alpha,4)
80  qc_orig.cx(3,1)
81 
82  qc_orig.rz(alpha,1)
83  qc_orig.cx(0,1)
84  qc_orig.cx(3,1)
85  qc_orig.cx(4,1)
86 
87  qc_orig.rz(-alpha,1)
88  qc_orig.cx(2,1)
89 
90  qc_orig.rz(alpha,1)
91  qc_orig.cx(3,1)
92  qc_orig.cx(4,1)
93 
94  qc_orig.rz(alpha,1)
95  qc_orig.cx(2,4)
96  qc_orig.cx(0,1)
97 
98  qc_orig.h(0)
99  qc_orig.cx(3,1)
100  qc_orig.cx(0,3)
101 
102  qc_orig.rx(-np.pi/2,0)
103  qc_orig.h(3)
104  qc_orig.cx(0,2)
105 
106  qc_orig.rx(-np.pi/2,0)
107  qc_orig.h(2)
108  qc_orig.cx(0,1)
109 
110  qc_orig.rx(-np.pi/2,0)
111  qc_orig.rx(-np.pi/2,1)
112  qc_orig.cx(2,4)
113  qc_orig.cx(1,2)
114  qc_orig.h(1)
115  return qc_orig
116 
117 
118 
119 
124 def get_unitary_distance( Umtx1, Umtx2 ):
125 
126  product_matrix = np.dot(Umtx1, Umtx2.conj().T)
127  phase = np.angle(product_matrix[0,0])
128  product_matrix = product_matrix*np.exp(-1j*phase)
129  product_matrix = np.eye(len(Umtx1))*2 - product_matrix - product_matrix.conj().T
130  distance = (np.real(np.trace(product_matrix)))/2
131 
132  return distance
133 
134 
135 
138 def get_optimized_circuit( alpha, optimizer='BFGS', optimized_parameters=None ):
139 
140  filename = 'data/19CNOT.qasm'
141  qc_trial = QuantumCircuit.from_qasm_file( filename )
142  #qc_trial = QuantumCircuit.from_qasm_file( filename )
143  qc_trial = transpile(qc_trial, optimization_level=1, basis_gates=['cz', 'cx', 'u3', 'rz', 'u', 'rx'], layout_method='sabre')
144 
145 
146  qc_orig = pauli_exponent( alpha )
147  #qc_orig = transpile(qc_orig, optimization_level=3, basis_gates=['cx', 'u3'], layout_method='sabre')
148  #print('global phase: ', qc_orig.global_phase)
149 
150  Umtx_orig = utils.get_unitary_from_qiskit_circuit( qc_orig )
151  #Umtx_trial = utils.get_unitary_from_qiskit_circuit( qc_trial )
152 
153  #print( np.abs(Umtx_orig @ Umtx_trial.conj().T) )
154 
155 
156  iteration_max = 1
157 
158  for jdx in range(iteration_max):
159 
160  cDecompose = N_Qubit_Decomposition_custom( Umtx_orig.conj().T )
161 
162  # setting the tolerance of the optimization process. The final error of the decomposition would scale with the square root of this value.
163  cDecompose.set_Optimization_Tolerance( 1e-5 )
164 
165  # importing the quantum circuit
166  cDecompose.import_Qiskit_Circuit( qc_trial )
167 
168  # set the number of successive identical blocks in the optimalization of disentanglement of the n-th qubits
169  cDecompose.set_Optimization_Blocks( 200 )
170 
171  # set the optimizer
172  cDecompose.set_Optimizer( optimizer )
173 
174  # set the initial parameters if given
175  if not ( optimized_parameters is None ):
176  cDecompose.set_Optimized_Parameters( optimized_parameters )
177 
178  # set verbosity
179  cDecompose.set_Verbose( 4 )
180 
181  # set the cost function to Hilbert-Schmidt test
182  cDecompose.set_Cost_Function_Variant( 3 )
183 
184  # starting the decomposition
185  cDecompose.Start_Decomposition()
186 
187  # getting the new optimized parameters
188  optimized_parameters_loc = cDecompose.get_Optimized_Parameters()
189 
190  qc_final = cDecompose.get_Qiskit_Circuit()
191 
192  # get the unitary of the final circuit
193  Umtx_recheck = utils.get_unitary_from_qiskit_circuit( qc_final )
194 
195  # get the decomposition error
196  decomposition_error = get_unitary_distance(Umtx_orig, Umtx_recheck)
197  print('recheck decomposition error: ', decomposition_error)
198 
199 
200  if decomposition_error < 1e-3:
201  break
202 
203 
204  if decomposition_error < 1e-3:
205  return qc_final, optimized_parameters_loc
206  else:
207  assert(decomposition_error < 1e-3)
208  return None, None
209 
210 
211 
212 
214  """This is a test class of the python iterface to the decompsition classes of the QGD package"""
215 
216 
218 
219 
220  qc_trial = pauli_exponent()
221 
222 
223 
224  # get the unitary of the final circuit
225  Umtx_check = utils.get_unitary_from_qiskit_circuit( qc_trial )
226 
227  qc_trial = pauli_exponent()
228  Circuit_Squander, parameters = Qiskit_IO.convert_Qiskit_to_Squander( qc_trial )
229 
230  cDecompose = N_Qubit_Decomposition_custom( Umtx_check.conj().T )
231 
232  cDecompose.set_Gate_Structure(Circuit_Squander)
233 
234  cDecompose.set_Cost_Function_Variant( 4 )
235 
236  # starting the decomposition
237  cost_function = cDecompose.Optimization_Problem( parameters )
238 
239  assert(cost_function < 1e-3)
240 
241 
242 
243  def test_optimizer(self):
244 
245 
246  # determine random parameter value alpha
247  alpha = 1.823631161607293
248 
249  # determine the quantum circuit at parameter value alpha with BFGS2 optimizer
250  qc, optimized_parameters = get_optimized_circuit( alpha, optimizer='BFGS2' )
251 
252  # determine the quantum circuit at parameter value alpha with BFGS optimizer
253  qc, optimized_parameters_tmp = get_optimized_circuit( alpha+0.005, optimizer='BFGS', optimized_parameters=optimized_parameters )
254 
255  # determine the quantum circuit at parameter value alpha with ADAM optimizer
256  qc, optimized_parameters_tmp = get_optimized_circuit( alpha+0.005, optimizer='GRAD_DESCEND', optimized_parameters=optimized_parameters )
257 
258  # determine the quantum circuit at parameter value alpha with ADAM optimizer
259  qc, optimized_parameters_tmp = get_optimized_circuit( alpha+0.005, optimizer='ADAM', optimized_parameters=optimized_parameters )
260 
261 
262 
263 
264 
def get_unitary_distance(Umtx1, Umtx2)
Calcuates the distance between two unitaries according to Eq.
def pauli_exponent(alpha=0.6217 *np.pi)
???????????
def get_optimized_circuit(alpha, optimizer='BFGS', optimized_parameters=None)
???????????