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 (RXRX) and rotation (RZRZ) 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∈C2n×2nA∈C2n×2n, a unitary U∈C2m×2mU∈C2m×2m is said to be a block encoding of AA if there exists two orthogonal projectors W,V∈C2m×2mW,V∈C2m×2m such that
WUV=[A000]=A⊕0I⊗(m−n).WUV=[A000]=A⊕0I⊗(m−n).(1)This could appear more familiar when W=V=|0⟩⊗(m−n)⟨0|⊗(m−n)⊗I⊗nW=V=|0⟩⊗(m−n)⟨0|⊗(m−n)⊗I⊗n, where above equation implies
U=[A.........].U=[A.........].(2)That is, AA is the left-upper block of matrix UU.
Encoding and Decoding of Matrices¶
There are numerous ways to find a block encoding unitary UU of a matrix AA. For example, when AA is diagonalizable, UU can be constructed as
U=[Ai√I⊗n−A2i√I⊗n−A2A];U=[Ai√I⊗n−A2i√I⊗n−A2A];(3)when AA is a unitary, we can simply select U=A⊕I⊗(m−n)U=A⊕I⊗(m−n) i.e. UU is a controlled AA. 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 A∈C2n×2nA∈C2n×2n can be extracted by applying its block encoding unitary U∈C2m×2mU∈C2m×2m on |0⟩⊗(m−n)⟨0|⊗(m−n)⊗ρ|0⟩⊗(m−n)⟨0|⊗(m−n)⊗ρ, where ρρ is an nn-qubit quantum state. Then the output state is
U(|0⟩⊗(m−n)⟨0|⊗(m−n)⊗ρ)U†=[AρA†.........]=|0⟩⊗(m−n)⟨0|⊗(m−n)⊗AρA†+(I⊗(m−n)−|0⟩⊗(m−n)⟨0|⊗(m−n))⊗...U(|0⟩⊗(m−n)⟨0|⊗(m−n)⊗ρ)U†=[AρA†.........]=|0⟩⊗(m−n)⟨0|⊗(m−n)⊗AρA†+(I⊗(m−n)−|0⟩⊗(m−n)⟨0|⊗(m−n))⊗...(4)Measuring the first quantum register with outcome 0⊗(m−n)0⊗(m−n) gives the desired information.
The success probability of decoding is the probability of measuring 0⊗(m−n)0⊗(m−n), which can be exponentially small as m−nm−n increases. However, this would no longer be a problem if we can control the size of first register. In this tutorial we will assume m=n+1m=n+1.
# 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 AρA†Tr(AρA†)AρA†Tr(AρA†).
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 m=1m=1 (and hence n=0n=0), AA is a scalar. In this case we can calculate the expectance of UU in terms of state |0⟩|0⟩ to receive AA.
Quantum Signal Processing¶
Theorem 3-4 in paper [2] state that, for a polynomial P∈C[x]P∈C[x] with deg(P)=kdeg(P)=k that satisfies
- if kk is even/odd, then PP is a even/odd polynomial;
- |P(x)|≤1,∀x∈[−1,1]|P(x)|≤1,∀x∈[−1,1];
- |P(x)|≥1,∀x∈(−∞,−1]∪[1,∞)|P(x)|≥1,∀x∈(−∞,−1]∪[1,∞);
- if kk is even, then P(ix)P∗(ix)≥1P(ix)P∗(ix)≥1,
there exists a polynomial Q∈C[x]Q∈C[x] and a vector of angles Φ=(ϕ0,...,ϕk)∈Rk+1Φ=(ϕ0,...,ϕk)∈Rk+1 such that ∀x∈[−1,1]∀x∈[−1,1],
WΦ(x):=RZ(−2ϕ0)k∏j=1RX(arccos(x))RZ(−2ϕj)=[P(x)iQ(x)√1−x2iQ∗(x)√1−x2P∗(x)].WΦ(x):=RZ(−2ϕ0)k∏j=1RX(arccos(x))RZ(−2ϕj)=[P(x)iQ(x)√1−x2iQ∗(x)√1−x2P∗(x)].(5)That is, we can use k+1k+1 rotation (RZRZ) gates and kk reflection (RXRX) gates to receive a block encoding unitary WΦ(x)WΦ(x) of P(x)P(x). After block decoding we have completed the simulation of P(x)P(x).
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 WΦWΦ naturally fits the family of Chebyshev polynomials. Indeed, for a Chebyshev polynomials of first kind with order kk, its corresponding ΦΦ is (0,π,...,π)(0,π,...,π) (if kk is even) or (π,...,π)(π,...,π) (if kk is odd).
PaddleQuantum has a built-in class ScalarQSP
for quantum signal processing and a function Phi_verification
to verify whether WΦWΦ is a block encoding of PP.
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 kk by finding a vector of angles Φ∈RkΦ∈Rk such that
RΦ(x):=k∏i=1eiϕiσzR(x)=[P(x).........],RΦ(x):=k∏i=1eiϕiσzR(x)=[P(x).........],(6)where
R(x)=[x√1−x2√1−x2−x].R(x)=[x√1−x2√1−x2−x].(7)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: R−Φ(x)R−Φ(x) is a block encoding of P∗(x)P∗(x).
Quantum Singular Value Transformation¶
To start with, we need to define the singular value transformation.
Let X∈C2m×2mX∈C2m×2m be a matrix. Suppose there exists a unitary UU and two orthogonal projectors WW and VV such that X=WUVX=WUV. Then there exists two orthonormal basis {|ωj⟩}2mj=1{|ωj⟩}2mj=1 and {|νj⟩}2mj=1{|νj⟩}2mj=1 extended by eigenstates of WW and VV respectively, such that
X=2m∑j=1λj|ωj⟩⟨νj|,X=2m∑j=1λj|ωj⟩⟨νj|,(8)where
λj={⟨ωj|U|νj⟩ if j≤rank(X)0 otherwise λj={⟨ωj|U|νj⟩ if j≤rank(X)0 otherwise (9)is a singular value of XX. The singular value transformation (SVT) of XX for a function ff is defined to be
f(SV)(X):=2m∑j=1f(λj){|νj⟩⟨νj| if f is even|ωj⟩⟨νj| if f is odd.f(SV)(X):=2m∑j=1f(λj){|νj⟩⟨νj| if f is even|ωj⟩⟨νj| if f is odd.(10)Decomposition and QSP¶
Paper [2] proves that under appropriate choices of basis, U,VU,V and WW can be decomposed by subspaces of λjλj such that
U=⨁j:λj∈[0,1)R(λj)⊕⨁j:λj=1[1]⊕⨁j:V|νj⟩=|νj⟩,WU|νj⟩=0[1]⊕⨁j:W|ωj⟩=|ωj⟩,VU†|ωj⟩=0[1]⊕⨁j:V|νj⟩=W|ωj⟩=0[...],V=⨁j:λj∈[0,1)|0⟩⟨0|⊕⨁j:λj=1[1]⊕⨁j:V|νj⟩=|νj⟩,WU|νj⟩=0[1]⊕⨁j:W|ωj⟩=|ωj⟩,VU†|ωj⟩=0[0]⊕⨁j:V|νj⟩=W|ωj⟩=0[0] andW=⨁j:λj∈[0,1)|0⟩⟨0|⊕⨁j:λj=1[1]⊕⨁j:V|νj⟩=|νj⟩,WU|νj⟩=0[0]⊕⨁j:W|ωj⟩=|ωj⟩,VU†|ωj⟩=0[1]⊕⨁j:V|νj⟩=W|ωj⟩=0[0].U=⨁j:λj∈[0,1)R(λj)⊕⨁j:λj=1[1]⊕⨁j:V|νj⟩=|νj⟩,WU|νj⟩=0[1]⊕⨁j:W|ωj⟩=|ωj⟩,VU†|ωj⟩=0[1]⊕⨁j:V|νj⟩=W|ωj⟩=0[...],V=⨁j:λj∈[0,1)|0⟩⟨0|⊕⨁j:λj=1[1]⊕⨁j:V|νj⟩=|νj⟩,WU|νj⟩=0[1]⊕⨁j:W|ωj⟩=|ωj⟩,VU†|ωj⟩=0[0]⊕⨁j:V|νj⟩=W|ωj⟩=0[0] andW=⨁j:λj∈[0,1)|0⟩⟨0|⊕⨁j:λj=1[1]⊕⨁j:V|νj⟩=|νj⟩,WU|νj⟩=0[0]⊕⨁j:W|ωj⟩=|ωj⟩,VU†|ωj⟩=0[1]⊕⨁j:V|νj⟩=W|ωj⟩=0[0].(11)We can decompose f(SV)f(SV) in a similar manner.
f(SV)(X)=∑j:λj∈[0,1)f(λj)...+∑j:λj=1f(1)...+∑j:λj=0,V|νj⟩≠0f(0)...+∑j:λj=0,W|ωj⟩≠0f(0)...+∑j:λj=0,otherwisef(0)...,f(SV)(X)=∑j:λj∈[0,1)f(λj)...+∑j:λj=1f(1)...+∑j:λj=0,V|νj⟩≠0f(0)...+∑j:λj=0,W|ωj⟩≠0f(0)...+∑j:λj=0,otherwisef(0)...,(12)Let's connect above decomposition with quantum signal processing. Suppose function ff is a polynomial P∈C[x]P∈C[x] with degree kk that satisfies the requirements in the last section. Then there exists Φ∈RkΦ∈Rk such that RΦRΦ is a block encoding of PP. Moreover, PP satisfies P(1)=ei∑kj=1ϕjP(1)=ei∑kj=1ϕj and P(0)=ei∑kj=1(−1)j+1ϕjP(0)=ei∑kj=1(−1)j+1ϕj. Now Define UΦUΦ as
UΦ:={∏k/2j=1eiϕ2j−1(2V−I)U†eiϕ2j(2W−I)U if k is eveneiϕ1(2W−I)U∏(k−1)/2j=1eiϕ2j(2V−I)U†eiϕ2j+1(2W−I)U if k is odd.UΦ:=⎧⎨⎩∏k/2j=1eiϕ2j−1(2V−I)U†eiϕ2j(2W−I)U if k is eveneiϕ1(2W−I)U∏(k−1)/2j=1eiϕ2j(2V−I)U†eiϕ2j+1(2W−I)U if k is odd.(13)Then
UΦ={⨁RΦ(λj)⊕⨁[ei∑kj=1ϕj]⊕⨁[ei∑kj=1(−1)j+1ϕj]⊕⨁[ei∑kj=1(−1)j+1ϕj]⊕⨁[...] if k is even⨁RΦ(λj)⊕⨁[ei∑kj=1ϕj]⊕⨁[ei∑kj=1(−1)j+1ϕj]⊕⨁[ei∑kj=1(−1)j+1ϕj]⊕⨁[...] if k is oddUΦ={⨁RΦ(λj)⊕⨁[ei∑kj=1ϕj]⊕⨁[ei∑kj=1(−1)j+1ϕj]⊕⨁[ei∑kj=1(−1)j+1ϕj]⊕⨁[...] if k is even⨁RΦ(λj)⊕⨁[ei∑kj=1ϕj]⊕⨁[ei∑kj=1(−1)j+1ϕj]⊕⨁[ei∑kj=1(−1)j+1ϕj]⊕⨁[...] if k is odd(14)and hence
P(SV)(X)={⨁[P(λj)000]⊕⨁[P(1)]⊕⨁[P(0)]⊕⨁[0]⊕⨁[0]=VUΦV if k is even⨁[P(λj)000]⊕⨁[P(1)]⊕⨁[0]⊕⨁[0]⊕⨁[0]=WUΦV if k is odd.P(SV)(X)=⎧⎪ ⎪ ⎪ ⎪⎨⎪ ⎪ ⎪ ⎪⎩⨁[P(λj)000]⊕⨁[P(1)]⊕⨁[P(0)]⊕⨁[0]⊕⨁[0]=VUΦV if k is even⨁[P(λj)000]⊕⨁[P(1)]⊕⨁[0]⊕⨁[0]⊕⨁[0]=WUΦV if k is odd.(15)If we can realize VUΦVVUΦV or WUΦVWUΦV by a quantum algorithm, then such transformation is the quantum singular value transformation of XX.
QSVT and Block Encoding¶
Suppose UU is a block encoding unitary of AA with respect to orthogonal projectors WW and VV. Then X=A⊕0I⊗(m−n)X=A⊕0I⊗(m−n) and P(SV)(X)=P(SV)(A)⊕0I⊗(m−n)P(SV)(X)=P(SV)(A)⊕0I⊗(m−n). Therefore, UΦUΦ is a block encoding of P(SV)(A)P(SV)(A).
When W=VW=V, {|ωj⟩}2mj=1={|νj⟩}2mj=1{|ωj⟩}2mj=1={|νj⟩}2mj=1. Denote P(x):=∑ki=0cixiP(x):=∑ki=0cixi. We can find that
P(SV)(X)=2m∑j=1P(λj)|νj⟩⟨νj|=P(X)=P(A)⊕0I⊗(m−n).P(SV)(X)=2m∑j=1P(λj)|νj⟩⟨νj|=P(X)=P(A)⊕0I⊗(m−n).(16)That is, when W=VW=V, QSVT maps a block encoding of AA to a block encoding of P(A)P(A).
PaddleQuantum has a built-in class QSVT
for quantum singular value transformation when W=V=|0⟩⊗(m−n)⟨0|⊗(m−n)⊗I⊗nW=V=|0⟩⊗(m−n)⟨0|⊗(m−n)⊗I⊗n. 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 UΦUΦ¶
One question is, how to realize UΦUΦ by a quantum circuit? Note that we assume UU is implementable. Then the remaining problem is the realization of unitaries eiϕ(2W−I)eiϕ(2W−I) and eiϕ(2V−I)eiϕ(2V−I), which are essentially RZ(−2ϕ)RZ(−2ϕ) operators performed in projection spaces of WW and VV 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 W=V=|0⟩⊗(m−n)⟨0|⊗(m−n)⊗I⊗nW=V=|0⟩⊗(m−n)⟨0|⊗(m−n)⊗I⊗n, 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 UΦUΦ, by comparison of the output state and (UΦ⊗I)|ψ⟩|0⟩(UΦ⊗I)|ψ⟩|0⟩ for |ψ⟩∈C2m|ψ⟩∈C2m.
# 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 mm-qubit quantum state such that |ψ⟩=sin(θ)|ψgood⟩+cos(θ)|ψbad⟩|ψ⟩=sin(θ)|ψgood⟩+cos(θ)|ψbad⟩. Amplitude amplification is a quantum algorithm that can increase the amplitude of |ψgood⟩ to approximately 1. This algorithm can be viewed in a different prospective. Let U be the unitary and W,V be two orthogonal projectors such that |ψ⟩=UV|0⟩⊗m and sin(θ)|ψgood⟩=W|ψ⟩, implying sin(θ)|ψgood⟩=WUV|0⟩⊗m. Therefore, amplitude amplification is essentially a singular value transformation of X=WUV. In this section we will show how to use QSVT to do fixed-point (i.e. θ=π2k for some k∈Z) amplitude amplification from Theorem 28 of paper [2].
Suppose sin(θ)|ψbad⟩=(|ψgood⟩⊗(m−n)⟨0|⊗(m−n)⊗I⊗n)|ψ⟩ and θ=π2k for some k∈Z. Then W=V=|0⟩⊗(m−n)⟨0|⊗(m−n)⊗I⊗n and X=A⊕0I⊗(m−n), where A is the left-upper block of U.
Observe that X|ψgood⟩⊗m=sin(π2k)|ψ⟩. Therefore, we aim to find a quantum circuit with unitary UΦ such that WUΦV=1sin(π2k)X, implying WUΦV|ψgood⟩⊗m=|0⟩. By X=WUV, the absolute values of singular values of X are sin(π2k). Moreover, choose P(x)=(−1)kTk(x), where Tk is the Chebyshev polynomial of the first kind of order k. Then
P(sin(π2k))=(−1)kTk(sin(π2k))=(−1)kcos(k−12π)=1,implies that the QSVT of X in terms of polynomial P is a block encoding of B:=1sin(π2k)A, as required.
We set k=3 so that sin(π2k)=12, and U is a randomly chosen unitary. We want to show that the left-upper block of U 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 X 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).