Sequential Quantum Gate Decomposer  v1.9.3
Powerful decomposition of general unitarias into one- and two-qubit gates 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 2020 Peter Rakyta, Ph.D.
4 
5 Licensed under the Apache License, Version 2.0 (the "License");
6 you may not use this file except in compliance with the License.
7 You may obtain a copy of the License at
8 
9  http://www.apache.org/licenses/LICENSE-2.0
10 
11 Unless required by applicable law or agreed to in writing, software
12 distributed under the License is distributed on an "AS IS" BASIS,
13 WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14 See the License for the specific language governing permissions and
15 limitations under the License.
16 
17 @author: Peter Rakyta, Ph.D.
18 */
25 
26 #include "BFGS_Powell.h"
27 
28 
29 #ifdef __MPI__
30 #include <mpi.h>
31 #endif // MPI
32 
33 //tbb::spin_mutex my_mutex;
34 
40 
41  // logical value. Set true if finding the minimum number of gate layers is required (default), or false when the maximal number of CNOT gates is used (ideal for general unitaries).
42  optimize_layer_num = false;
43 
44  // A string describing the type of the class
46 
47  // The global minimum of the optimization problem
49 
50  // logical value indicating whether the quasi-unitarization of the submatrices was done or not
51  subdisentaglement_done = false;
52 
53  // custom gate structure used in the decomposition
54  unit_gate_structure = NULL;
55 
57 
58 }
59 
60 
69 Sub_Matrix_Decomposition::Sub_Matrix_Decomposition( Matrix Umtx_in, int qbit_num_in, bool optimize_layer_num_in, std::map<std::string, Config_Element>& config_in, guess_type initial_guess_in= CLOSE_TO_ZERO ) : Decomposition_Base(Umtx_in, qbit_num_in, config_in, initial_guess_in) {
70 
71  // logical value. Set true if finding the minimum number of gate layers is required (default), or false when the maximal number of CNOT gates is used (ideal for general unitaries).
72  optimize_layer_num = optimize_layer_num_in;
73 
74  // A string describing the type of the class
76 
77  // The global minimum of the optimization problem
79 
80  // number of iteratrion loops in the optimization
81  iteration_loops[2] = 3;
82 
83  // logical value indicating whether the quasi-unitarization of the submatrices was done or not
84  subdisentaglement_done = false;
85 
86  // filling in numbers that were not given in the input
87  for ( std::map<int,int>::iterator it = max_layer_num_def.begin(); it!=max_layer_num_def.end(); it++) {
88  if ( max_layer_num.count( it->first ) == 0 ) {
89  max_layer_num.insert( std::pair<int, int>(it->first, it->second) );
90  }
91  }
92 
93  // custom gate structure used in the decomposition
94  unit_gate_structure = NULL;
95 
97 
98 
99 }
100 
101 
106 
107  if ( unit_gate_structure != NULL ) {
108  delete unit_gate_structure;
109  }
110  unit_gate_structure = NULL;
111 
112 }
113 
114 
115 
116 
117 
118 
123 
124  if (subdisentaglement_done) {
125  std::stringstream sstream;
126  sstream << "Sub-disentaglement already done." << std::endl;
127  print(sstream, 2);
128  return;
129  }
130 
131  std::stringstream sstream;
132  sstream << std::endl << "Disentagling submatrices." << std::endl;
133  print(sstream, 2);
134 
135 
136  // setting the global target minimum
139 
140  // check if it needed to do the subunitarization
142  std::stringstream sstream;
143  sstream << "Disentanglig not needed" << std::endl;
144  print(sstream, 2);
146  subdisentaglement_done = true;
147  return;
148  }
149 
150 
151  if ( !check_optimization_solution() ) {
152  // Adding the gates of the successive layers
153 
154  //measure the time for the decompositin
155  tbb::tick_count start_time = tbb::tick_count::now();
156 
157  // the maximal number of layers in the subdeconposition
158  int max_layer_num_loc;
159  try {
160  max_layer_num_loc = max_layer_num[qbit_num];
161  }
162  catch (...) {
163  std::stringstream sstream;
164  sstream << "Layer number not given" << std::endl;
165  print(sstream, 0);
166  throw sstream.str();
167  }
168 
169  while ( layer_num < max_layer_num_loc ) {
170 
171  // add another gate layers to the gate structure used in the decomposition
172  add_gate_layers();
173 
174  // get the number of blocks
175  layer_num = gates.size();
176 
177  // Do the optimization
178  if (optimize_layer_num || layer_num >= max_layer_num_loc ) {
179 
180  // solve the optimization problem to find the correct minimum
181  solve_optimization_problem( NULL, 0);
182 
184  break;
185  }
186  }
187 
188  }
189 
190 
191  tbb::tick_count current_time = tbb::tick_count::now();
192 
193  std::stringstream sstream;
194  sstream << "--- " << (current_time - start_time).seconds() << " seconds elapsed during the decomposition ---" << std::endl << std::endl;
195  print(sstream, 2);
196 
197  }
198 
199 
200 
202  std::stringstream sstream;
203  sstream << "Sub-disentaglement was succesfull." << std::endl << std::endl;
204  print(sstream, 1);
205  }
206  else {
207  std::stringstream sstream;
208  sstream << "Sub-disentaglement did not reach the tolerance limit." << std::endl << std::endl;
209  print(sstream, 1);
210  }
211 
212 
213  // indicate that the unitarization of the sumbatrices was done
214  subdisentaglement_done = true;
215 
216  // The subunitarized matrix
219 
220 }
221 
226 
227  if ( unit_gate_structure == NULL ) {
228  // add default gate structure if custom is not given
230  return;
231  }
232 
234  // add custom gate structure
235 
236  // the number of succeeding identical layers in the subdecomposition
237  int identical_blocks_loc;
238  try {
239  identical_blocks_loc = identical_blocks[qbit_num];
240  if (identical_blocks_loc==0) {
241  identical_blocks_loc = 1;
242  }
243  }
244  catch (...) {
245  identical_blocks_loc=1;
246  }
247 
248  // get the list of sub-blocks in the custom gate structure
249  std::vector<Gate*> gates = unit_gate_structure->get_gates();
250 
251  for (std::vector<Gate*>::iterator gates_it = gates.begin(); gates_it != gates.end(); gates_it++ ) {
252 
253  // The current gate
254  Gate* gate = *gates_it;
255 
256  for (int idx=0; idx<identical_blocks_loc; idx++) {
257  add_gate( (Gate*)(gate->clone()) );
258  }
259 
260  }
261 
262 
263 }
264 
265 
266 
271 
272  int control_qbit_loc = qbit_num-1;
273 
274  // the number of succeeding identical layers in the subdecomposition
275  int identical_blocks_loc;
276  try {
277  identical_blocks_loc = identical_blocks[qbit_num];
278  if (identical_blocks_loc==0) {
279  identical_blocks_loc = 1;
280  }
281  }
282  catch (...) {
283  identical_blocks_loc=1;
284  }
285 
286 
287  for (int target_qbit_loc = 0; target_qbit_loc<control_qbit_loc; target_qbit_loc++ ) {
288 
289  for (int idx=0; idx<identical_blocks_loc; idx++) {
290 
291  // creating block of gates
292  Gates_block* block = new Gates_block( qbit_num );
293 
294  // adding U3 gate to the block
295  block->add_u3(target_qbit_loc);
296  block->add_u3(control_qbit_loc);
297 
298 
299  // add CNOT gate to the block
300  block->add_cnot(target_qbit_loc, control_qbit_loc);
301 
302  // adding the opeartion block to the gates
303  add_gate( block );
304 
305  }
306  }
307 
308 
309 }
310 
311 
312 
318 
319  unit_gate_structure = block_in->clone();
320 
321 }
322 
323 
324 
325 
326 
333 
334 
335  if (gates.size() == 0 ) {
336  return;
337  }
338 
339 
340  if (solution_guess.size() == 0 ) {
341  solution_guess = Matrix_real(num_of_parameters,1);
342  }
343 
344 
345  if (optimized_parameters_mtx.size() == 0) {
346  optimized_parameters_mtx = Matrix_real(1, num_of_parameters);
347  memcpy(optimized_parameters_mtx.get_data(), solution_guess.get_data(), num_of_parameters*sizeof(double) );
348  }
349 
350 
351  // maximal number of iteration loops
352  int iteration_loops_max;
353  try {
354  iteration_loops_max = std::max(iteration_loops[qbit_num], 1);
355  }
356  catch (...) {
357  iteration_loops_max = 1;
358  }
359 
360  // maximal number of inner iterations overriden by config
361  long long max_inner_iterations_loc;
362  if ( config.count("max_inner_iterations") > 0 ) {
363  config["max_inner_iterations"].get_property( max_inner_iterations_loc );
364  }
365  else {
366  max_inner_iterations_loc =max_inner_iterations;
367  }
368 
369  // random generator of real numbers
370  std::uniform_real_distribution<> distrib_real(0.0, 2*M_PI);
371 
372  // do the optimization loops
373  for (long long idx=0; idx<iteration_loops_max; idx++) {
374 
375  BFGS_Powell cBFGS_Powell(optimization_problem_combined, this);
376  double f = cBFGS_Powell.Start_Optimization(solution_guess, max_inner_iterations);
377 
378 
379  if (current_minimum > f) {
380  current_minimum = f;
381  memcpy( optimized_parameters_mtx.get_data(), solution_guess.get_data(), num_of_parameters*sizeof(double) );
382  for ( int jdx=0; jdx<num_of_parameters; jdx++) {
383  solution_guess[jdx] = solution_guess[jdx] + distrib_real(gen)/100;
384  }
385  }
386  else {
387  for ( int jdx=0; jdx<num_of_parameters; jdx++) {
388  solution_guess[jdx] = solution_guess[jdx] + distrib_real(gen);
389  }
390  }
391 
392 #ifdef __MPI__
393  MPI_Bcast( (void*)solution_guess.get_data(), num_of_parameters, MPI_DOUBLE, 0, MPI_COMM_WORLD);
394 #endif
395 
396 
397 
398  }
399 
400 
401 }
402 
403 
404 
405 
412 
413  // get the transformed matrix with the gates in the list
414  Matrix_real parameters_mtx(parameters, 1, parameter_num );
415  Matrix matrix_new = Umtx.copy();
416  apply_to( parameters_mtx, matrix_new );
417 
418 #ifdef DEBUG
419  if (matrix_new.isnan()) {
420  std::stringstream sstream;
421  sstream << "Sub_Matrix_Decomposition::optimization_problem: matrix_new contains NaN a. Exiting" << std::endl;
422  print(sstream, 1);
423  }
424 #endif
425 
426  double cost_function = get_submatrix_cost_function(matrix_new);
427 
428 
429  return cost_function;
430 }
431 
432 
440 
441  Sub_Matrix_Decomposition* instance = reinterpret_cast<Sub_Matrix_Decomposition*>(void_instance);
442  std::vector<Gate*> gates_loc = instance->get_gates();
443 
444  Matrix Umtx_loc, matrix_new;
445 
446 //{
447 //tbb::spin_mutex::scoped_lock my_lock{my_mutex};
448 
449  Umtx_loc = instance->get_Umtx();
450  matrix_new = Umtx_loc.copy();
451  instance->apply_to( parameters, matrix_new );
452 
453 #ifdef DEBUG
454  if (matrix_new.isnan()) {
455  std::stringstream sstream;
456  sstream << "Sub_Matrix_Decomposition::optimization_problem matrix_new contains NaN b." << std::endl;
457  print(sstream, 1);
458  }
459 #endif
460 
461 //}
462 
463 
464  double cost_function = get_submatrix_cost_function(matrix_new); //NEW METHOD
465 
466 
467  return cost_function;
468 }
469 
470 
471 
472 
473 
481 
482  // The function value at x0
483  double f0;
484 
485  // calculate the approximate gradient
486  optimization_problem_combined( parameters, void_instance, &f0, grad);
487 
488 }
489 
490 
491 
500 
501  Sub_Matrix_Decomposition* instance = reinterpret_cast<Sub_Matrix_Decomposition*>(void_instance);
502 
503  int parameter_num_loc = instance->get_parameter_num();
504 
505  // storage for the function values calculated at the displaced points x
506  Matrix_real f(1, grad.size());
507 
508  // the difference in one direction in the parameter for the gradient calculation
509  double dparam = 1e-8;
510 
511  // calculate the function values at displaced x and the central x0 points through TBB parallel for
512  tbb::parallel_for(0, parameter_num_loc+1, 1, [&](int i) {
513 
514  if (i == (int)parameters.size()) {
515  // calculate function value at x0
516  *f0 = instance->optimization_problem(parameters, reinterpret_cast<void*>(instance));
517  }
518  else {
519 
520  Matrix_real parameters_d = parameters.copy();
521  parameters_d[i] = parameters_d[i] + dparam;
522 
523  // calculate the cost function at the displaced point
524  f[i] = instance->optimization_problem(parameters_d, reinterpret_cast<void*>(instance));
525 
526  }
527  });
528 
529 
530 /*
531  // sequential version
532  functor_sub_optimization_grad<Sub_Matrix_Decomposition> tmp = functor_grad<Sub_Matrix_Decomposition>( parameters, instance, f, f0, dparam );
533  #pragma omp parallel for
534  for (int idx=0; idx<parameter_num_loc+1; idx++) {
535  tmp(idx);
536  }
537 */
538 
539 
540  for (int idx=0; idx<parameter_num_loc; idx++) {
541  // set the gradient
542 #ifdef DEBUG
543  if (isnan(f->data[idx])) {
544  std::string err("Sub_Matrix_Decomposition::optimization_problem_combined: f->data[i] is NaN.");
545  throw err;
546  }
547 #endif // DEBUG
548  grad[idx] = (f[idx]-(*f0))/dparam;
549  }
550 
551 
552 
553 }
554 
555 
556 
557 
558 
559 
560 
561 
562 
563 
570 int Sub_Matrix_Decomposition::set_identical_blocks( int qbit, int identical_blocks_in ) {
571 
572  std::map<int,int>::iterator key_it = identical_blocks.find( qbit );
573 
574  if ( key_it != identical_blocks.end() ) {
575  identical_blocks.erase( key_it );
576  }
577 
578  identical_blocks.insert( std::pair<int, int>(qbit, identical_blocks_in) );
579 
580  return 0;
581 
582 }
583 
584 
590 int Sub_Matrix_Decomposition::set_identical_blocks( std::map<int, int> identical_blocks_in ) {
591 
592  for ( std::map<int,int>::iterator it=identical_blocks_in.begin(); it!= identical_blocks_in.end(); it++ ) {
593  set_identical_blocks( it->first, it->second );
594  }
595 
596  return 0;
597 
598 }
599 
600 
606 
607 
609 
610  // setting computational parameters
616 
617  if ( extract_gates(static_cast<Gates_block*>(ret)) != 0 ) {
618  std::string err("Sub_Matrix_Decomposition::clone(): extracting gates was not succesfull.");
619  throw err;
620  }
621 
622  return ret;
623 
624 }
625 
626 
bool optimize_layer_num
logical value. Set true to optimize the minimum number of gate layers required in the decomposition...
bool isnan()
Call to check the array for NaN entries.
Definition: matrix.cpp:128
Sub_Matrix_Decomposition()
Nullary constructor of the class.
static void optimization_problem_grad(Matrix_real parameters, void *void_instance, Matrix_real &grad)
Calculate the approximate derivative (f-f0)/(x-x0) of the cost function with respect to the free para...
void set_custom_gate_layers(Gates_block *block_in)
Call to set custom layers to the gate structure that are intended to be used in the subdecomposition...
void print(const std::stringstream &sstream, int verbose_level=1) const
Call to print output messages in the function of the verbosity level.
Definition: logging.cpp:55
Gates_block * unit_gate_structure
Custom gate structure describing the gate structure used in the decomposition. The gate structure is ...
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...
Matrix_real copy() const
Call to create a copy of the matrix.
double current_minimum
The current minimum of the optimization problem.
virtual Gate * clone()
Call to create a clone of the present class.
Definition: Gate.cpp:513
void add_gate(Gate *gate)
Append a general gate to the list of gates.
long long max_inner_iterations
maximal number of inner iterations
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::vector< Gate * > get_gates()
Call to get the gates stored in the class.
int layer_num
number of gate layers
Definition: Gates_block.h:48
void add_cnot(int target_qbit, int control_qbit)
Append a CNOT gate gate to the list of gates.
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...
virtual void apply_to(Matrix_real &parameters_mtx, Matrix &input, int parallel=0)
Call to apply the gate on the input array/matrix Gates_block*input.
scalar * get_data() const
Call to get the pointer to the stored data.
guess_type initial_guess
type to guess the initial values for the optimization. Possible values: ZEROS=0, RANDOM=1, CLOSE_TO_ZERO=2
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...
double get_submatrix_cost_function(Matrix &matrix)
Call to calculate the cost function of a given matrix during the submatrix decomposition process...
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...
std::vector< Gate * > gates
The list of stored gates.
Definition: Gates_block.h:46
int max_outer_iterations
Maximal number of iterations allowed in the optimization process.
A class implementing the BFGS optimizer based on conjugate gradient direction method of M...
Definition: BFGS_Powell.h:26
gate_type type
The type of the operation (see enumeration gate_type)
Definition: Gate.h:83
int get_parameter_num()
Call to get the number of free parameters.
double Start_Optimization(Matrix_real &x, long maximal_iterations_in=5001)
Call this method to start the optimization.
void add_u3(int target_qbit)
Append a U3 gate to the list of gates.
bool subdisentaglement_done
logical value indicating whether the disentamglement of a qubit from the othetrs was done or not ...
virtual void add_gate_layers()
Call to add further layer to the gate structure used in the subdecomposition.
Sub_Matrix_Decomposition * clone()
Call to create a clone of the present class.
std::map< int, int > iteration_loops
A map of <int n: int num> indicating the number of iteration in each step of the decomposition.
~Sub_Matrix_Decomposition()
Destructor of the class.
void add_default_gate_layers()
Call to add default gate layers to the gate structure used in the subdecomposition.
A class containing basic methods for the decomposition process.
int extract_gates(Gates_block *op_block)
Call to extract the gates stored in the class.
virtual Gates_block * clone()
Create a clone of the present class.
Class to store data of complex arrays and its properties.
Definition: matrix.h:38
int size() const
Call to get the number of the allocated elements.
Gates_block()
Default constructor of the class.
Definition: Gates_block.cpp:64
A class responsible for grouping two-qubit (CNOT,CZ,CH) and one-qubit gates into layers.
Definition: Gates_block.h:41
Matrix subdecomposed_mtx
The subdecomposed matrix.
guess_type
Type definition of the types of the initial guess.
std::map< std::string, Config_Element > config
config metadata utilized during the optimization
Base class for the representation of general gate operations.
Definition: Gate.h:73
static void optimization_problem_combined(Matrix_real parameters, void *void_instance, double *f0, Matrix_real &grad)
Call to calculate both the cost function and the its gradient components.
Header file for a class responsible for the disentanglement of one qubit from the others...
int parameter_num
the number of free parameters of the operation
Definition: Gate.h:91
double optimization_problem(double *parameters)
The optimization problem of the final optimization.
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.
Matrix copy()
Call to create a copy of the matrix.
Definition: matrix.cpp:105
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...
int qbit_num
number of qubits spanning the matrix of the operation
Definition: Gate.h:81
void set_max_iteration(int max_outer_iterations_in)
Call to set the maximal number of the iterations in the optimization process.
bool check_optimization_solution()
check_optimization_solution
A class responsible for the disentanglement of one qubit from the others.
int optimization_block
number of gate blocks used in one shot of the optimization process
int set_iteration_loops(int n, int iteration_loops_in)
Set the number of iteration loops during the subdecomposition of the n-th qubit.
Matrix_real optimized_parameters_mtx
The optimized parameters for the gates.
double global_target_minimum
The global target minimum of the optimization problem.
Header file for the paralleized calculation of the cost function of the subdecomposition (supporting ...
void set_optimization_blocks(int optimization_block_in)
Call to set the number of gate blocks to be optimized in one shot.
Class to store data of complex arrays and its properties.
Definition: matrix_real.h:39
std::mt19937 gen
Standard mersenne_twister_engine seeded with rd()
void solve_layer_optimization_problem(int num_of_parameters, Matrix_real solution_guess_gsl)
Call to solve layer by layer the optimization problem.