Generalized Quantum Signal Processing (GQSP)#

GQSP(anc: QuantumBool, *qargs: QuantumVariable, unitary: Callable[[...], None], p: ArrayLike | None = None, angles: Tuple[ArrayLike, ArrayLike, ArrayLike] | None = None, k: int = 0, kwargs: Dict[str, Any] = {}) None[source]#

Performs Generalized Quantum Signal Processing.

Generalized Quantum Signal Processing was introduced by Motlagh and Wiebe. Its equivalence to the non-linear Fourier transform (NLFT) was subsequently established by Laneve. By adopting the conventions from Laneve, we treat Generalized Quantum Signal Processing (GQSP) as the overarching framework for single-qubit signal processing. In this setting, previously known QSP variants—such as those constrained to \(X\) or \(Y\) rotations—emerge naturally as special cases where the underlying complex sequences of the NLFT are restricted to purely imaginary or real values, respectively.

GQSP is described as follows:

  • Given a unitary \(U\),

  • and a complex degree \(d\) polynomial \(p(z)\in\mathbb C[z]\) such that \(|p(e^{ix})|^2\leq 1\) for all \(x\in\mathbb R\),

this algorithm implements the unitary

\[\begin{split}\begin{pmatrix} p(U) & * \\ * & * \end{pmatrix} &=R(\theta_0,\phi_0,\lambda)\left(\prod_{j=1}^{d}AR(\theta_j,\phi_j,0)\right)\\ &=R(\theta_0,\phi_0,\lambda)AR(\theta_1,\phi_1,0)A\dotsc AR(\theta_d,\phi_d,0)\end{split}\]

where

\[R(\theta,\phi,\lambda) = e^{i\lambda Z}e^{i\phi X}e^{i\theta Z} \in SU(2)\]

and the angles \(\theta,\phi\in\mathbb R^{d+1}\), \(\lambda\in\mathbb R\) are calculated from the polynomial \(p\), \(A=\begin{pmatrix}U & 0\\ 0 & I\end{pmatrix}\) is the signal operator.

Parameters:
ancQuantumBool

Auxiliary variable in state \(\ket{0}\) for applying the GQSP protocol. Must be measured in state \(\ket{0}\) for the GQSP protocol to be successful.

*qargsQuantumVariable

QuantumVariables serving as operands for the unitary.

unitaryCallable

A function applying a unitary to the variables *qargs. Typically, \(U=e^{iH}\) for a Hermitian operator \(H\) and GQSP applies a function of \(H\).

pArrayLike, optional

1-D array containing the polynomial coefficients, ordered from lowest order term to highest. Either the polynomial p or angles must be specified.

anglestuple(ArrayLike, ArrayLike, ArrayLike), optional

A tuple of angles \((\theta,\phi,\lambda)\) where \(\theta,\phi\in\mathbb R^{d+1}\) are 1-D arrays and \(\lambda\in\mathbb R\) is a scalar.

kint

If specified, the Laurent polynomial \(\tilde p(x)=x^{-k}p(x)\) is applied. The default is 0.

kwargsdict

A dictionary of keyword arguments to pass to unitary. The default is {}.

Notes

  • The polynomial \(p\) is rescaled automatically to satisfy \(|p(e^{ix})|^2\leq 1\) for all \(x\in\mathbb R\).

Examples

Applying a transformation in Fourier basis

We apply the operator

\[\cos(H) = \frac{e^{iH}+e^{-iH}}{2}\]

for some Hermitian operator \(H\) to the input state \(\ket{\psi}=\ket{0}\).

First, we define an operator \(H\) and the unitary performing the Hamiltonian evolution \(e^{iH}\). (In this case, Trotterization will perform Hamiltonian evolution exactly since the individual terms commute.)

from qrisp import *
from qrisp.gqsp import *
from qrisp.operators import X,Y,Z
import jax.numpy as jnp

H = Z(0)*Z(1) + X(0)*X(1)

def U(operand):
    H.trotterization(forward_evolution=False)(operand)

Next, we define the operand_prep function that prepares a QuantumVariable is state \(\ket{\psi}=\ket{0}\).

def operand_prep():
    operand = QuantumVariable(2)
    return operand

The transformation \(\cos(H)\) is achieved by applying \(\tilde p(x)=0.5x^{-1} + 0.5x^1\) to the unitary \(e^{iH}\). This corresponds to the polynomial \(p(x)=0.5+0.5x^2\) (i.e., p=[0.5,0,0.5]) and k=1.

Finally, we apply QSP within a Repeat-Until-Success protocol.

@RUS
def inner():

    p = jnp.array([0.5,0,0.5])

    operand = operand_prep()
    anc = QuantumBool()
    GQSP(anc, operand, unitary=U, p=p, k=1)

    success_bool = measure(anc) == 0
    reset(anc)
    anc.delete()
    return success_bool, operand

@terminal_sampling
def main(): 

    qv = inner()
    return qv

and simulate

>>> main()
{3: 0.85471756539818, 0: 0.14528243460182003}

Let’s compare to the classically calculated result:

>>> A = H.to_array()
>>> from scipy.linalg import cosm
>>> print(cosm(A))
[[ 0.29192658+0.j  0.        +0.j  0.        +0.j -0.70807342+0.j]
[ 0.        +0.j  0.29192658+0.j  0.70807342+0.j  0.        +0.j]
[ 0.        +0.j  0.70807342+0.j  0.29192658+0.j  0.        +0.j]
[-0.70807342+0.j  0.        +0.j  0.        +0.j  0.29192658+0.j]]

That is, starting in state \(\ket{\psi}=\ket{0}=(1,0,0,0)\), we obtain

>>> result = cosm(A)@(np.array([1,0,0,0]).transpose())
>>> result = result/np.linalg.norm(result) # normalise
>>> result = result**2 # compute measurement probabilities
>>> print(result)
[0.1452825+0.j 0.       +0.j 0.       +0.j 0.8547175-0.j]

which are exactly the probabilities we observed in the quantum simulation.

Note

While GQSP allows you to apply arbitrary polynomials to operators, applying abitrary polynomials to BlockEncodings requires an additional step. This is because raising the operator

\[\begin{split}U = \begin{pmatrix} \frac{A}{\alpha} & * \\ * & * \end{pmatrix}\end{split}\]

to a given power \(k\) does not necessarily give you

\[\begin{split}\tilde{U} = \begin{pmatrix} \left(\frac{A}{\alpha}\right)^k & * \\ * & * \end{pmatrix}\end{split}\]

In order to still apply polynomials also to them, we need to call the qubitization method and transform the polynomial into Chebychev basis. More to that in the GQSP Tutorials.

Algorithms & Applications#

GQET(H, p[, kind, rescale])

Returns a BlockEncoding representing a polynomial transformation of the operator via Generalized Quantum Eigenvalue Transform.

QET(H, p[, kind, rescale])

Returns a BlockEncoding representing a polynomial transformation of the operator via Quantum Eigenvalue Transform.

hamiltonian_simulation(H[, t, N])

Returns a BlockEncoding approximating Hamiltonian simulation of the operator.

inversion(A, eps, kappa)

Quantum Linear System Solver via Quantum Eigenvalue Transformation (QET).

convolve(qarg, weights)

Performs cyclic convolution of a quantum state with a filter.

fourier_series_loader(qarg[, signal, ...])

Performs quantum state preparation for the first \(k\) Fourier modes.

Utilities#

poly2cheb(poly)

Convert a polynomial from monomial to Chebyshev basis.

cheb2poly(cheb)

Convert a polynomial from Chebyshev to monomial basis.

gqsp_angles(p)

Computes the GQSP angles for a given polynomial.