Quantum Signal Processing and Quantum Singular Value Transformation¶
Copyright (c) 2022 Institute for Quantum Computing, Baidu Inc. All Rights Reserved.
Overview¶
Quantum signal processing (QSP), originally purposed by Guang Hao Low and Issac L. Chuang [1], is an algorithm for simulation of polynomial transformation of scalars using reflection () and rotation () operators. QSP states that, for a specific group of polynomials like the Chebyshev polynomials of the first kind, such simulation can be done by a single qubit. Moreover, using a technique called "linear combination of unitaries", we can use two more ancilla qubit to simulate all complex polynomials and hence approximate any functions that can be expanded by Taylor series.
The idea of QSP is further generalized by paper [2], so that it can simulate polynomial transformation of matrices using high-dimensional rotation operators. Since the polynomial transformation of a diagonalizable matrix is intrinsically the polynomial transformation of its singular values, this algorithm is called the quantum singular value transformation (QSVT).
Moreover, a technique called qubitization is further employed to substantially reduce the complexity of high-dimensional rotation operators, so that one (ancilla) qubit suffices to perform the rotation in larger space.
Based on paper [2], this tutorial will illustrate quantum signal processing, quantum singular value transformation and qubitization, and how to implement these algorithms in PaddleQuantum. But first of all, we need to present the idea of block encoding.
Here are some necessary libraries and packages.
import numpy as np
from numpy.polynomial.polynomial import Polynomial
from numpy.polynomial import Chebyshev
import paddle
# PaddleQuantum related packages
import paddle_quantum
from paddle_quantum.ansatz import Circuit
from paddle_quantum.qinfo import dagger
from paddle_quantum.linalg import unitary_random_with_hermitian_block, density_matrix_random
from paddle_quantum.qsvt import poly_matrix, ScalarQSP, Phi_verification, reflection_based_quantum_signal_processing, QSVT
Block Encoding¶
Block encoding is a method to encode matrices of data. In quantum computing, for a matrix , a unitary is said to be a block encoding of if there exists two orthogonal projectors such that
This could appear more familiar when , where above equation implies
That is, is the left-upper block of matrix .
Encoding and Decoding of Matrices¶
There are numerous ways to find a block encoding unitary of a matrix . For example, when is diagonalizable, can be constructed as
when is a unitary, we can simply select i.e. is a controlled . Despite the fact that we do not hold a mature method for doing this in a quantum system, several studies [3] [4] have been conducted on the quantum implementation of block encodings. In this tutorial, we will assume that there exists an efficient algorithm to block encode a matrix (or at least a Hermitian matrix) in a quantum circuit.
Realization in PaddleQuantum.¶
num_qubits = 3
num_block_qubits = 2
# create a 3-qubit block encoding unitary U of a random 2-qubit Hermitian matrix A
U = unitary_random_with_hermitian_block(num_qubits)
A = U[0: 2 ** num_block_qubits, 0: 2 ** num_block_qubits]
As for the decoding part, the information of matrix can be extracted by applying its block encoding unitary on , where is an -qubit quantum state. Then the output state is
Measuring the first quantum register with outcome gives the desired information.
The success probability of decoding is the probability of measuring , which can be exponentially small as increases. However, this would no longer be a problem if we can control the size of first register. In this tutorial we will assume .
# set density matrix backend
paddle_quantum.set_backend('density_matrix')
# define input state
rho = density_matrix_random(num_block_qubits)
zero_state = paddle.eye(2 ** (num_qubits - num_block_qubits), 1) # create a 0 state (in state_vector form)
input_state = paddle_quantum.State(paddle.kron(zero_state @ dagger(zero_state), rho))
# define auxiliary register
aux_register = list(range(num_qubits - num_block_qubits))
# construct the circuit
cir = Circuit(num_qubits)
cir.oracle(U, list(range(num_qubits)))
cir.collapse(aux_register, desired_result='0', if_print = True) # call Collapse operator
# get output_state and actual output rho
output_state = cir(input_state)
output_rho = output_state.data[0: 2 ** num_block_qubits, 0: 2 ** num_block_qubits]
qubits [0] collapse to the state |0> with probability 0.25072506070137024
We can verify that the output_rho
is the same as .
expect_rho = A @ rho @ dagger(A)
expect_rho /= paddle.trace(expect_rho)
print(f"the difference between input and output is {paddle.norm(paddle.abs(expect_rho - output_rho)).item()}")
the difference between input and output is 0.0
When (and hence ), is a scalar. In this case we can calculate the expectance of in terms of state to receive .
Quantum Signal Processing¶
Theorem 3-4 in paper [2] state that, for a polynomial with that satisfies
- if is even/odd, then is a even/odd polynomial;
- ;
- ;
- if is even, then ,
there exists a polynomial and a vector of angles such that ,
That is, we can use rotation () gates and reflection () gates to receive a block encoding unitary of . After block decoding we have completed the simulation of .
Here are some polynomials that satisfy above requirements.
# simple odd polynomials
# P = Polynomial([0, 0, 0, 0, 0, 1])
# P = Polynomial([0, 0.5, 0, 0, 0, 0.5])
P = Polynomial([0, 1 / 3, 0, 0, 0, 1 / 3, 0, 1 / 3])
# Chebyshev polynomials of first kind with degree 10
# P = Chebyshev([0 for _ in range(10)] + [1]).convert(kind=Polynomial)
# Chebyshev polynomials of first kind with degree 11
# P = Chebyshev([0 for _ in range(11)] + [1]).convert(kind=Polynomial)
Note that such structure of naturally fits the family of Chebyshev polynomials. Indeed, for a Chebyshev polynomials of first kind with order , its corresponding is (if is even) or (if is odd).
PaddleQuantum has a built-in class ScalarQSP
for quantum signal processing and a function Phi_verification
to verify whether is a block encoding of .
qsp = ScalarQSP(P)
print(qsp.Phi)
assert Phi_verification(qsp.Phi.numpy(), P)
Tensor(shape=[8], dtype=float64, place=Place(cpu), stop_gradient=True, [ 3.14159265, 1.39460474, -0.44313783, 1.09757975, -1.09757975, 0.44313783, -1.39460474, 3.14159265])
We can also verify the correctness of from the corresponding circuit.
x = np.random.rand() * 2 - 1 # random x in [-1, 1]
cir = qsp.block_encoding(x)
print(cir)
print(f"the error of simulating P(x) is {np.abs(cir.unitary_matrix()[0, 0].item() - P(x))}")
--Rz(-6.28)----Rx(-1.74)----Rz(2.789)----Rx(-1.74)----Rz(-0.88)----Rx(-1.74)----Rz(2.195)----Rx(-1.74)----Rz(-2.19)----Rx(-1.74)----Rz(0.886)----Rx(-1.74)----Rz(-2.78)----Rx(-1.74)----Rz(-6.28)-- the error of simulating P(x) is 6.378721218542117e-08
Variant of QSP¶
Moreover, the number of parameters can be further simplified to by finding a vector of angles such that
where
In PaddleQuantum, we can compute such by function reflection_based_quantum_signal_processing
.
print(reflection_based_quantum_signal_processing(P))
[15.70796327 -0.17619159 -2.01393416 -0.47321657 -2.66837608 -1.12765849 -2.96540107]
Note: is a block encoding of .
Quantum Singular Value Transformation¶
To start with, we need to define the singular value transformation.
Let be a matrix. Suppose there exists a unitary and two orthogonal projectors and such that . Then there exists two orthonormal basis and extended by eigenstates of and respectively, such that
where
is a singular value of . The singular value transformation (SVT) of for a function is defined to be
Decomposition and QSP¶
Paper [2] proves that under appropriate choices of basis, and can be decomposed by subspaces of such that
We can decompose in a similar manner.
Let's connect above decomposition with quantum signal processing. Suppose function is a polynomial with degree that satisfies the requirements in the last section. Then there exists such that is a block encoding of . Moreover, satisfies and . Now Define as
Then
and hence
If we can realize or by a quantum algorithm, then such transformation is the quantum singular value transformation of .
QSVT and Block Encoding¶
Suppose is a block encoding unitary of with respect to orthogonal projectors and . Then and . Therefore, is a block encoding of .
When , . Denote . We can find that
That is, when , QSVT maps a block encoding of to a block encoding of .
PaddleQuantum has a built-in class QSVT
for quantum singular value transformation when . We can call its member function QSVT.block_encoding_matrix()
to verify the correctness of above theories, from entries to entries and from eigenvalues to eigenvalues.
qsvt = QSVT(P, U, num_block_qubits)
# find P(A) and its expected eigenvalues, note that they are computed in different ways
expect_PX = poly_matrix(P, A).numpy()
expect_eig_result = np.sort(list(map(lambda x: P(x), np.linalg.eigvals(A.numpy()))))
# Calculate U_\Phi and extract eigenvalues of block encoded matrix
U_Phi = qsvt.block_encoding_matrix().numpy()
actual_PX = U_Phi[0:2 ** num_block_qubits, 0:2 ** num_block_qubits]
actual_eig_result = np.sort(np.linalg.eigvals(actual_PX))
print("error of simulating P(X)")
print(" maximum absolute, ", np.amax(abs(expect_PX - actual_PX)))
print(" percentage, ", np.linalg.norm(expect_PX - actual_PX) / np.linalg.norm(expect_PX))
error of simulating P(X) maximum absolute, 1.6997650495437753e-07 percentage, 3.4201093195057237e-07
print("error of eigenvalues of simulating P(X)")
print(" maximum absolute, ", np.amax(abs(expect_eig_result - actual_eig_result)))
print(" percentage, ", np.linalg.norm(expect_eig_result - actual_eig_result) / np.linalg.norm(expect_eig_result))
error of eigenvalues of simulating P(X) maximum absolute, 2.3308557146996258e-07 percentage, 2.2962108806419593e-07
Qubitization: Quantum Realization of ¶
One question is, how to realize by a quantum circuit? Note that we assume is implementable. Then the remaining problem is the realization of unitaries and , which are essentially operators performed in projection spaces of and respectively.
To achieve this, Lemma 10 in paper [1] projects the image space of projector into an ancilla qubit and perform the rotation operation on that qubit. Then by entanglement such rotation is simultaneously applied to the main register.
The results are summarized to the following circuit:
When , PaddleQuantum can create a circuit in Figure 2 by member function QSVT.block_encoding_circuit
. We can show that the constructed quantum circuit can simulate , by comparison of the output state and for .
# set state vector backend
paddle_quantum.set_backend('state_vector')
# arbitrary state such that last qubit is in state |0\rangle
ket_0 = [1.0, 0.0]
psi = np.array([np.random.rand() + 1j * np.random.rand() for _ in range(2 ** num_qubits)])
psi = np.kron(psi / np.linalg.norm(psi), ket_0)
expect_state = np.kron(U_Phi, np.eye(2)) @ psi
cir = qsvt.block_encoding_circuit()
actual_state = cir(paddle_quantum.State(psi, dtype='complex64')).data.numpy()
print(f"the different between expected and actual state is {np.linalg.norm(expect_state - actual_state)}")
the different between expected and actual state is 2.383485634049416e-07
Note: if we add Hadamard gates at the beginning and the end of the ancilla register in Figure 2, then the quantum circuit turns to simulate the real part of the polynomial. Combined with the technique of linear combination of unitaries, theocratically we can simulate all complex polynomials with norm smaller than 1 using quantum circuit.
Application: Amplitude Amplification¶
Let be a -qubit quantum state such that . Amplitude amplification is a quantum algorithm that can increase the amplitude of to approximately 1. This algorithm can be viewed in a different prospective. Let be the unitary and be two orthogonal projectors such that and , implying . Therefore, amplitude amplification is essentially a singular value transformation of . In this section we will show how to use QSVT to do fixed-point (i.e. for some ) amplitude amplification from Theorem 28 of paper [2].
Suppose and for some . Then and , where is the left-upper block of .
Observe that . Therefore, we aim to find a quantum circuit with unitary such that , implying . By , the absolute values of singular values of are . Moreover, choose , where is the Chebyshev polynomial of the first kind of order . Then
implies that the QSVT of in terms of polynomial is a block encoding of , as required.
We set so that , and is a randomly chosen unitary. We want to show that the left-upper block of is amplified by 2 after QSVT.
P = -1 * Chebyshev([0 for _ in range(3)] + [1]).convert(kind=Polynomial)
U = unitary_random_with_hermitian_block(num_qubits, is_unitary=True)
A = U[0:2 ** num_block_qubits, 0:2 ** num_block_qubits].numpy()
amplifier = QSVT(P, U, num_block_qubits)
U_Phi = amplifier.block_encoding_unitary()
B = U_Phi[0:2 ** num_block_qubits, 0:2 ** num_block_qubits].numpy()
print(f"the accuracy of quantum singular value transformation is {np.linalg.norm(B - 2 * A)}")
the accuracy of quantum singular value transformation is 3.029211939065135e-07
As shown above, we successfully amplify by 2.
References¶
[1] Low, Guang Hao, and Isaac L. Chuang. "Hamiltonian simulation by qubitization." Quantum 3 (2019): 163.
[2] Gilyén, András, et al. "Quantum singular value transformation and beyond: exponential improvements for quantum matrix arithmetics." Proceedings of the 51st Annual ACM SIGACT Symposium on Theory of Computing. 2019.
[3] Camps, Daan, et al. "Explicit Quantum Circuits for Block Encodings of Certain Sparse Matrices." arXiv preprint arXiv:2203.10236 (2022).
[4] Clader, B. David, et al. "Quantum Resources Required to Block-Encode a Matrix of Classical Data." arXiv preprint arXiv:2206.03505 (2022).