QQuantLib.AE
ae_class
This module contains a general class for solving AE problems using the algorithm classes from QQuantLib.AE library package Authors: Alberto Pedro Manzano Herrero & Gonzalo Ferro
- class tnbs.BTC_02_AE.QQuantLib.AE.ae_class.AE(oracle=None, target=None, index=None, ae_type=None, **kwargs)
- Bases: - object- Class for creating and solving an AE problem - Parameters:
- oracle (QLM gate) – QLM gate with the Oracle for implementing the Grover operator 
- target (list of ints) – python list with the target for the amplitude estimation 
- index (list of ints) – qubits which mark the register to do the amplitude estimation 
- ae_type (string) – string with the desired AE algorithm: MLAE, CQPEAE, IQPEAE, IQAE, RQAE 
- kwars (dictionary) – dictionary that allows the configuration of the AE algorithm. The different configration keys of the different AE algorithms can be provided. 
 
 - property ae_type
- creating ae_type property 
 - create_ae_solver()
- Method for instantiate the AE algorithm class. 
 - run()
- Method for running an AE problem 
 
ae_classical_qpe
This module contains the CQPEAE class. Given a quantum oracle operator, this class estimates the amplitude of a given target state using the classical QPE algorithm (with QFT). This module uses the QQuantLib.PE.classical_qpe one.
The following references were used:
Brassard, G., Hoyer, P., Mosca, M., & Tapp, A. (2000). Quantum amplitude amplification and estimation. AMS Contemporary Mathematics Series, 305. https://arxiv.org/abs/quant-ph/0005055v1
NEASQC deliverable: D5.1: Review of state-of-the-art for Pricing and Computation of VaR
Author: Gonzalo Ferro Costas & Alberto Manzano Herrero
- class tnbs.BTC_02_AE.QQuantLib.AE.ae_classical_qpe.CQPEAE(oracle: qat.lang.AQASM.QRoutine, target: list, index: list, **kwargs)
- Bases: - object- Class for doing Amplitude Estimation (AE) using classical Quantum Amplitude Estimation (with QFT) algorithm - Parameters:
- oracle (QLM gate) – QLM gate with the Oracle for implementing the Grover operator 
- target (list of ints) – python list with the target for the amplitude estimation 
- index (list of ints) – qubits which mark the register to do the amplitude estimation 
- kwars (dictionary) – - dictionary that allows the configuration of the CQPEAE algorithm: Implemented keys: - qpuQLM solver
- solver for simulating the resulting circuits 
- shotsint
- number of measurements 
- mcz_qlmbool
- for using or not QLM implementation of the multi controlled Z gate 
 
 
 - property index
- creating index property 
 - property oracle
- creating oracle property 
 - run()
- run method for the class. - Returns:
- the estimation of a 
- Return type:
- result 
 - Notes \[a = \cos^2(\theta)\]- Where \(\theta\) is: \[\mathcal{Q}|\Psi\rangle = e^{2i\theta}|\Psi\rangle\]- And \(\mathcal{Q}\) the Grover Operator 
 - property target
- creating target property 
 
ae_iterative_quantum_pe
This module contains necessary functions and classes to implement amplitude estimation algorithm using Iterative Quantum Phase Estimation (IQPE). The implementation is based on following paper:
Dobšíček, Miroslav and Johansson, Göran and Shumeiko, Vitaly and Wendin, Göran*. Arbitrary accuracy iterative quantum phase estimation algorithm using a single ancillary qubit: A two-qubit benchmark. Physical Review A 3(76), 2007. https://arxiv.org/abs/quant-ph/0610214
Author: Gonzalo Ferro Costas & Alberto Manzano Herrero
- class tnbs.BTC_02_AE.QQuantLib.AE.ae_iterative_quantum_pe.IQPEAE(oracle: qat.lang.AQASM.QRoutine, target: list, index: list, **kwargs)
- Bases: - object- Class for using Iterative Quantum Phase Estimation (IQPE) class for doing Amplitude Estimation (AE) - Parameters:
- oracle (QLM gate) – QLM gate with the Oracle for implementing the Grover operator 
- target (list of ints) – python list with the target for the amplitude estimation 
- index (list of ints) – qubits which mark the register to do the amplitude estimation 
- kwars (dictionary) – - dictionary that allows the configuration of the IQPEAE algorithm: Implemented keys: - qpuQLM solver
- solver for simulating the resulting circutis 
- shotsint
- number of measurements 
- mcz_qlmbool
- for using or not QLM implementation of the multi controlled Z gate 
 
 
 - property index
- creating index property 
 - property oracle
- creating oracle property 
 - run()
- run method for the class. - Returns:
- the estimation of a 
- Return type:
- result 
 - Notes \[a = \cos^2(\theta)\]- Where \(\theta\) is: \[\mathcal{Q}|\Psi\rangle = e^{2i\theta}|\Psi\rangle\]- And \(\mathcal{Q}\) the Grover Operator 
 - property target
- creating target property 
 
iterative_quantum_ae
This module contains necessary functions and classes to implement Iterative Quantum Amplitude Estimation based on the paper:
Grinko, D., Gacon, J., Zoufal, C. et al. Iterative Quantum Amplitude Estimation npj Quantum Inf 7, 52 (2021). https://doi.org/10.1038/s41534-021-00379-1
Author: Gonzalo Ferro Costas & Alberto Manzano Herrero
- class tnbs.BTC_02_AE.QQuantLib.AE.iterative_quantum_ae.IQAE(oracle: qat.lang.AQASM.QRoutine, target: list, index: list, **kwargs)
- Bases: - object- Class for Iterative Quantum Amplitude Estimation (IQAE) algorithm - Parameters:
- oracle (QLM gate) – QLM gate with the Oracle for implementing the Grover operator 
- target (list of ints) – python list with the target for the amplitude estimation 
- index (list of ints) – qubits which mark the register to do the amplitude estimation 
- kwars (dictionary) – - dictionary that allows the configuration of the IQAE algorithm: Implemented keys: - qpuQLM solver
- solver for simulating the resulting circuits 
- epsilonfloat
- precision 
- alphafloat
- accuracy 
- shotsint
- number of measurements on each iteration 
- mcz_qlmbool
- for using or not QLM implementation of the multi controlled Z gate 
 
 
 - static chebysev_bound(n_samples: int, gamma: float)
- Computes the length of the confidence interval for a given number of samples n_samples and an accuracy gamma: \[\epsilon = \dfrac{1}{\sqrt{2N}}\log\left(\dfrac{2}{\gamma} \ \right)\]- Parameters:
- n_samples (int) – number of samples 
- gamma (float) – accuracy 
 
- Returns:
- length of the confidence interval 
- Return type:
- float 
 
 - static compute_info(epsilon: float = 0.01, shots: int = 100, alpha: float = 0.05)
- This function computes theoretical values of the IQAE algorithm. - Parameters:
- epsilon (float) – precision 
- alpha (float) – accuracy 
- shots (int) – number of measurements on each iteration 
 
- Returns:
- info – python dictionary with the computed information 
- Return type:
- dict 
 
 - static display_information(epsilon: float = 0.01, shots: int = 100, alpha: float = 0.05)
- This function displays information of the properties of the method for a given set of parameters - Parameters:
- epsilon (float) – precision 
- alpha (float) – accuracy 
- shots (int) – number of measurements on each iteration 
 
 
 - static find_next_k(k: int, theta_lower: float, theta_upper: float, flag: bool, ratio: float = 2)
- This is an implementation of Algorithm 2 from the IQAE paper. This function computes the next suitable k. - Parameters:
- k (int) – number of times to apply the grover operator to the quantum circuit 
- theta_lower (float) – lower bound for the estimation of the angle 
- theta_upper (float) – upper bound for the estimation of the angle 
- flag (bool) – flag to keep track of weather we are in the upper or lower half pane 
- ratio (float) – ratio of amplifications between consecutive iterations 
 
- Returns:
- k (int) – number of times to apply the grover operator to the quantum circuit 
- flag (bool) – flag to keep track of weather we are in the upper or lower half pane 
 
 
 - property index
- creating index property 
 - static invert_sector(a_min: float, a_max: float, flag: bool = True)
- This function inverts the expression: \[a = \dfrac{1-\cos(\theta)}{2}\]- for a pair of bounds (a_min,a_max). The result belongs to the domain (0,2pi) - Parameters:
- a_min (float) – lower bound 
- a_max (float) – upper bound 
- flag (bool) – flag to keep track of weather we are in the upper or lower half pane 
 
- Returns:
- theta_min (float) – lower bound for the associated angle 
- theta_max (float) – upper bound for the associated angle 
 
 
 - iqae(epsilon: float = 0.01, shots: int = 100, alpha: float = 0.05)
- This function implements Algorithm 1 from the IQAE paper. The result is an estimation of the desired probability with precision at least epsilon and accuracy at least alpha. - Parameters:
- epsilon (float) – precision 
- alpha (float) – accuracy 
- shots (int) – number of measurements on each iteration 
 
- Returns:
- a_l (float) – lower bound for the probability to be estimated 
- a_u (float) – upper bound for the probability to be estimated 
 
 
 - property oracle
- creating oracle property 
 - quantum_step(k)
- Create the quantum routine needed for the iqae step - Parameters:
- k (int) – number of Grover operator applications 
- Returns:
- routine – qlm routine for the iqae step 
- Return type:
- qlm routine 
 
 - run()
- run method for the class. - Returns:
- amplitude estimation parameter 
- Return type:
- self.ae 
 
 - property target
- creating target property 
 
maximum_likelihood_ae
This module contains necessary functions and classes to implement Maximum Likelihood Amplitude Estimation based on the paper:
Suzuki, Y., Uno, S., Raymond, R., Tanaka, T., Onodera, T., & Yamamoto, N. Amplitude estimation without phase estimation Quantum Information Processing, 19(2), 2020 arXiv: quant-ph/1904.10246v2
Author: Gonzalo Ferro Costas & Alberto Manzano Herrero
- class tnbs.BTC_02_AE.QQuantLib.AE.maximum_likelihood_ae.MLAE(oracle: qat.lang.AQASM.QRoutine, target: list, index: list, **kwargs)
- Bases: - object- Class for using Maximum Likelihood Quantum Amplitude Estimation (MLAE) algorithm - Parameters:
- oracle (QLM gate) – QLM gate with the Oracle for implementing the Grover operator: init_q_prog and q_gate will be interpreted as None 
- target (list of ints) – python list with the target for the amplitude estimation 
- index (list of ints) – qubits which mark the register to do the amplitude estimation 
- kwars (dictionary) – - dictionary that allows the configuration of the MLAE algorithm: Implemented keys: - qpuQLM solver
- solver for simulating the resulting circuits 
- schedulelist of two lists
- the schedule for the algorithm 
- optimizer :
- an optimizer with just one possible entry 
- deltafloat
- tolerance to avoid division by zero warnings 
- nsint
- number of grid points for brute scipy optimizer 
- mcz_qlmbool
- for using or not QLM implementation of the multi controlled Z gate 
 
 
 - static cost_function(angle: float, m_k: list, n_k: list, h_k: list) float
- This method calculates the -Likelihood of angle theta for a given schedule m_k,n_k - Notes \[L(\theta,\mathbf{h}) = -\sum_{k = 0}^M\log{l_k(\theta|h_k)}\]- Parameters:
- angle (float) – Angle (radians) for calculating the probability of measure a positive event. 
- m_k (list of ints) – number of times the grover operator was applied. 
- n_k (list of ints) – number of total events measured for the specific m_k 
- h_k (list of ints) – number of positive events measured for each m_k 
 
- Returns:
- cost – the aggregation of the individual likelihoods 
- Return type:
- float 
 
 - property index
- creating index property 
 - static likelihood(theta: float, m_k: int, n_k: int, h_k: int) float
- Calculates Likelihood from Suzuki paper. For h_k positive events of n_k total events, this function calculates the probability of this taking into account that the probability of a positive event is given by theta and by m_k The idea is use this function to minimize it for this reason it gives minus Likelihood - Notes \[l_k(\theta|h_k) = \sin^2\left((2m_k+1)\theta\right)^{h_k} \ \cos^2 \left((2m_k+1)\theta\right)^{n_k-h_k}\]- Parameters:
- theta (float) – Angle (radians) for calculating the probability of measure a positive event. 
- m_k (int) – number of times the grover operator was applied. 
- n_k (int) – number of total events measured for the specific m_k 
- h_k (int) – number of positive events measured for each m_k 
 
- Returns:
- Gives the Likelihood p(h_k with m_k amplifications|theta) 
- Return type:
- float 
 
 - static log_likelihood(theta: float, m_k: int, n_k: int, h_k: int) float
- Calculates log of the likelihood from Suzuki paper. - Notes \[\log{l_k(\theta|h_k)} = 2h_k\log\big[\sin\left((2m_k+1) \ \theta\right)\big] +2(n_k-h_k)\log\big[\cos\left((2m_k+1) \ \theta\right)\big]\]- Parameters:
- theta (float) – Angle (radians) for calculating the probability of measure a positive event. 
- m_k (int) – number of times the grover operator was applied. 
- n_k (int) – number of total events measured for the specific m_k 
- h_k (int) – number of positive events measured for each m_k 
 
- Returns:
- Gives the log Likelihood p(h_k with m_k amplifications|theta) 
- Return type:
- float 
 
 - mlae(schedule, optimizer)
- This method executes a complete Maximum Likelihood Algorithm, including executing schedule, defining the correspondent cost function and optimizing it. - Parameters:
- schedule (list of two lists) – the schedule for the algorithm 
- optimizer (optimization routine.) – the optimizer should receive a function of one variable the angle to be optimized. Using lambda functions is the recommended way. 
 
- Returns:
- result (optimizer results) – the type of the result is the type of the result of the optimizer 
- h_k (list) – list with number of positive outcomes from quantum circuit for each pair element of the input schedule 
- cost_function_partial (function) – partial cost function with the m_k, n_k and h_k fixed to the obtained values of the different experiments. 
 
 
 - property oracle
- creating oracle property 
 - run() float
- run method for the class. - Returns:
- list with the estimation of a 
- Return type:
- result 
 - Notes \[a^* = \sin^2(\theta^*) \; where \; \theta^* = \arg \ \min_{\theta} L(\theta,\mathbf{h})\]
 - run_schedule(schedule)
- This method execute the run_step method for each pair of values of a given schedule. - Parameters:
- schedule (list of two lists) – the schedule for the algorithm 
- Returns:
- h_k – list with the h_k result of each pair of the input schedule 
- Return type:
- list 
 
 - run_step(m_k: int, n_k: int) int
- This method executes on step of the MLAE algorithm - Parameters:
- m_k (int) – number of times to apply the self.q_gate to the quantum circuit 
- n_k (int) – number of shots 
 
- Returns:
- h_k (int) – number of positive events 
- routine (QLM Routine object) 
 
 
 - property schedule
- creating schedule property 
 - set_exponential_schedule(n_t: int, n_s: int)
- Creates a scheduler of exponential increasing of m_ks. - Parameters:
- n_t (int) – number of maximum applications of the grover operator 
- n_s (int) – number of shots for each m_k grover applications 
 
 
 - set_linear_schedule(n_t: int, n_s: int)
- Creates a scheduler of linear increasing of m_ks. - Parameters:
- n_t (int) – number of maximum applications of the grover operator 
- n_s (int) – number of shots for each m_k grover applications 
 
 
 - property target
- creating target property 
 
montecarlo_ae
This module contains necessary functions and classes to implement a MonterCarlo Amplitude Estimation. In this case not amplification is used. The probability of the target stat of the oracle is measured.
Author: Gonzalo Ferro Costas & Alberto Manzano Herrero
- class tnbs.BTC_02_AE.QQuantLib.AE.montecarlo_ae.MCAE(oracle: qat.lang.AQASM.QRoutine, target: list, index: list, **kwargs)
- Bases: - object- Class for MonteCarlo Amplitude Estimation (MCAE). - Parameters:
- oracle (QLM gate) – QLM gate with the Oracle for implementing the Grover operator 
- target (list of ints) – python list with the target for the amplitude estimation 
- index (list of ints) – qubits which mark the register to do the amplitude estimation 
- kwars (dictionary) – - dictionary that allows the configuration of the IQAE algorithm: Implemented keys: - qpuQLM solver
- solver for simulating the resulting circuits 
- shotsint
- number of measurements 
- mcz_qlmbool
- for using or not QLM implementation of the multi controlled Z gate 
 
 
 - property index
- creating index property 
 - property oracle
- creating oracle property 
 - run()
- run method for the class. - Returns:
- amplitude estimation parameter 
- Return type:
- self.ae 
 
 - property target
- creating target property 
 
real_quantum_ae
This module contains necessary functions and classes to implement Real Quantum Amplitude Estimation based on the paper:
Manzano, A., Musso, D., Leitao, A. et al. Real Quantum Amplitude Estimation Preprint
Author: Gonzalo Ferro Costas & Alberto Manzano Herrero
- class tnbs.BTC_02_AE.QQuantLib.AE.real_quantum_ae.RQAE(oracle: qat.lang.AQASM.QRoutine, target: list, index: list, **kwargs)
- Bases: - object- Class for Real Quantum Amplitude Estimation (RQAE) algorithm - Parameters:
- oracle (QLM gate) – QLM gate with the Oracle for implementing the Grover operator 
- target (list of ints) – python list with the target for the amplitude estimation 
- index (list of ints) – qubits which mark the register to do the amplitude estimation 
- kwars (dictionary) – - dictionary that allows the configuration of the IQAE algorithm: Implemented keys: - qpuQLM solver
- solver for simulating the resulting circuits 
- qint
- amplification ratio 
- epsilonint
- precision 
- gammafloat
- accuracy 
- mcz_qlmbool
- for using or not QLM implementation of the multi controlled Z gate 
 
 
 - static chebysev_bound(n_samples: int, gamma: float)
- Computes the length of the confidence interval for a given number of samples n_samples and an accuracy gamma. - Parameters:
- n_samples (int) – number of samples 
- gamma (float) – accuracy 
 
- Return type:
- length of the confidence interval 
 
 - static compute_info(ratio: float = 2, epsilon: float = 0.01, gamma: float = 0.05, **kwargs)
- This function computes theoretical values of the IQAE algorithm. - Parameters:
- ratio (float) – amplification ratio/policy 
- epsilon (float) – precision 
- gamma (float) – accuracy 
 
- Returns:
- info – python dictionary with the computed information 
- Return type:
- dict 
 
 - static display_information(ratio: float = 2, epsilon: float = 0.01, gamma: float = 0.05, **kwargs)
- This function displays information of the properties of the method for a given set of parameters - Parameters:
- ratio (float) – amplification ratio/policy 
- epsilon (float) – precision 
- gamma (float) – accuracy 
 
 
 - first_step(shift: float, shots: int, gamma: float)
- This function implements the first step of the RQAE paper. The result is a first estimation of the desired amplitude. - Parameters:
- shift (float) – shift for the first iteration 
- shots (int) – number of measurements 
- gamma (float) – accuracy 
 
- Returns:
- amplitude_min (float) – lower bound for the amplitude to be estimated 
- amplitude_max (float) – upper bound for the amplitude to be estimated 
 
 
 - property index
- creating index property 
 - property oracle
- creating oracle property 
 - rqae(ratio: float = 2, epsilon: float = 0.01, gamma: float = 0.05)
- This function implements the first step of the RQAE paper. The result is an estimation of the desired amplitude with precision epsilon and accuracy gamma. - Parameters:
- ratio (int) – amplification ratio 
- epsilon (int) – precision 
- gamma (float) – accuracy 
 
- Returns:
- amplitude_min (float) – lower bound for the amplitude to be estimated 
- amplitude_max (float) – upper bound for the amplitude to be estimated 
 
 
 - run()
- run method for the class. - Returns:
- amplitude estimation parameter 
- Return type:
- self.ae 
 
 - run_step(shift: float, shots: int, gamma: float, k: int)
- This function implements a step of the RQAE paper. The result is a refined estimation of the desired amplitude. - Parameters:
- shift (float) – shift for the first iteration 
- shots (int) – number of measurements 
- gamma (float) – accuracy 
- k (int) – number of amplifications 
 
- Returns:
- amplitude_min (float) – lower bound for the amplitude to be estimated 
- amplitude_max (float) – upper bound for the amplitude to be estimated 
 
 
 - property shifted_oracle
- creating shifted_oracle property 
 - property target
- creating target property