Sequential Quantum Gate Decomposer  v1.9.3
Powerful decomposition of general unitarias into one- and two-qubit gates gates
COSINE.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 */
24 #include "Optimization_Interface.h"
26 
27 
28 #include <fstream>
29 
30 
31 #ifdef __DFE__
32 #include "common_DFE.h"
33 #endif
34 
35 
36 
37 
44 
45 
46  if ( cost_fnc != FROBENIUS_NORM && cost_fnc != VQE ) {
47  std::string err("Optimization_Interface::solve_layer_optimization_problem_COSINE: Only cost functions FROBENIUS_NORM and VQE are implemented for this strategy");
48  throw err;
49  }
50 
51 
52 #ifdef __DFE__
53  if ( qbit_num >= 5 && get_accelerator_num() > 0 ) {
54  upload_Umtx_to_DFE();
55  }
56 #endif
57 
58  tbb::tick_count t0_CPU = tbb::tick_count::now();
59 
60  if (gates.size() == 0 ) {
61  return;
62  }
63 
64 
65  double M_PI_quarter = M_PI/4;
66  double M_PI_half = M_PI/2;
67  double M_PI_double = M_PI*2;
68 
69  if (solution_guess.size() == 0 ) {
70  solution_guess = Matrix_real(num_of_parameters,1);
71  std::uniform_real_distribution<> distrib_real(0, M_PI_double);
72  for ( int idx=0; idx<num_of_parameters; idx++) {
73  solution_guess[idx] = distrib_real(gen);
74  }
75 
76  }
77 
78 
79 
80  if (optimized_parameters_mtx.size() == 0) {
81  optimized_parameters_mtx = Matrix_real(1, num_of_parameters);
82  memcpy(optimized_parameters_mtx.get_data(), solution_guess.get_data(), num_of_parameters*sizeof(double) );
83  }
84 
85 
86 
87 
88  std::stringstream sstream;
89  double optimization_time = 0.0;
90  tbb::tick_count optimization_start = tbb::tick_count::now();
91 
92 
93  // the current result
95 
96 
97  // the array storing the optimized parameters
98  Matrix_real solution_guess_tmp_mtx = Matrix_real( num_of_parameters, 1 );
99  memcpy(solution_guess_tmp_mtx.get_data(), optimized_parameters_mtx.get_data(), num_of_parameters*sizeof(double) );
100 
101  int batch_size;
102  if ( config.count("batch_size_cosine") > 0 ) {
103  long long value;
104  config["batch_size_cosine"].get_property( value );
105  batch_size = (int) value;
106  }
107  else if ( config.count("batch_size") > 0 ) {
108  long long value;
109  config["batch_size"].get_property( value );
110  batch_size = (int) value;
111  }
112  else {
113  batch_size = 64 <= num_of_parameters ? 64 : num_of_parameters;
114 
115  }
116 
117  if( batch_size > num_of_parameters ) {
118  std::string err("Optimization_Interface::solve_layer_optimization_problem_COSINE: batch size should be lower or equal to the number of free parameters");
119  throw err;
120  }
121 
122 
123  long long max_inner_iterations_loc;
124  if ( config.count("max_inner_iterations_cosine") > 0 ) {
125  config["max_inner_iterations_cosine"].get_property( max_inner_iterations_loc );
126  }
127  else if ( config.count("max_inner_iterations") > 0 ) {
128  config["max_inner_iterations"].get_property( max_inner_iterations_loc );
129  }
130  else {
131  max_inner_iterations_loc = max_inner_iterations;
132  }
133 
134 
135  long long export_circuit_2_binary_loc;
136  if ( config.count("export_circuit_2_binary_cosine") > 0 ) {
137  config["export_circuit_2_binary_cosine"].get_property( export_circuit_2_binary_loc );
138  }
139  else if ( config.count("export_circuit_2_binary") > 0 ) {
140  config["export_circuit_2_binary"].get_property( export_circuit_2_binary_loc );
141  }
142  else {
143  export_circuit_2_binary_loc = 0;
144  }
145 
146 
147  double optimization_tolerance_loc;
148  if ( config.count("optimization_tolerance_cosine") > 0 ) {
149  config["optimization_tolerance_cosine"].get_property( optimization_tolerance_loc );
150  }
151  else if ( config.count("optimization_tolerance") > 0 ) {
152  double value;
153  config["optimization_tolerance"].get_property( optimization_tolerance_loc );
154  }
155  else {
156  optimization_tolerance_loc = optimization_tolerance;
157  }
158 
159 
160  // The number if iterations after which the current results are displed/exported
161  int output_periodicity;
162  if ( config.count("output_periodicity_cosine") > 0 ) {
163  long long value = 1;
164  config["output_periodicity_cosine"].get_property( value );
165  output_periodicity = (int) value;
166  }
167  if ( config.count("output_periodicity") > 0 ) {
168  long long value = 1;
169  config["output_periodicity"].get_property( value );
170  output_periodicity = (int) value;
171  }
172  else {
173  output_periodicity = 0;
174  }
175 
176 
177 
178  if ( output_periodicity>0 ) {
180  }
181 
182  // vector stroing the lates values of current minimums to identify convergence
183  Matrix_real f0_vec(1, 100);
184  memset( f0_vec.get_data(), 0.0, f0_vec.size()*sizeof(double) );
185  double f0_mean = 0.0;
186  int f0_idx = 0;
187 
188 
189 
190  Matrix_real param_update_mtx( batch_size, 1 );
191  matrix_base<int> param_idx_agents( batch_size, 1 );
192 
193 
194  std::vector<Matrix_real> parameters_mtx_vec(batch_size);
195  parameters_mtx_vec.reserve(batch_size);
196 
197 
198 
199 
200  bool three_point_line_search = cost_fnc == FROBENIUS_NORM;
201  bool three_point_line_search_double_period = cost_fnc == VQE;
202 
203 
204  for (unsigned long long iter_idx=0; iter_idx<max_inner_iterations_loc; iter_idx++) {
205 
206  // build up a vector of indices providing a set from which we can draw random (but unique) choices
207  std::vector<int> indices(num_of_parameters);
208  indices.reserve(num_of_parameters);
209  for( int idx=0; idx<num_of_parameters; idx++ ) {
210  indices[idx] = idx;
211  }
212 
213  for( int idx=0; idx<batch_size; idx++ ) {
214  parameters_mtx_vec[idx] = solution_guess_tmp_mtx.copy();
215 
216  // random generator of integers
217  std::uniform_int_distribution<> distrib_int(0, indices.size()-1);
218 
219  // The index array of the chosen parameters
220  int chosen_idx = distrib_int(gen);
221  param_idx_agents[ idx ] = indices[ chosen_idx ];
222  indices.erase( indices.begin()+chosen_idx );
223  }
224 
225 
226 #ifdef __MPI__
227  MPI_Bcast( (void*)param_idx_agents.get_data(), batch_size, MPI_INT, 0, MPI_COMM_WORLD);
228 #endif
229 
230 
231 
232  if ( three_point_line_search ) {
233 
234  for(int idx=0; idx<batch_size; idx++) {
235  Matrix_real& solution_guess_mtx_idx = parameters_mtx_vec[ idx ];
236  solution_guess_mtx_idx[ param_idx_agents[idx] ] += M_PI_half;
237  }
238 
239  Matrix_real f0_shifted_pi2_agents = optimization_problem_batched( parameters_mtx_vec );
240 
241 
242  for(int idx=0; idx<batch_size; idx++) {
243  Matrix_real& solution_guess_mtx_idx = parameters_mtx_vec[ idx ];
244  solution_guess_mtx_idx[ param_idx_agents[idx] ] += M_PI_half;
245  }
246 
247  Matrix_real f0_shifted_pi_agents = optimization_problem_batched( parameters_mtx_vec );
248 
249  for( int idx=0; idx<batch_size; idx++ ) {
250 
251  double f0_shifted_pi = f0_shifted_pi_agents[idx];
252  double f0_shifted_pi2 = f0_shifted_pi2_agents[idx];
253 
254 
255  double A_times_cos = (current_minimum-f0_shifted_pi)/2;
256  double offset = (current_minimum+f0_shifted_pi)/2;
257 
258  double A_times_sin = offset - f0_shifted_pi2;
259 
260  double phi0 = atan2( A_times_sin, A_times_cos);
261 
262 
263  double parameter_shift = phi0 > 0 ? M_PI-phi0 : -phi0-M_PI;
264 
265 
266  param_update_mtx[ idx ] = parameter_shift;
267 
268  //revert the changed parameters
269  Matrix_real& solution_guess_mtx_idx = parameters_mtx_vec[idx];
270  solution_guess_mtx_idx[ param_idx_agents[idx] ] = solution_guess_tmp_mtx[ param_idx_agents[idx] ];
271 
272  }
273 
274  }
275  else if ( three_point_line_search_double_period ) {
276 
277  for(int idx=0; idx<batch_size; idx++) {
278  Matrix_real& solution_guess_mtx_idx = parameters_mtx_vec[ idx ];
279  solution_guess_mtx_idx[ param_idx_agents[idx] ] += M_PI_quarter;
280  }
281 
282  Matrix_real f0_shifted_pi4_agents = optimization_problem_batched( parameters_mtx_vec );
283 
284 
285  for(int idx=0; idx<batch_size; idx++) {
286  Matrix_real& solution_guess_mtx_idx = parameters_mtx_vec[ idx ];
287  solution_guess_mtx_idx[ param_idx_agents[idx] ] += M_PI_quarter;
288  }
289 
290  Matrix_real f0_shifted_pi2_agents = optimization_problem_batched( parameters_mtx_vec );
291 
292  for( int idx=0; idx<batch_size; idx++ ) {
293 
294  double f0_shifted_pi = f0_shifted_pi2_agents[idx];
295  double f0_shifted_pi2 = f0_shifted_pi4_agents[idx];
296 
297 
298  double A_times_cos = (current_minimum-f0_shifted_pi)/2;
299  double offset = (current_minimum+f0_shifted_pi)/2;
300 
301  double A_times_sin = offset - f0_shifted_pi2;
302 
303  double phi0 = atan2( A_times_sin, A_times_cos);
304 
305 
306  double parameter_shift = phi0 > 0 ? M_PI_half-phi0/2 : -phi0/2-M_PI_half;
307 
308 
309  param_update_mtx[ idx ] = parameter_shift;
310 
311  //revert the changed parameters
312  Matrix_real& solution_guess_mtx_idx = parameters_mtx_vec[idx];
313  solution_guess_mtx_idx[ param_idx_agents[idx] ] = solution_guess_tmp_mtx[ param_idx_agents[idx] ];
314 
315  }
316 
317 
318 
319  }
320  else {
321  std::string err("solve_layer_optimization_problem_COSINE: Not implemented method.");
322  throw err;
323  }
324 
325 
327  // the line search is converted onto a onediemnsional search between x0+a*param_update_mtx and x0+b*param_update_mtx
328  double interval_coeff = 2.0/(sqrt(5.0) + 1); // = 1/tau in Fig 1 of https://ieeexplore.ieee.org/stamp/stamp.jsp?tp=&arnumber=1606817
329 
330 
331  // lowest point on the search interval
332  double point_a = 0.0;
333  double point_b = 1.0;
334  double epsilon = 1e-2;
335 
336  double interval_length = point_b - point_a;
337 
338  double x1 = point_a + interval_length*interval_coeff*interval_coeff;
339  double x2; // = point_a + interval_length*interval_coeff;
340 
341  Matrix_real parameters_x1 = solution_guess_tmp_mtx.copy();
342  Matrix_real parameters_x2;
343  for( int idx=0; idx<batch_size; idx++ ) {
344  parameters_x1[ param_idx_agents[idx] ] += param_update_mtx[ idx ]*x1;
345  }
346 
347 
348  double val_x1 = optimization_problem( parameters_x1 );
349  double val_x2;
350 
351 
352 
353  int iter_max = 50;
354  int iter = 0;
355 
356  double current_best_point = point_a;
357  double current_best_value = current_minimum;
358 
359 
360  while( true ) {
361 
362  interval_length = point_b - point_a;
363 
364  if ( interval_length < epsilon) {
365  break;
366  }
367 
368  if ( (x1-point_a) < (point_b-x1) ) {
369  // x1 is closer to "a" than to "b" --> x2 should be closer to "b"
370  x2 = point_a + interval_length*interval_coeff;
371 
372  parameters_x2 = solution_guess_tmp_mtx.copy();
373 
374  for( int idx=0; idx<batch_size; idx++ ) {
375  parameters_x2[ param_idx_agents[idx] ] += param_update_mtx[ idx ]*x2;
376  }
377 
378  val_x2 = optimization_problem( parameters_x2 );
379  }
380  else {
381  // x1 should be always closer to "a"
382  x2 = x1;
383  val_x2 = val_x1;
384 
385  x1 = point_a + interval_length*interval_coeff*interval_coeff;
386 
387  parameters_x1 = solution_guess_tmp_mtx.copy();
388 
389  for( int idx=0; idx<batch_size; idx++ ) {
390  parameters_x1[ param_idx_agents[idx] ] += param_update_mtx[ idx ]*x1;
391  }
392 
393  val_x1 = optimization_problem( parameters_x1 );
394 
395  }
396 
397 
398  //std::cout << point_a << " " << x1 << " " << x2 << " " << point_b << std::endl;
399  //std::cout << val_x1 << " " << val_x2 << " " << current_minimum << " " << std::endl;
400 
401  if ( val_x1 < val_x2 ) {
402  point_b = x2;
403  if ( current_best_value > val_x1 ) {
404  current_best_point = x1;
405  current_best_value = val_x1;
406  }
407  }
408  else {
409  point_a = x1;
410  if ( current_best_value > val_x2 ) {
411  current_best_point = x2;
412  current_best_value = val_x2;
413  }
414 
415  x1 = x2;
416  val_x1 = val_x2;
417  }
418 
419 
420  iter = iter + 1;
421 
422  if ( iter > iter_max) {
423  std::cout << "line search not converged: interval length: " << interval_length << " " << interval_coeff << std::endl;
424  break;
425  }
426  }
427 
428  //std::cout << "number of costfunction evaluations : " << iter+1 << " " << interval_coeff << std::endl;
429 
430  current_minimum = current_best_value;
431 
432  // update parameters
433  for (int param_idx=0; param_idx<batch_size; param_idx++) {
434  solution_guess_tmp_mtx[ param_idx_agents[param_idx] ] += param_update_mtx[ param_idx ]*current_best_point;
435  }
436 
438 /*
439  // parameters for line search
440  int line_points = 128;
441 
442  std::vector<Matrix_real> parameters_line_search_mtx_vec(line_points);
443  parameters_line_search_mtx_vec.reserve(line_points);
444 
445  // perform line search over the deriction determined previously
446  for( int line_idx=0; line_idx<line_points; line_idx++ ) {
447 
448  Matrix_real parameters_line_idx = solution_guess_tmp_mtx.copy();
449 
450  for( int idx=0; idx<batch_size; idx++ ) {
451  parameters_line_idx[ param_idx_agents[idx] ] += param_update_mtx[ idx ]*(double)line_idx/line_points;
452  }
453 
454  parameters_line_search_mtx_vec[ line_idx] = parameters_line_idx;
455 
456  }
457 
458 
459  Matrix_real line_values = optimization_problem_batched( parameters_line_search_mtx_vec );
460 
461 
462  // find the smallest value
463  double f0_min = line_values[0];
464  int idx_min = 0;
465  for (int idx=1; idx<line_points; idx++) {
466  if ( line_values[idx] < f0_min ) {
467  idx_min = idx;
468  f0_min = line_values[idx];
469  }
470  }
471 
472  current_minimum = f0_min;
473 
474  // update parameters
475  for (int param_idx=0; param_idx<batch_size; param_idx++) {
476  solution_guess_tmp_mtx[ param_idx_agents[param_idx] ] += param_update_mtx[ param_idx ]*(double)idx_min/line_points;
477  }
478 */
479 #ifdef __MPI__
480  MPI_Bcast( solution_guess_tmp_mtx.get_data(), num_of_parameters, MPI_DOUBLE, 0, MPI_COMM_WORLD);
481 #endif
482 
483 
484  // update the current cost function
485  //current_minimum = optimization_problem( solution_guess_tmp_mtx );
486 
487  if ( output_periodicity>0 && iter_idx % output_periodicity == 0 ) {
488  std::stringstream sstream;
489  sstream << "COSINE: processed iterations " << (double)iter_idx/max_inner_iterations_loc*100 << "\%, current minimum:" << current_minimum;
490  sstream << " circuit simulation time: " << circuit_simulation_time << std::endl;
491  print(sstream, 0);
492  if ( export_circuit_2_binary_loc > 0 ) {
493  std::string filename("initial_circuit_iteration.binary");
494  if (project_name != "") {
495  filename=project_name+ "_" +filename;
496  }
497  export_gate_list_to_binary(solution_guess_tmp_mtx, this, filename, verbose);
498  }
499 
500  memcpy( optimized_parameters_mtx.get_data(), solution_guess_tmp_mtx.get_data(), num_of_parameters*sizeof(double) );
501 
502  if ( output_periodicity>0 && iter_idx % output_periodicity == 0 ) {
504  }
505 
506 
507  }
508 
509  if (current_minimum < optimization_tolerance_loc ) {
510  break;
511  }
512 
513 
514  // test local minimum convergence
515  f0_mean = f0_mean + (current_minimum - f0_vec[ f0_idx ])/f0_vec.size();
516  f0_vec[ f0_idx ] = current_minimum;
517  f0_idx = (f0_idx + 1) % f0_vec.size();
518 
519  double var_f0 = 0.0;
520  for (int idx=0; idx<f0_vec.size(); idx++) {
521  var_f0 = var_f0 + (f0_vec[idx]-f0_mean)*(f0_vec[idx]-f0_mean);
522  }
523  var_f0 = std::sqrt(var_f0)/f0_vec.size();
524 
525 
526 
527  if ( std::abs( f0_mean - current_minimum) < 1e-7 && var_f0/f0_mean < 1e-7 ) {
528  std::stringstream sstream;
529  sstream << "COSINE: converged to minimum at iterations " << (double)iter_idx/max_inner_iterations_loc*100 << "\%, current minimum:" << current_minimum;
530  sstream << " circuit simulation time: " << circuit_simulation_time << std::endl;
531  print(sstream, 0);
532  if ( export_circuit_2_binary_loc > 0 ) {
533  std::string filename("initial_circuit_iteration.binary");
534  if (project_name != "") {
535  filename=project_name+ "_" +filename;
536  }
537  export_gate_list_to_binary(solution_guess_tmp_mtx, this, filename, verbose);
538  }
539 
540 
541  break;
542  }
543 
544  }
545 
546 
547 
548  memcpy( optimized_parameters_mtx.get_data(), solution_guess_tmp_mtx.get_data(), num_of_parameters*sizeof(double) );
549 
550  // CPU time
551  CPU_time += (tbb::tick_count::now() - t0_CPU).seconds();
552 
553  sstream.str("");
554  sstream << "obtained minimum: " << current_minimum << std::endl;
555 
556 
557  tbb::tick_count optimization_end = tbb::tick_count::now();
558  optimization_time = optimization_time + (optimization_end-optimization_start).seconds();
559  sstream << "COSINE time: " << CPU_time << " seconds, obtained minimum: " << current_minimum << std::endl;
560 
561  print(sstream, 0);
562 
563 }
564 
565 
void export_current_cost_fnc(double current_minimum)
Call to print out into a file the current cost function and the second Rényi entropy on the subsyste...
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
Matrix_real copy() const
Call to create a copy of the matrix.
double current_minimum
The current minimum of the optimization problem.
cost_function_type cost_fnc
The chosen variant of the cost function.
int get_accelerator_num()
Get the number of accelerators to be reserved on DFEs on users demand.
double optimization_problem(double *parameters)
Evaluate the optimization problem of the optimization.
scalar * get_data() const
Call to get the pointer to the stored data.
std::vector< Gate * > gates
The list of stored gates.
Definition: Gates_block.h:46
std::string project_name
the name of the project
double optimization_tolerance
The maximal allowed error of the optimization problem (The error of the decomposition would scale wit...
double CPU_time
time spent on optimization
int verbose
Set the verbosity level of the output messages.
Definition: logging.h:50
double circuit_simulation_time
Time spent on circuit simulation/cost function evaluation.
int size() const
Call to get the number of the allocated elements.
std::map< std::string, Config_Element > config
config metadata utilized during the optimization
Header file for the paralleized calculation of the cost function of the final optimization problem (s...
void solve_layer_optimization_problem_COSINE(int num_of_parameters, Matrix_real &solution_guess)
Call to solve layer by layer the optimization problem via the COSINE algorithm.
Definition: COSINE.cpp:43
void export_gate_list_to_binary(Matrix_real &parameters, Gates_block *gates_block, const std::string &filename, int verbosity)
?????????
int qbit_num
number of qubits spanning the matrix of the operation
Definition: Gate.h:81
Header file for DFE support in unitary simulation.
int max_inner_iterations
the maximal number of iterations for which an optimization engine tries to solve the optimization pro...
Matrix_real optimized_parameters_mtx
The optimized parameters for the gates.
Matrix_real optimization_problem_batched(std::vector< Matrix_real > &parameters_vec)
The cost function of the optimization with batched input (implemented only for the Frobenius norm cos...
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()