Sequential Quantum Gate Decomposer  v1.9.3
Powerful decomposition of general unitarias into one- and two-qubit gates gates
Optimization_Interface.h
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 */
23 #ifndef Optimization_Interface_H
24 #define Optimization_Interface_H
25 
26 #include "Decomposition_Base.h"
27 #include "BFGS_Powell.h"
28 #include "Bayes_Opt.h"
29 #include "Powells_method.h"
30 
31 
32 #if defined __DFE__
33  #include "common_DFE.h"
34 #elif defined __GROQ__
35  #include "common_GROQ.h"
36 #endif
37 
38 
39 
44 
45 
46 
49 
50 
51 #ifdef __cplusplus
52 extern "C"
53 {
54 #endif
55 
57 int LAPACKE_zggev ( int matrix_layout,
58  char jobvl,
59  char jobvr,
60  int n,
61  QGD_Complex16 * a,
62  int lda,
63  QGD_Complex16 * b,
64  int ldb,
66  QGD_Complex16 * beta,
67  QGD_Complex16 * vl,
68  int ldvl,
69  QGD_Complex16 * vr,
70  int ldvr
71  );
72 
73 #ifdef __cplusplus
74 }
75 #endif
76 
77 
78 
79 
80 
86 
87 
88 public:
89 
95  int id;
96 protected:
97 
98 
101 
103  std::map<int,int> identical_blocks;
114 
115 
118 
122  double radius;
127 
130 
131 
135  double CPU_time;
136 
137 
138 
139 
140 public:
141 
147 
148 
149 
159 Optimization_Interface( Matrix Umtx_in, int qbit_num_in, bool optimize_layer_num_in, std::map<std::string, Config_Element>& config, guess_type initial_guess_in, int accelerator_num_in=0 );
160 
161 
162 
166 virtual ~Optimization_Interface();
167 
174 
175 
179 virtual void add_finalyzing_layer();
180 
181 
185 void final_optimization();
186 
187 
188 
195 
196 
203 
204 
211 
212 
219 
220 
227 
228 
229 
235 
236 
243 
244 
251 static void export_current_cost_fnc(double current_minimum, Matrix_real& parameters, void* void_instance);
252 
253 
260 
261 
268 
269 
276 
283 
284 
291 
292 
293 
300 
307 
314 void randomize_parameters( Matrix_real& input, Matrix_real& output, const double& f0 );
315 
316 
322 double optimization_problem( double* parameters);
323 
324 
331 
332 
333 
334 #ifdef __DFE__
335 
340 Matrix_real optimization_problem_batched_DFE( std::vector<Matrix_real>& parameters_vec);
341 #endif
342 
343 #ifdef __GROQ__
344 
351 virtual double optimization_problem_Groq( Matrix_real& parameters, int chosen_device);
352 
353 
359 Matrix_real optimization_problem_batched_Groq( std::vector<Matrix_real>& parameters_vec);
360 #endif
361 
362 
368 Matrix_real optimization_problem_batched( std::vector<Matrix_real>& parameters_vec);
369 
370 
371 
379 double optimization_problem( Matrix_real parameters, void* void_instance, Matrix ret_temp);
380 
381 
388 static double optimization_problem( Matrix_real parameters, void* void_instance);
389 
390 
396 virtual double optimization_problem_non_static( Matrix_real parameters, void* void_instance);
397 
398 
405 static void optimization_problem_grad( Matrix_real parameters, void* void_instance, Matrix_real& grad );
406 
407 
408 
409 
417 static void optimization_problem_combined( Matrix_real parameters, void* void_instance, double* f0, Matrix_real& grad );
418 
419 
427 virtual void optimization_problem_combined_non_static( Matrix_real parameters, void* void_instance, double* f0, Matrix_real& grad );
428 
429 
437 virtual void optimization_problem_combined( Matrix_real parameters, double* f0, Matrix_real grad );
438 
439 
447 static void optimization_problem_combined_unitary( Matrix_real parameters, void* void_instance, Matrix& Umtx, std::vector<Matrix>& Umtx_deriv );
448 
455 void optimization_problem_combined_unitary( Matrix_real parameters, Matrix& Umtx, std::vector<Matrix>& Umtx_deriv );
456 
457 
464 double optimization_problem_panelty( double* parameters, Gates_block* gates_block );
465 
466 
471 
472 
477 
478 
484 
485 
490 double get_correction1_scale();
491 
492 
493 
494 
495 
500 double get_correction2_scale();
501 
502 
503 
504 
505 
506 
511 void set_max_inner_iterations( int max_inner_iterations_in );
512 
513 
518 void set_random_shift_count_max( int random_shift_count_max_in );
519 
520 
526 
527 
531 int get_trace_offset();
532 
536 void set_trace_offset(int trace_offset_in);
537 
538 
542 int get_num_iters();
543 
548 void set_custom_gate_structure( Gates_block* gate_structure_in );
549 
550 
551 #ifdef __DFE__
552 
556 void upload_Umtx_to_DFE();
557 
558 
559 #endif
560 
561 
565 int get_accelerator_num();
566 
567 
568 
569 
570 };
571 
572 
573 #endif
optimization_aglorithms alg
The optimization algorithm to be used in the optimization.
bool adaptive_eta
logical variable indicating whether adaptive learning reate is used in the ADAM algorithm ...
int id
unique id indentifying the instance of the class
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...
int get_num_iters()
Get the number of processed iterations during the optimization process.
void set_optimizer(optimization_aglorithms alg_in)
Call to set the optimizer engine to be used in solving the optimization problem.
void solve_layer_optimization_problem_AGENTS(int num_of_parameters, Matrix_real &solution_guess)
Call to solve layer by layer the optimization problem via the AGENT algorithm.
Definition: AGENTS.cpp:42
void set_custom_gate_structure(Gates_block *gate_structure_in)
Call to set custom layers to the gate structure that are intended to be used in the subdecomposition...
std::map< int, int > identical_blocks
A map of <int n: int num> indicating that how many identical successive blocks should be used in the ...
double current_minimum
The current minimum of the optimization problem.
static void optimization_problem_combined_unitary(Matrix_real parameters, void *void_instance, Matrix &Umtx, std::vector< Matrix > &Umtx_deriv)
Call to calculate both the effect of the circuit on th eunitary and it&#39;s gradient componets...
double correction1_scale
prefactor of the single-bitflip errors in the cost function. (see Eq. (21) in arXiv:2210.09191)
double get_correction2_scale()
Call to get the prefactor of the two-bitflip errors in the cost function.
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.
void solve_layer_optimization_problem_GRAD_DESCEND(int num_of_parameters, Matrix_real &solution_guess)
Call to solve layer by layer the optimization problem via the GRAD_DESCEND (line search in the direct...
void set_trace_offset(int trace_offset_in)
Set the trace offset used in the evaluation of the cost function.
double prev_cost_fnv_val
the previous value of the cost funtion to be used to evaluate bitflip errors in the cost funtion (see...
int trace_offset
The offset in the first columns from which the "trace" is calculated. In this case Tr(A) = sum_(i-off...
void set_random_shift_count_max(int random_shift_count_max_in)
Call to set the maximal number of parameter randomization tries to escape a local minimum...
virtual double optimization_problem_non_static(Matrix_real parameters, void *void_instance)
The optimization problem of the final optimization.
double correction2_scale
prefactor of the double-bitflip errors in the cost function. (see Eq. (21) in arXiv:2210.09191)
cost_function_type get_cost_function_variant()
Call to get the variant of the cost function used in the calculations.
int get_trace_offset()
Get the trace ffset used in the evaluation of the cost function.
optimization_aglorithms
implemented optimization strategies
void set_max_inner_iterations(int max_inner_iterations_in)
Call to set the maximal number of iterations for which an optimization engine tries to solve the opti...
void solve_layer_optimization_problem_BAYES_AGENTS(int num_of_parameters, Matrix_real &solution_guess)
Call to solve layer by layer the optimization problem via Bayes & Agents algorithm.
int accelerator_num
number of utilized accelerators
int number_of_iters
number of iterations
virtual void add_finalyzing_layer()
Call to add further layer to the gate structure used in the subdecomposition.
A base class to determine the decomposition of an N-qubit unitary into a sequence of CNOT and U3 gate...
double randomization_rate
randomization rate
void solve_layer_optimization_problem_BFGS2(int num_of_parameters, Matrix_real solution_guess)
Call to solve layer by layer the optimization problem via BBFG algorithm.
Definition: BFGS2.cpp:42
Structure type representing complex numbers in the SQUANDER package.
Definition: QGDTypes.h:38
double optimization_problem_panelty(double *parameters, Gates_block *gates_block)
The optimization problem of the final optimization.
double CPU_time
time spent on optimization
A class containing basic methods for the decomposition process.
double circuit_simulation_time
Time spent on circuit simulation/cost function evaluation.
decomposed_matrix
the unitary matrix from the result object
Definition: example.py:90
void solve_layer_optimization_problem_BFGS(int num_of_parameters, Matrix_real &solution_guess)
Call to solve layer by layer the optimization problem via BBFG algorithm.
Definition: BFGS.cpp:42
Class to store data of complex arrays and its properties.
Definition: matrix.h:38
Header file for a class containing basic methods for the decomposition process.
double get_previous_cost_function_value()
Call to retrieve the previous value of the cost funtion to be used to evaluate bitflip errors in the ...
void calc_decomposition_error(Matrix &decomposed_matrix)
Calculate the error of the decomposition according to the spectral norm of , where is the unitary pr...
cost_function_type
Type definition of the fifferent types of the cost function.
A class responsible for grouping two-qubit (CNOT,CZ,CH) and one-qubit gates into layers.
Definition: Gates_block.h:41
int LAPACKE_zggev(int matrix_layout, char jobvl, char jobvr, int n, QGD_Complex16 *a, int lda, QGD_Complex16 *b, int ldb, QGD_Complex16 *alpha, QGD_Complex16 *beta, QGD_Complex16 *vl, int ldvl, QGD_Complex16 *vr, int ldvr)
Definition of the zggev function from Lapacke to calculate the eigenvalues of a complex matrix...
guess_type
Type definition of the types of the initial guess.
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.
std::map< std::string, Config_Element > config
config metadata utilized during the optimization
void solve_layer_optimization_problem_ADAM(int num_of_parameters, Matrix_real &solution_guess)
Call to solve layer by layer the optimization problem via ADAM algorithm.
static void optimization_problem_grad(Matrix_real parameters, void *void_instance, Matrix_real &grad)
Calculate the derivative of the cost function with respect to the free parameters.
void solve_layer_optimization_problem(int num_of_parameters, Matrix_real solution_guess)
Call to solve layer by layer the optimization problem via calling one of the implemented algorithms...
void solve_layer_optimization_problem_GRAD_DESCEND_PARAMETER_SHIFT_RULE(int num_of_parameters, Matrix_real &solution_guess)
Call to solve layer by layer the optimization problem via the GRAD_DESCEND_PARAMETER_SHIFT_RULE algor...
Matrix Umtx
The unitary to be decomposed.
virtual ~Optimization_Interface()
Destructor of the class.
double get_correction1_scale()
Call to get the prefactor of the single-bitflip errors in the cost function.
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 final_optimization()
final optimization procedure improving the accuracy of the decompositin when all the qubits were alre...
bool optimize_layer_num
logical value. Set true to optimize the minimum number of gate layers required in the decomposition...
Header file for DFE support in unitary simulation.
void randomize_parameters(Matrix_real &input, Matrix_real &output, const double &f0)
Call to randomize the parameter.
int max_inner_iterations
the maximal number of iterations for which an optimization engine tries to solve the optimization pro...
void solve_layer_optimization_problem_AGENTS_COMBINED(int num_of_parameters, Matrix_real &solution_guess)
Call to solve layer by layer the optimization problem via the AGENT COMBINED algorithm.
Definition: AGENTS.cpp:920
int random_shift_count_max
the maximal number of parameter randomization tries to escape a local minimum.
Optimization_Interface()
Nullary constructor of the class.
void solve_layer_optimization_problem_ADAM_BATCHED(int num_of_parameters, Matrix_real &solution_guess_)
Call to solve layer by layer the optimization problem via batched ADAM algorithm. ...
void set_cost_function_variant(cost_function_type variant)
Call to set the variant of the cost function used in the calculations.
double radius
parameter to contron the radius of parameter randomization around the curren tminimum ...
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
virtual void optimization_problem_combined_non_static(Matrix_real parameters, void *void_instance, double *f0, Matrix_real &grad)
Call to calculate both the cost function and the its gradient components.
void solve_layer_optimization_problem_BAYES_OPT(int num_of_parameters, Matrix_real &solution_guess)
Call to solve layer by layer the optimization problem via Bayes algorithm.