quil-lang / quilc

The optimizing Quil compiler.
Apache License 2.0
454 stars 73 forks source link

seeking to simplify `find-diagonalizer-in-e-basis` #841

Open stylewarning opened 2 years ago

stylewarning commented 2 years ago

EDIT This issue was previously about using generalized Schur to re-implement F-D-I-E-B, but it turned out to not be e a working approach. As such, I'm changing this issue into a discussion about how we might approach doing so.

I'm mostly interpreting some things @kilimanjaro told me (so credit goes to him), though errors I may have below are my own.

find-diagonalizer-in-e-basis aims to diagonalize a symmetric unitary matrix gammag = u u^T in terms of an orthogonal matrix of eigenvectors. Normally, all that the spectral theorem gives you (in this context) is that you can do this in terms of a unitary matrix of eigenvectors. In this case, however, the real and imaginary parts of gammag commute, so it's sufficient to simultaneously diagonalize them. They are real, symmetric matrices, their eigenvectors will be real and will give an orthogonal matrix

To simultaneously diagonalize a pair of commuting matrices, it's not quite enough to just compute eigenvectors and eigenvalues of one and then hope that this works for the other (consider that the identity matrix, which commutes with anything, but we need to write its eigenspaces as spanned by eigenvectors of some other matrix.) The current QUILC approach is to try to diagonalize a linear combination of the real and imaginary parts. Since there could be some relationships between them that we don't know, the current code tries to pick a random combination, and just repeats until we get something that works

This problem has a standard solution that linear algebra libraries implement, called the generalized Schur decomposition. LAPACK documents it here

A Python implementation (which uses NumPy's qz) supplied by @kilimanjaro is here:

def orthogonal_decomposition(U):
    """
    Given unitary U, decompose UU^T into O @ np.diag(d) @ O.T where O is special orthogonal.
    """
    # Generalized Schur decomposition writes
    #    Re[U] = L D_r R^T,  Im[U] = L D_i R^T
    # In general, D_r, D_i are upper triangular, but for unitary U they end up diagonal
    diag_r, diag_i, left, right = scipy.linalg.qz(np.real(U), np.imag(U))
    diag = np.diagonal(diag_r) + np.diagonal(diag_i) * 1j
    if np.linalg.det(left) < 0:
        left[:,0] *= -1
    return left, diag*diag
genos commented 2 years ago

After my myriad false starts I hesitate to be even cautiously optimistic, but: using the Schur decomposition of $U U^T$ to find a diagonalizing $X$ and then normalizing by its eigenvalues to ensure $\det(X) = 1$ might do a thing.

import numpy as np
from rich import print
from rich.progress import track
import scipy.linalg as la
from scipy.stats import unitary_group
import typer

def decomp(u: np.ndarray) -> np.ndarray:
    _, z = la.schur(u @ u.T)
    return z / la.eigvals(z)

def _test(u: np.ndarray):
    x = decomp(u)
    uut = u @ u.T
    m = x.T @ uut @ x
    _, n = u.shape
    np.testing.assert_allclose(la.det(x), 1, err_msg="det(X) ≠ 1")
    np.testing.assert_equal(
        np.any(np.not_equal(x, u)), True, err_msg="Just returned X = U"
    )
    assert np.abs(m - np.diag(np.diag(m))).max() < 1e-10, "X^T U U^T X not diagonal"

def test_specifics():
    for unitary in track(
        [
            np.array([[1, 0, 0, 0], [0, 0, 0, 1], [0, 0, 1, 0], [0, 1, 0, 0]]),
            np.diag([1, 1, 1, 1j]),
            np.kron(np.eye(2), 1 / np.sqrt(2) * np.array([[1, -1j], [-1j, 1]])),
            np.array(
                [
                    [-0.965 + 0.133j, 0.093 + 0.000j, -0.000 + 0.000j, -0.000 + 0.206j],
                    [0.093 + 0.000j, 0.965 + 0.133j, 0.000 - 0.206j, -0.000 + 0.000j],
                    [-0.000 + 0.000j, 0.000 - 0.206j, 0.965 - 0.133j, -0.093 + 0.000j],
                    [
                        -0.000 + 0.206j,
                        -0.000 + 0.000j,
                        -0.093 + 0.000j,
                        -0.965 - 0.133j,
                    ],
                ]
            ),
        ],
        description="[cyan]Testing specific examples…[/cyan]",
    ):
        _test(unitary)

def test_su4(iterations: int, seed: int):
    rng = unitary_group(dim=4, seed=seed)
    for _ in track(range(iterations), description="[cyan]Testing SU(4)…[/cyan]"):
        _test(rng.rvs())

def test_su2_x_su2(iterations: int, seed: int):
    rng = unitary_group(dim=2, seed=seed)
    for _ in track(range(iterations), description="[cyan]Testing SU(2)⊗SU(2)…[/cyan]"):
        _test(np.kron(rng.rvs(), rng.rvs()))

def main(
    iterations: int = typer.Option(100_000, help="Number of tests to run"),
    su4_seed: int = typer.Option(96692877, help="Seed for SU(4) tests (random.org)"),
    su2_seed: int = typer.Option(29676226, help="Seed for SU(2) tests (random.org)"),
):
    test_specifics()
    test_su4(iterations, su4_seed)
    test_su2_x_su2(iterations, su2_seed)
    print("[bold][green]Passed![/green][/bold]")

if __name__ == "__main__":
    typer.run(main)
Screen Shot 2022-09-12 at 11 38 42 AM

Though I run into problems when trying to make a similar change in the genos-diagonalos branch. With this diff:

diff --git a/src/compilers/approx.lisp b/src/compilers/approx.lisp
index 53abb86..f1faee1 100644
--- a/src/compilers/approx.lisp
+++ b/src/compilers/approx.lisp
@@ -212,46 +212,30 @@
     (magicl::map-to #'complex m cm)
     cm))

-(defun takagi-decomposition-of-uu^t (u)
+(defun decomposition-of-uu^t (u)
   "Given a unitary U, finds an X such that

     X^T (UU^T) X

 is a diagonal matrix. Return (VALUES X UU^T)."
-  (let* ((uut (magicl:@ u (magicl:transpose u)))
-         (n (magicl:nrows uut))
-         (a (magicl:.realpart uut))
-         (b (magicl:.imagpart uut))
-         (m (magicl:block-matrix (list (magicl:map #'- a) b b a) '(2 2))))
-    (multiple-value-bind (p d) (magicl:schur m)
-      (let* ((positive-eigs (loop :for i :below (magicl:nrows d)
-                                  :for e := (magicl:tref d i i)
-                                  :when (and (double~ 0.0d0 (imagpart e))
-                                             (plusp (realpart e)))
-                                    :collect i))
-             (diagonalizer (magicl:zeros (list n n) :type '(complex double-float))))
-        (assert (= 4 (length positive-eigs))
-            ()
-            "Expected 4 positive eigenvalues. Got ~D." (length positive-eigs))
-        (dotimes (row n)
-          (loop :for col :below n
-                :for from-col :in positive-eigs :do
-                  (setf (magicl:tref diagonalizer row col)
-                        (complex (magicl:tref p (+ n row) from-col)
-                                 (magicl:tref p row from-col)))))
-        (setf diagonalizer (magicl:transpose (magicl:inv diagonalizer)))
-        (when *check-math*
-          (assert (diagonal-matrix-p (magicl:@ (magicl:transpose diagonalizer)
-                                               uut
-                                               diagonalizer))))
-        (values diagonalizer uut)))))
+  (let* ((u-u^t (magicl:@ u (magicl:transpose u)))
+         (z (nth-value 0 (magicl:schur u-u^t)))
+         (eig-vals-z (nth-value 0 (magicl:eig z)))
+         (n (first (magicl:shape u)))
+         (eig-vals-mat (magicl:from-list
+                         (loop :with m = '()
+                               :for i :below n
+                               :do (setf m (append m eig-vals-z)) :finally (return m))
+                         (list n n)))
+         (x (magicl:./ z eig-vals-mat)))
+    (values x u-u^t)))

 (defun find-diagonalizer-in-e-basis (m)
   "For M in SU(4), compute an SO(4) column matrix of eigenvectors of E^* M E (E^* M E)^T."
   (check-type m magicl:matrix)
   (assert (magicl:unitary-matrix-p m))
   (let ((u (magicl:@ +edag-basis+ m +e-basis+)))
-    (multiple-value-bind (evecs gammag) (takagi-decomposition-of-uu^t u)
+    (multiple-value-bind (evecs gammag) (decomposition-of-uu^t u)

I receive the following complaint with make test:

CL-QUIL-TESTS (Suite)
  TEST-LOGICAL-MATRIX-SANITY
    I 0
    X 0
    Y 0
    Z 0
    X 0;X 0;
Unhandled SIMPLE-ERROR in thread #<SB-THREAD:THREAD "main thread" RUNNING
                                    {1001670003}>:
  The calculated eigenvectors were not found to be orthonormal. EE^T =
#<MATRIX/COMPLEX-DOUBLE-FLOAT (4x4):
  -0.712 + 0.702j     0.000 + 0.000j     0.000 - 0.000j     0.000 + 0.000j
   0.000 + 0.000j    -0.921 - 0.391j    -0.000 + 0.000j     0.000 + 0.000j
   0.000 - 0.000j    -0.000 + 0.000j     0.852 - 0.524j     0.000 - 0.000j
   0.000 + 0.000j     0.000 + 0.000j     0.000 - 0.000j     0.599 + 0.801j>
genos commented 2 years ago

The Gram-Schmidt process in orthonormalize-matrix! looks correct, so I suspect the check should be about $E E^\dagger$ instead of $EE^T$, because we want evecs to be unitary, not necessarily orthogonal (i.e. all real). However, adding the change

-        (assert (magicl:every #'double~
-                              (eye 4 :type 'double-float)
-                              (magicl:@ (magicl:transpose evecs)
-                                        evecs))
+        (assert (magicl:unitary-matrix-p evecs)
             (evecs)
-            "The calculated eigenvectors were not found to be orthonormal. ~
-               EE^T =~%~A"
-            (magicl:@ (magicl:transpose evecs)
-                      evecs)))
+            "The calculated eigenvectors were not found to be unitary. ~
+               EE^† =~%~A"
+            (magicl:@ (magicl:dagger evecs) evecs)))

leads to a different failure later on

CL-QUIL-TESTS (Suite)
  TEST-LOGICAL-MATRIX-SANITY
    I 0
    X 0
    Y 0
    Z 0
    X 0;X 0;
Unhandled SIMPLE-ERROR in thread #<SB-THREAD:THREAD "main thread" RUNNING
                                    {1001698003}>:
  Matrices do not lie in the same projective class.
 #<MATRIX/COMPLEX-DOUBLE-FLOAT (4x4):
    1.000 + 0.000j     0.000 + 0.000j     0.000 + 0.000j     0.000 + 0.000j
    0.000 + 0.000j     1.000 + 0.000j     0.000 + 0.000j     0.000 + 0.000j
    0.000 + 0.000j     0.000 + 0.000j     1.000 + 0.000j     0.000 + 0.000j
    0.000 + 0.000j     0.000 + 0.000j     0.000 + 0.000j     1.000 + 0.000j>
 #<MATRIX/COMPLEX-DOUBLE-FLOAT (4x4):
   -0.668 + 0.627j     0.155 + 0.085j     0.225 - 0.173j     0.066 + 0.211j
    0.155 + 0.085j     0.570 + 0.711j    -0.049 + 0.296j    -0.194 - 0.109j
    0.225 - 0.173j    -0.049 + 0.296j     0.759 - 0.415j     0.272 - 0.090j
    0.066 + 0.211j    -0.194 - 0.109j     0.272 - 0.090j    -0.841 - 0.334j>

So perhaps we do in fact want evecs to be a (real) orthogonal matrix rather than simply a (complex) unitary one.