Quantum Gate Decomposer  v1.3
Powerful decomposition of almost any unitary into U3 and CNOT gates
Decomposition_Base.h
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 */
24 #ifndef DECOMPOSITION_BASE_H
25 #define DECOMPOSITION_BASE_H
26 
27 
28 #include "qgd/Operation_block.h"
29 #include "qgd/CNOT.h"
30 #include "qgd/U3.h"
31 #include <map>
32 #include <cstdlib>
33 #include <time.h>
34 #include <ctime>
35 #include "gsl/gsl_multimin.h"
36 #include "gsl/gsl_statistics.h"
37 
38 
41 
42 
43 
48 
49 
50 public:
52  bool verbose;
53 
56 
58  static std::map<int,int> max_layer_num_def;
59 
62 
63 protected:
64 
66  std::map<int,int> max_layer_num;
67 
69  std::map<int,int> iteration_loops;
70 
73 
76 
79 
82 
85 
88 
91 
94 
97 
100 
103 
106 
107 public:
108 
116 Decomposition_Base( Matrix Umtx_in, int qbit_num_in, guess_type initial_guess_in);
117 
121 virtual ~Decomposition_Base();
122 
123 
128 void set_optimization_blocks( int optimization_block_in );
129 
134 void set_max_iteration( int max_iterations_in);
135 
136 
141 
142 
147 void list_operations( int start_index );
148 
156 Matrix get_finalizing_operations( Matrix& mtx, Operation_block* finalizing_operations, double* finalizing_parameters);
157 
158 
164 void solve_optimization_problem( double* solution_guess, int solution_guess_num );
165 
171 virtual void solve_layer_optimization_problem( int num_of_parameters, gsl_vector *solution_guess_gsl);
172 
173 
174 
179 virtual double optimization_problem( const double* parameters );
180 
186 
187 
195 std::vector<Matrix, tbb::cache_aligned_allocator<Matrix>> get_operation_products(double* parameters, std::vector<Operation*>::iterator operations_it, int num_of_operations);
196 
197 
202 Matrix get_Umtx();
203 
208 int get_Umtx_size();
209 
214 double* get_optimized_parameters();
215 
220 void get_optimized_parameters( double* ret );
221 
229 Matrix get_transformed_matrix( const double* parameters, std::vector<Operation*>::iterator operations_it, int num_of_operations );
230 
231 
232 
241 Matrix get_transformed_matrix( const double* parameters, std::vector<Operation*>::iterator operations_it, int num_of_operations, Matrix& initial_matrix );
242 
243 
249 
250 
257 Matrix apply_operation( Matrix& operation_mtx, Matrix& input_matrix );
258 
263 void reorder_qubits( std::vector<int> qbit_list);
264 
271 int set_max_layer_num( int n, int max_layer_num_in );
272 
278 int set_max_layer_num( std::map<int, int> max_layer_num_in );
279 
280 
287 int set_iteration_loops( int n, int iteration_loops_in );
288 
294 int set_iteration_loops( std::map<int, int> iteration_loops_in );
295 
296 
300 static void Init_max_layer_num();
301 
302 
307 
314 std::vector<Operation*> prepare_operations_to_export( std::vector<Operation*> ops, const double* parameters );
315 
316 
317 
324 std::vector<Operation*> prepare_operations_to_export( Operation_block* block_op, const double* parameters );
325 
335 int get_operation( unsigned int n, operation_type &type, int &target_qbit, int &control_qbit, double* parameters );
336 
337 
342 void set_verbose( bool verbose_in );
343 
344 
349 double get_decomposition_error( );
350 
351 
352 
353 };
354 
355 #endif //DECOMPOSITION_BASE
Matrix get_finalizing_operations(Matrix &mtx, Operation_block *finalizing_operations, double *finalizing_parameters)
This method determine the operations needed to rotate the indepent qubits into the state |0>
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.
void finalize_decomposition()
After the main optimization problem is solved, the indepent qubits can be rotated into state |0> by t...
void prepare_operations_to_export()
Call to prepare the optimized operations to export.
Matrix get_decomposed_matrix()
Calculate the decomposed matrix resulted by the effect of the optimized operations on the unitary Umt...
int target_qbit
The index of the qubit on which the operation acts (target_qbit >= 0)
Definition: Operation.h:50
int control_qbit
The index of the qubit which acts as a control qubit (control_qbit >= 0) in controlled operations.
Definition: Operation.h:52
bool decomposition_finalized
logical value describing whether the decomposition was finalized or not (i.e. whether the decomposed ...
Matrix get_Umtx()
Call to retrive a pointer to the unitary to be transformed.
double current_minimum
The current minimum of the optimization problem.
bool optimization_problem_solved
logical value describing whether the optimization problem was solved or not
Matrix apply_operation(Matrix &operation_mtx, Matrix &input_matrix)
Apply an operations on the input matrix.
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...
guess_type initial_guess
type to guess the initial values for the optimization. Possible values: ZEROS=0, RANDOM=1,...
static void Init_max_layer_num()
Initializes default layer numbers.
std::vector< Matrix, tbb::cache_aligned_allocator< Matrix > > get_operation_products(double *parameters, std::vector< Operation * >::iterator operations_it, int num_of_operations)
Calculate the list of gate operation matrices such that the i>0-th element in the result list is the ...
int finalizing_operations_num
number of finalizing (deterministic) opertaions rotating the disentangled qubits into state |0>.
double * get_optimized_parameters()
Call to get the optimized parameters.
Header file for a class representing a CNOT operation.
Decomposition_Base(Matrix Umtx_in, int qbit_num_in, guess_type initial_guess_in)
Contructor of the class.
double optimization_tolerance
The maximal allowed error of the optimization problem.
operation_type type
The type of the operation (see enumeration operation_type)
Definition: Operation.h:48
double get_decomposition_error()
Call to get the error of the decomposition.
int get_Umtx_size()
Call to get the size of the unitary to be transformed.
virtual void solve_layer_optimization_problem(int num_of_parameters, gsl_vector *solution_guess_gsl)
Abstarct function to be used to solve a single sub-layer optimization problem.
A class responsible for grouping CNOT and U3 operations into layers.
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.
int num_threads
Number of outer OpenMP threads. (During the calculations OpenMP multithreading is turned off....
A class containing basic methods for the decomposition process.
Class to store data of complex arrays and its properties.
Definition: matrix.h:12
virtual ~Decomposition_Base()
Destructor of the class.
Header file for a class representing a U3 operation.
void reorder_qubits(std::vector< int > qbit_list)
Call to reorder the qubits in the matrix of the operation.
bool verbose
Logical variable. Set true for verbose mode, or to false to suppress output messages.
guess_type
Type definition of the types of the initial guess.
double * optimized_parameters
The optimized parameters for the operations.
int get_operation(unsigned int n, operation_type &type, int &target_qbit, int &control_qbit, double *parameters)
Call to prepare the optimized operations to export.
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...
void set_verbose(bool verbose_in)
Call to set the verbose attribute to true or false.
int finalizing_parameter_num
the number of the finalizing (deterministic) parameters of operations rotating the disentangled qubit...
bool check_optimization_solution()
check_optimization_solution
int optimization_block
number of operation blocks used in one shot of the optimization process
Header file for a class responsible for grouping CNOT and U3 operations into layers.
double decomposition_error
error of the unitarity of the final decomposition
int set_iteration_loops(int n, int iteration_loops_in)
Set the number of iteration loops during the subdecomposition of the n-th qubit.
operation_type
Type definition of operation types (also generalized for decomposition classes derived from the class...
Definition: Operation.h:33
void list_operations(int start_index)
Call to print the operations decomposing the initial unitary.
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.
virtual double optimization_problem(const double *parameters)
This is an abstact definition of function giving the cost functions measuring the entaglement of the ...