Quantum Gate Decomposer  v1.3
Powerful decomposition of almost any unitary into U3 and CNOT gates
Sub_Matrix_Decomposition.cpp
Go to the documentation of this file.
1 /*
2 Created on Fri Jun 26 14:13:26 2020
3 Copyright (C) 2020 Peter Rakyta, Ph.D.
4 
5 This program is free software: you can redistribute it and/or modify
6 it under the terms of the GNU General Public License as published by
7 the Free Software Foundation, either version 3 of the License, or
8 (at your option) any later version.
9 
10 This program is distributed in the hope that it will be useful,
11 but WITHOUT ANY WARRANTY; without even the implied warranty of
12 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 GNU General Public License for more details.
14 
15 You should have received a copy of the GNU General Public License
16 along with this program. If not, see http://www.gnu.org/licenses/.
17 
18 @author: Peter Rakyta, Ph.D.
19 */
25 
26 
27 
28 
29 
38 Sub_Matrix_Decomposition::Sub_Matrix_Decomposition( Matrix Umtx_in, int qbit_num_in, bool optimize_layer_num_in=false, guess_type initial_guess_in= CLOSE_TO_ZERO ) : Decomposition_Base(Umtx_in, qbit_num_in, initial_guess_in) {
39 
40  // logical value. Set true if finding the minimum number of operation layers is required (default), or false when the maximal number of CNOT gates is used (ideal for general unitaries).
41  optimize_layer_num = optimize_layer_num_in;
42 
43  // A string describing the type of the class
45 
46  // The global minimum of the optimization problem
48 
49  // number of iteratrion loops in the optimization
50  iteration_loops[2] = 3;
51 
52  // logical value indicating whether the quasi-unitarization of the submatrices was done or not
53  subdisentaglement_done = false;
54 
55  // filling in numbers that were not given in the input
56  for ( std::map<int,int>::iterator it = max_layer_num_def.begin(); it!=max_layer_num_def.end(); it++) {
57  if ( max_layer_num.count( it->first ) == 0 ) {
58  max_layer_num.insert( std::pair<int, int>(it->first, it->second) );
59  }
60  }
61 
62 
63 
64 
65 
66 
67 }
68 
69 
74 
75 }
76 
77 
78 
79 
80 
81 
86 
88  if (verbose) {
89  printf("Sub-disentaglement already done.\n");
90  }
91  return;
92  }
93 
94  if (verbose) {
95  printf("\nDisentagling submatrices.\n");
96  }
97 
98  // setting the global target minimum
101 
102  // check if it needed to do the subunitarization
104  if (verbose) {
105  printf("Disentanglig not needed\n");
106  }
108  subdisentaglement_done = true;
109  return;
110  }
111 
112 
113 
114  if ( !check_optimization_solution() ) {
115  // Adding the operations of the successive layers
116 
117  //measure the time for the decompositin
118  clock_t start_time = time(NULL);
119 
120  // the maximal number of layers in the subdeconposition
121  int max_layer_num_loc;
122  try {
123  max_layer_num_loc = max_layer_num[qbit_num];
124  }
125  catch (...) {
126  throw "Layer number not given";
127  }
128 
129 
130  // the number of succeeding identicallayers in the subdeconposition
131  int identical_blocks_loc;
132  try {
133  identical_blocks_loc = identical_blocks[qbit_num];
134  if (identical_blocks_loc==0) {
135  identical_blocks_loc = 1;
136  }
137  }
138  catch (...) {
139  identical_blocks_loc=1;
140  }
141 
142  while ( layer_num < max_layer_num_loc ) {
143 
144  int control_qbit_loc = qbit_num-1;
145  int solution_guess_num = parameter_num;
146 
147  for (int target_qbit_loc = 0; target_qbit_loc<control_qbit_loc; target_qbit_loc++ ) {
148 
149  for (int idx=0; idx<identical_blocks_loc; idx++) {
150 
151  // creating block of operations
152  Operation_block* block = new Operation_block( qbit_num );
153 
154  // add CNOT gate to the block
155  block->add_cnot_to_end(control_qbit_loc, target_qbit_loc);
156 
157  // adding U3 operation to the block
158  bool Theta = true;
159  bool Phi = false;
160  bool Lambda = true;
161  block->add_u3_to_end(target_qbit_loc, Theta, Phi, Lambda);
162  block->add_u3_to_end(control_qbit_loc, Theta, Phi, Lambda);
163 
164  // adding the opeartion block to the operations
165  add_operation_to_end( block );
166 
167  }
168  }
169 
170  // get the number of blocks
171  layer_num = operations.size();
172 
173  // Do the optimization
174  if (optimize_layer_num || layer_num >= max_layer_num_loc ) {
175 
176  // solve the optzimalization problem to find the correct mninimum
177  if ( optimized_parameters == NULL ) {
179  }
180  else {
181  solve_optimization_problem( optimized_parameters, solution_guess_num);
182  }
183 
185  break;
186  }
187  }
188 
189  }
190 
191  if (verbose) {
192  printf("--- %f seconds elapsed during the decomposition ---\n\n", float(time(NULL) - start_time));
193  }
194 
195 
196  }
197 
198 
199 
201  if (verbose) {
202  printf("Sub-disentaglement was succesfull.\n\n");
203  }
204  }
205  else {
206  if (verbose) {
207  printf("Sub-disentaglement did not reach the tolerance limit.\n\n");
208  }
209  }
210 
211 
212  // indicate that the unitarization of the sumbatrices was done
213  subdisentaglement_done = true;
214 
215  // The subunitarized matrix
217 }
218 
219 
220 
221 
222 
228 void Sub_Matrix_Decomposition::solve_layer_optimization_problem( int num_of_parameters, gsl_vector *solution_guess_gsl) {
229 
230 
231  if (operations.size() == 0 ) {
232  return;
233  }
234 
235 
236  if (solution_guess_gsl == NULL) {
237  solution_guess_gsl = gsl_vector_alloc(num_of_parameters);
238  }
239 
240  if (optimized_parameters == NULL) {
241  optimized_parameters = (double*)qgd_calloc(num_of_parameters,sizeof(double), 64);
242  memcpy(optimized_parameters, solution_guess_gsl->data, num_of_parameters*sizeof(double) );
243  }
244 
245  // maximal number of iteration loops
246  int iteration_loops_max;
247  try {
248  iteration_loops_max = std::max(iteration_loops[qbit_num], 1);
249  }
250  catch (...) {
251  iteration_loops_max = 1;
252  }
253 
254 
255  // do the optimization loops
256  for (int idx=0; idx<iteration_loops_max; idx++) {
257 
258  size_t iter = 0;
259  int status;
260 
261  const gsl_multimin_fdfminimizer_type *T;
262  gsl_multimin_fdfminimizer *s;
263 
264  Sub_Matrix_Decomposition* par = this;
265 
266 
267  gsl_multimin_function_fdf my_func;
268 
269 
270  my_func.n = num_of_parameters;
271  my_func.f = optimization_problem;
272  my_func.df = optimization_problem_grad;
273  my_func.fdf = optimization_problem_combined;
274  my_func.params = par;
275 
276 
277  T = gsl_multimin_fdfminimizer_vector_bfgs2;
278  s = gsl_multimin_fdfminimizer_alloc (T, num_of_parameters);
279 
280 
281  gsl_multimin_fdfminimizer_set (s, &my_func, solution_guess_gsl, 0.1, 0.1);
282 
283  do {
284  iter++;
285 
286  status = gsl_multimin_fdfminimizer_iterate (s);
287 
288  if (status) {
289  break;
290  }
291 
292  status = gsl_multimin_test_gradient (s->gradient, 1e-1);
293  /*if (status == GSL_SUCCESS) {
294  printf ("Minimum found\n");
295  }*/
296 
297  } while (status == GSL_CONTINUE && iter < 100);
298 
299  if (current_minimum > s->f) {
300  current_minimum = s->f;
301  memcpy( optimized_parameters, s->x->data, num_of_parameters*sizeof(double) );
302  gsl_multimin_fdfminimizer_free (s);
303 
304  for ( int jdx=0; jdx<num_of_parameters; jdx++) {
305  solution_guess_gsl->data[jdx] = solution_guess_gsl->data[jdx] + (2*double(rand())/double(RAND_MAX)-1)*2*M_PI/100;
306  }
307  }
308  else {
309  for ( int jdx=0; jdx<num_of_parameters; jdx++) {
310  solution_guess_gsl->data[jdx] = solution_guess_gsl->data[jdx] + (2*double(rand())/double(RAND_MAX)-1)*2*M_PI;
311  }
312 
313  gsl_multimin_fdfminimizer_free (s);
314  }
315 
316 
317 
318  }
319 
320 
321 
322 }
323 
324 
325 
326 
332 double Sub_Matrix_Decomposition::optimization_problem( const double* parameters ) {
333 
334  // get the transformed matrix with the operations in the list
335  Matrix matrix_new = get_transformed_matrix( parameters, operations.begin(), operations.size(), Umtx );
336 
337  double cost_function = get_submatrix_cost_function(matrix_new); //NEW METHOD
338 
339 
340  return cost_function;
341 }
342 
343 
350 double Sub_Matrix_Decomposition::optimization_problem( const gsl_vector* parameters, void* void_instance ) {
351 
352  Sub_Matrix_Decomposition* instance = reinterpret_cast<Sub_Matrix_Decomposition*>(void_instance);
353  std::vector<Operation*> operations_loc = instance->get_operations();
354 
355  Matrix Umtx_loc = instance->get_Umtx();
356  Matrix matrix_new = instance->get_transformed_matrix( parameters->data, operations_loc.begin(), operations_loc.size(), Umtx_loc );
357 
358 
359  double cost_function = get_submatrix_cost_function(matrix_new); //NEW METHOD
360 
361 
362  return cost_function;
363 }
364 
365 
366 
367 
368 
369 
376 int Sub_Matrix_Decomposition::set_identical_blocks( int qbit, int identical_blocks_in ) {
377 
378  std::map<int,int>::iterator key_it = identical_blocks.find( qbit );
379 
380  if ( key_it != identical_blocks.end() ) {
381  identical_blocks.erase( key_it );
382  }
383 
384  identical_blocks.insert( std::pair<int, int>(qbit, identical_blocks_in) );
385 
386  return 0;
387 
388 }
389 
390 
396 int Sub_Matrix_Decomposition::set_identical_blocks( std::map<int, int> identical_blocks_in ) {
397 
398  for ( std::map<int,int>::iterator it=identical_blocks_in.begin(); it!= identical_blocks_in.end(); it++ ) {
399  set_identical_blocks( it->first, it->second );
400  }
401 
402  return 0;
403 
404 }
405 
406 
412 
414 
415  // setting computational parameters
421 
422  if ( extract_operations(static_cast<Operation_block*>(ret)) != 0 ) {
423  printf("Sub_Matrix_Decomposition::clone(): extracting operations was not succesfull\n");
424  exit(-1);
425  }
426 
427  return ret;
428 
429 }
430 
431 
432 
bool optimize_layer_num
logical value. Set true to optimize the minimum number of operation layers required in the decomposit...
std::vector< Operation * > get_operations()
Call to get the operations stored in the class.
Matrix get_transformed_matrix(const double *parameters, std::vector< Operation * >::iterator operations_it, int num_of_operations)
Calculate the transformed matrix resulting by an array of operations on the matrix Umtx.
int extract_operations(Operation_block *op_block)
Call to extract the operations stored in the class.
void solve_layer_optimization_problem(int num_of_parameters, gsl_vector *solution_guess_gsl)
Call to solve layer by layer the optimization problem.
Sub_Matrix_Decomposition(Matrix Umtx_in, int qbit_num_in, bool optimize_layer_num_in, guess_type initial_guess_in)
Constructor of the class.
void * qgd_calloc(size_t element_num, size_t size, size_t alignment)
custom defined memory allocation function.
Definition: common.cpp:38
Matrix get_Umtx()
Call to retrive a pointer to the unitary to be transformed.
void disentangle_submatrices()
Start the optimization process to disentangle the most significant qubit from the others.
double current_minimum
The current minimum of the optimization problem.
int set_max_layer_num(int n, int max_layer_num_in)
Set the maximal number of layers used in the subdecomposition of the n-th qubit.
std::map< int, int > max_layer_num
A map of <int n: int num> indicating that how many layers should be used in the subdecomposition proc...
static void optimization_problem_combined(const gsl_vector *parameters, void *void_instance, double *f0, gsl_vector *grad)
Call to calculate both the cost function and the its gradient components.
guess_type initial_guess
type to guess the initial values for the optimization. Possible values: ZEROS=0, RANDOM=1,...
std::map< int, int > identical_blocks
A map of <int n: int num> indicating that how many identical succesive blocks should be used in the d...
int set_identical_blocks(int qbit, int identical_blocks_in)
Set the number of identical successive blocks during the subdecomposition of the qbit-th qubit.
operation_type type
The type of the operation (see enumeration operation_type)
Definition: Operation.h:48
double get_submatrix_cost_function(Matrix &matrix)
Call to calculate the cost function of a given matrix during the submatrix decomposition process.
bool subdisentaglement_done
logical value indicating whether the disentamglement of a qubit from the othetrs was done or not
int qbit_num
number of qubits spanning the matrix of the operation
Definition: Operation.h:46
A class responsible for grouping CNOT and U3 operations into layers.
Operation_block()
Deafult constructor of the class.
Sub_Matrix_Decomposition * clone()
Call to create a clone of the present class.
void set_max_iteration(int max_iterations_in)
Call to set the maximal number of the iterations in the optimization process.
std::map< int, int > iteration_loops
A map of <int n: int num> indicating the number of iteration in each step of the decomposition.
double optimization_problem(const double *parameters)
The optimization problem of the final optimization.
~Sub_Matrix_Decomposition()
Destructor of the class.
unsigned int parameter_num
the number of free parameters of the operation
Definition: Operation.h:56
A class containing basic methods for the decomposition process.
int layer_num
number of operation layers
Class to store data of complex arrays and its properties.
Definition: matrix.h:12
std::vector< Operation * > operations
The list of stored operations.
bool verbose
Logical variable. Set true for verbose mode, or to false to suppress output messages.
Matrix subdecomposed_mtx
The subdecomposed matrix.
guess_type
Type definition of the types of the initial guess.
void add_u3_to_end(int target_qbit, bool Theta, bool Phi, bool Lambda)
Append a U3 gate to the list of operations.
double * optimized_parameters
The optimized parameters for the operations.
void add_operation_to_end(Operation *operation)
Append a general operation to the list of operations.
Header file for a class responsible for the disentanglement of one qubit from the others.
static void optimization_problem_grad(const gsl_vector *parameters, void *void_instance, gsl_vector *grad)
Calculate the approximate derivative (f-f0)/(x-x0) of the cost function with respect to the free para...
void solve_optimization_problem(double *solution_guess, int solution_guess_num)
This method can be used to solve the main optimization problem which is devidid into sub-layer optimi...
Matrix Umtx
The unitary to be decomposed.
static std::map< int, int > max_layer_num_def
A map of <int n: int num> indicating that how many layers should be used in the subdecomposition proc...
bool check_optimization_solution()
check_optimization_solution
A class responsible for the disentanglement of one qubit from the others.
int optimization_block
number of operation blocks used in one shot of the optimization process
void add_cnot_to_end(int control_qbit, int target_qbit)
Append a C_NOT gate operation to the list of operations.
int set_iteration_loops(int n, int iteration_loops_in)
Set the number of iteration loops during the subdecomposition of the n-th qubit.
double global_target_minimum
The global target minimum of the optimization problem.
void set_optimization_blocks(int optimization_block_in)
Call to set the number of operation blocks to be optimized in one shot.
int max_iterations
Maximal number of iterations allowed in the optimization process.