sagemath / sage

Main repository of SageMath
https://www.sagemath.org
Other
1.41k stars 475 forks source link

unitary DFT for the symmetric group algebra #38456

Open jacksonwalters opened 2 months ago

jacksonwalters commented 2 months ago

Problem Description

At present the only two DFT's available for the symmetric group algebra are the seminormal basis and modular DFT. Neither of these are unitary transformations.

Proposed Solution

In Beals ['97], Quantum computation of Fourier transforms over symmetric groups he remarks that the DFT can be made unitary with two modifications:

  1. include a factor of \sqrt{d_\rho/|G|} in front of each Fourier coefficient at \rho
  2. the rep'ns being summed over must be themselves unitary

Alternatives Considered

The seminormal basis is unitary in the sense that A*A^T is diagonal meaning the columns are orthogonal. There may be some modification of it to yield orthonormal columns.

Additional Information

No response

Is there an existing issue for this?

jacksonwalters commented 1 week ago

Weyl's unitary trick effectively solves this for the complex numbers and number fields.

For finite fields, the situation is more delicate.

There is: "orthogonal" – for Young’s orthogonal representation for sage.combinat.symmetric_group_representations.SymmetricGroupRepresentation.

However, when ring=GF(7) is set, we obtain TypeErrors for partition=[2,1] corresponding to n=3. What appears to be happening is sqrt(3) is not able to be converted to an element of the finite field, nor are extensions being created. Perhaps it is possible to work over an extension so that these matrices become "orthogonal".

TypeError: self must be a numeric expression

jacksonwalters commented 1 week ago

PR: https://github.com/sagemath/sage/pull/38455

The error with square roots has been resolved by first using self._ring on beta before passing to 2x2 matrix element in symmetric_group_representations.py on line 662.

Regarding unitary representations of finite groups over finite fields, see:

https://mathoverflow.net/questions/327823/unitary-representations-of-finite-groups-over-finite-fields

One must check that the action of the q^th power of Frobenius on the Brauer character (which is just the complex character in the non-modular case) is complex conjugation. Even then, it is not immediately clear to me that Young's orthogonal representations will be unitary, or which sesquilinear form they will be unitary with respect to, exactly.

jacksonwalters commented 1 week ago

With the fix implemented locally, we do get orthogonal representations over GF(q^2). However, these representations are not unitary w.r.t. the form <x,y> = \sum_i x_i*y_i^q since x |--> x^q does not fix elements outside of GF(q). It may be that the matrices are unitary w.r.t. a different sesquilinear form.

dimpase commented 3 days ago

Perhaps the invariant form can be obtained by summing $g\overline{g}^\top$ over the group. (At least this works in characteristic 0 case, here it's a question of luck)

jacksonwalters commented 3 days ago

@dimpase I tried this after doing the characteristic zero case (complex numbers, number fields) and it didn't seem to work. This led to this question:

https://math.stackexchange.com/questions/4952613/is-there-a-version-of-weyls-unitary-trick-for-positive-characteristic

It does seem natural to try to just sum. One similar thought I haven't tried is averaging the sesquilinear form $<x,y>=\sum_i x_i y_i^q$ to get a $G$-invariant form.

Updates: representations of the symmetric group have Frobenius-Schur indicator +1, so Lemma 4.4.1 in the linked MO question doesn't apply. However, this question seems of relevance:

https://mathoverflow.net/questions/271932/formula-for-the-frobenius-schur-indicator-of-a-finite-group

This appears to imply (please check) that since $\chi$ is real, and irreducible rep'ns of the symmetric group are self-dual (theorem 1.5 in Gordon James' book), and we know the F-S indicator is +1, then as long as $d(\chi,\phi)$ is odd (it is often +1, see tables in Appendix of James' book), then there should exist a $G$-invariant symmetric bilinear form.

I just don't know how to find it. Perhaps I am overlooking something obvious.

dimpase commented 3 days ago

an obvious way to find invariant forms is by linear algebra in the representation $\rho(g)$. Take a matrix of unknowns $U$, and for each generator $g$ of $G$, write $\rho (g)U=\lambda_g U\overline{\rho( g)}^{-1^\top}$, where $\lambda_g:=\det \rho (g) \det \overline{\rho( g)}^{\top}$. In the solution space pick a nonzero element.

EDIT: nontrivial $\lambda_g$ arise when one has sesquilinear forms. (Corrected above). Also note that a convenient way to form the system of equations is to take Kronecker products $\rho (g)U\otimes \lambda_g U\overline{\rho( g)}^{-1^\top}$. Then $U$ can be recovered from a common, for all the generators, eigenvector with eigenvalue 1, of these Kronecker products. such an eigenvector gives you an embedding into the full unitary group of the corresponding form - I don't recall off top of my head whether it will be unique for the irreducible representation $\rho$ in positive characteristic. @jacksonwalters

jacksonwalters commented 2 days ago

@dimpase Ah, okay, great. Thanks very much! This is very concrete. I will attempt to implement this, ideally by the end of the week, and see what I can say.

dimpase commented 2 days ago

There is also a GAP package forms which does this in GAP, for matrix groups over finite fields. We probably should add it to the list of GAP packages we can easily install. For now, you can install gap_packages SPKG and then use InstallPackage to install forms.

jacksonwalters commented 2 days ago

I tried installing the GAP package forms via extracting the .tar to the appropriate directory ($SAGE_ROOT/local/share/gap/pkg/). It appears to load the package successfully, then I get all kinds of errors and claims that functions don't exist:

sage: gap('LoadPackage("forms");')
true
sage: gap('Packages();')
---------------------------------------------------------------------------
RuntimeError                              Traceback (most recent call last)
File /private/var/tmp/sage-10.4-current/local/var/lib/sage/venv-python3.12.4/lib/python3.12/site-packages/sage/interfaces/gap.py:700, in Gap_generic._eval_line(self, line, allow_use_file, wait_for_prompt, restart_if_needed)
    699     error = error.replace('\r', '')
--> 700     raise RuntimeError("%s produced error output\n%s\n   executing %s" % (self, error, line))
    701 if not normal:

RuntimeError: Gap produced error output
Error, Variable: 'Packages' must have a value

   executing \$sage2:=Packages();;;

During handling of the above exception, another exception occurred:

RuntimeError                              Traceback (most recent call last)
File /private/var/tmp/sage-10.4-current/local/var/lib/sage/venv-python3.12.4/lib/python3.12/site-packages/sage/interfaces/expect.py:1519, in ExpectElement.__init__(self, parent, value, is_name, name)
   1518 try:
-> 1519     self._name = parent._create(value, name=name)
   1520 # Convert ValueError and RuntimeError to TypeError for
   1521 # coercion to work properly.

File /private/var/tmp/sage-10.4-current/local/var/lib/sage/venv-python3.12.4/lib/python3.12/site-packages/sage/interfaces/interface.py:517, in Interface._create(self, value, name)
    516 name = self._next_var_name() if name is None else name
--> 517 self.set(name, value)
    518 return name

File /private/var/tmp/sage-10.4-current/local/var/lib/sage/venv-python3.12.4/lib/python3.12/site-packages/sage/interfaces/gap.py:1354, in Gap.set(self, var, value)
   1353 cmd = ('%s:=%s;;' % (var, value)).replace('\n', '')
-> 1354 self._eval_line(cmd, allow_use_file=True)

File /private/var/tmp/sage-10.4-current/local/var/lib/sage/venv-python3.12.4/lib/python3.12/site-packages/sage/interfaces/gap.py:734, in Gap_generic._eval_line(self, line, allow_use_file, wait_for_prompt, restart_if_needed)
    733     else:
--> 734         raise RuntimeError(exc)
    736 except KeyboardInterrupt:

RuntimeError: Gap produced error output
Error, Variable: 'Packages' must have a value

   executing \$sage2:=Packages();;;

During handling of the above exception, another exception occurred:

TypeError                                 Traceback (most recent call last)
Cell In[2], line 1
----> 1 gap('Packages();')

File /private/var/tmp/sage-10.4-current/local/var/lib/sage/venv-python3.12.4/lib/python3.12/site-packages/sage/interfaces/interface.py:299, in Interface.__call__(self, x, name)
    296         pass
    298 if isinstance(x, str):
--> 299     return cls(self, x, name=name)
    300 try:
    301     # Special methods do not and should not have an option to
    302     # set the name directly, as the identifier assigned by the
    303     # interface should stay consistent. An identifier with a
    304     # user-assigned name might change its value, so we return a
    305     # new element.
    306     result = self._coerce_from_special_method(x)

File /private/var/tmp/sage-10.4-current/local/var/lib/sage/venv-python3.12.4/lib/python3.12/site-packages/sage/interfaces/expect.py:1524, in ExpectElement.__init__(self, parent, value, is_name, name)
   1522 except (RuntimeError, ValueError) as x:
   1523     self._session_number = -1
-> 1524     raise TypeError(*x.args)
   1525 except BaseException:
   1526     self._session_number = -1

TypeError: Gap produced error output
Error, Variable: 'Packages' must have a value

   executing \$sage2:=Packages();;;
dimpase commented 2 days ago

Does it work in a separate GAP session? (start Sage shell: ./sage -sh and then type gap)

jacksonwalters commented 2 days ago

@dimpase Ah, that fixed it! Thank you.

gap> LoadPackage("forms");
---------------------------------------------------------------------
Loading 'Forms' 1.2.12 (30/08/2024)
by John Bamberg (http://school.maths.uwa.edu.au/~bamberg/)
   Jan De Beule (http://www.debeule.eu)
For help, type: ?Forms 
---------------------------------------------------------------------
true
gap> IsPackageLoaded("forms");
true
gap> gf := GF(8);
GF(2^3)
gap> r := PolynomialRing(gf,3);
GF(2^3)[x_1,x_2,x_3]
gap> poly := r.1^2 + r.2 * r.3;
x_1^2+x_2*x_3
gap> form := QuadraticFormByPolynomial(poly,r);
< quadratic form >

I'm still parsing the linear algebra approach. I generated some (non-functional) code but I need to tweak it. The main difficulty right now is just which rings to define these matrices over, in particular U which is symbolic. I want to solve these systems over GF(q^2) but I'm not quite sure how to do symbolic linear algebra in Sage.

TypeError: positive characteristic not allowed in symbolic computations

q = 7
n = 3
#partitions should be such that decomposition numbers d(\chi,\phi) is odd
partition = [2,1]

# Define the group G and its representation as a Specht module
SGA = SymmetricGroupAlgebra(GF(q^2), n)
SM = SGA.specht_module(partition)

# Get representation from Specht module
rho = SM.representation_matrix

# Define the matrix of unknowns U (2x2 matrix in this case)
n = 2  # Dimension of U
# Create symbolic variables in the symbolic ring (SR)
u_11, u_12, u_21, u_22 = var('u_11 u_12 u_21 u_22')
U = Matrix(SR, n, n, [[u_11, u_12], [u_21, u_22]])  # Use SR for symbolic ring

# Create a list to hold the Kronecker products
kron_products = []

# Iterate over each generator g in G
for g in G.gens():
    # Compute rho(g) - you need to replace this with your actual representation logic
    rho_g = rho(Permutation(g))
    rho_g_inv_T = rho_g.inverse().transpose()

    # Compute lambda_g
    det_rho_g = det(rho_g)  # This will be an element in GF(q)
    lambda_g = det_rho_g * (det_rho_g ** q)  # Correct conjugation operation

    # Convert lambda_g to a matrix for the Kronecker product
    lambda_g_matrix = lambda_g * identity_matrix(GF(q^2), n)

    # Form the Kronecker product
    # Note: Use rho_g and its inverse correctly in the tensor product
    kron_prod = (rho_g * U).tensor_product(lambda_g_matrix * U * rho_g_inv_T)
    kron_products.append(kron_prod)

# Now, we need to find the common eigenvector with eigenvalue 1
common_eigenvectors = []
for k in kron_products:
    eigenspaces = k.eigenvectors_right(1)  # Find eigenvectors for eigenvalue 1
    common_eigenvectors.extend(eigenspaces)

# Check if there are non-zero solutions
non_zero_solutions = [vec for eigens in common_eigenvectors for vec in eigens if vec[0] != 0]

# Display results
if non_zero_solutions:
    print("Found non-zero common eigenvectors:")
    for vec in non_zero_solutions:
        print(vec)
else:
    print("No non-zero common eigenvectors found.")
dimpase commented 2 days ago

you don't need to make U a variable! Just write down these Kronecker products, stack them up in a matrix, etc

anyway, the symbolic ring is the last resort, and it doesn't work in positive char. So all you have are polynomial rings.

dimpase commented 1 day ago

Computing an eigenvector $v$ of $A$ with a known eigenvalue $\lambda$ is the same as finding a vector in the nullspace of $A-\lambda I$, and the nullspace (or kernel) are methods of Sage matrices. So you can stack up these matrices $A-\lambda I$ and call the (right?) nullspace for this non-square matrix. This way avoids dealing with explicit variables etc.

jacksonwalters commented 1 day ago

@dimpase Yes, one can find eigenvalues by computing $\text{null}(A-\lambda I)$ and we just want the common $\lambda=1$ eigenspace for each $g \in G$, so we can stack the $A-\lambda I$ where $A=\rho(g)U\otimes \lambda_g U\overline{\rho(g)}^{-1 \top}$.

Can you just write in code how you would instantiate $U$? You seem to be saying two contradictory things:

1) $U$ is a matrix of unknowns 2) $U$ is not a variable

So I don't know how to create $U$, therefore any of these matrices because U=matrix(GF(q^2),d_rho,d_rho) is just the zero matrix. I mean, it's clear that this results in a system of linear equations, but I don't see how it wouldn't be in the variables $u_{ij}$, so I'm confused.

Also, you write $\det(\overline{\rho(g)}^{\top})$. Isn't the determinant of the transpose just the determinant?

EDIT: Okay, beginning to see. We are just setting $\lambda=1$, so this system can be arranged so that $U$ is the unknown vector of a linear system which is non-square, so we just need to find the matrix itself. $U$ won't appear explicitly, e.g. when solving $Ax = b$ we use the augmented system $[A|b]$ with no mention of $x$. I think I do need to work over a polynomial ring $GF(q^2)[u0,...,u{d\rho^2}]$ though.

dimpase commented 1 day ago

yeah, right, matrix of unknowns was there for the purpose of explanation. By the way, I am not sure that $\lambda=1$ always holds in sesquilinear case. There perhaps could be different $\lambda$ for different groups generators.

In the end you can fill in the entries of $U$ with the solution entries of the system, no need for variables

jacksonwalters commented 1 day ago

@dimpase Since $U$ occurs in both factors of the Kroenecker product, don't we get terms quadratic in the $u_{ij}$?

dimpase commented 1 day ago

no, you "vectorise" $U$ and get its $n^2$ entries as a vector. So you need to fill them back in as a matrix.

jacksonwalters commented 1 day ago

Yes, we can write $U=(u0, \ldots, u{n^2})$ as a vector where $n=d_\rho$ , but when I do $\rho(g)U \otimes \lambdag U \rho(g)$ I am forming the block matrix with first block $(\rho(g)U){11}\lambda_g U \rho(g)$ so I will have terms like $u_0^2$.

dimpase commented 23 hours ago

Yes, we can write U = ( u 0 , … , u n 2 ) as a vector where n = d ρ , but when I do ρ ( g ) U ⊗ λ g U ρ ( g ) I am forming the block matrix with first block ( ρ ( g ) U ) 11 λ g U ρ ( g ) so I will have terms like u 0 2 .

Why do you need to do ρ ( g ) U ⊗ λ g U ρ ( g ) ? All I was saying is that the action of the group on the space of $n\times n$ matrices defining sesquilinear forms is the tensor product of the suitable $n$-dimensional representations. So for fixed $A$, $B$ the system of linear equations $AU=UB$ in $(u{ij})=U$ can be written as something like $(A\otimes B) u=0$, with $u=(u{11},u{12},\dots, u{nn})$. You solve the latter system to obtain a solution of the former. $u_{ij}$ as explicit variables don't appear anywhere.

jacksonwalters commented 22 hours ago

Here is a version without using the Kroenecker products which should work (at least in cases where the decomposition number is odd) for n=3, q=3, partition = [2,1].

EDIT: Corrected swapped indices in conjugate matrix.

# Setup
q = 3
n = 3
partition = [2, 1]

# Define the group G and its representation as a Specht module
SGA = SymmetricGroupAlgebra(GF(q^2), n)
SM = SGA.specht_module(partition)
G = SGA.group()

# Get representation from Specht module
rho = SM.representation_matrix
d_rho = SM.dimension()

# Initialize U as a matrix of variables over GF(q^2)
R = PolynomialRing(GF(q^2), 'u', d_rho^2)
U_vars = R.gens()  # List of variable generators for U
U = Matrix(R, d_rho, d_rho, U_vars)  # U is a d_rho x d_rho matrix of variables

#define conjugation as x |--> x**q, an order two automorphism of F_q^2. note x**q == x for x \in F_q.
def conjugate_pos_char(A):
    assert A.nrows() == A.ncols()
    field_size = A.base_ring().order()
    q = sqrt(field_size) if field_size.is_square() else field_size
    return matrix(GF(q**2),[[A[i][j]**q for j in range(A.nrows())] for i in range(A.nrows())])

# for each generator of G, form the augmented system 
def augmented_matrix(g):

    rho_g = rho(Permutation(g))
    rho_g_conj_inv_T = conjugate_pos_char(rho_g.inverse().transpose())

    # Compute lambda_g
    det_rho_g = det(rho_g)
    lambda_g = det_rho_g * (det_rho_g ** q)

    # Form the matrix equation
    equation_matrix = rho_g * U - lambda_g * U * rho_g_conj_inv_T

    # Initialize a list to hold rows of the augmented system
    augmented_system = []

    # Extract coefficients for each linear equation in the matrix
    for i in range(d_rho):
        for j in range(d_rho):
            # Get the (i, j) entry of the equation matrix, which is a linear combination of the u variables
            linear_expression = equation_matrix[i, j]

            # Extract the coefficients of each u_k in the linear expression
            row = [linear_expression.coefficient(u) for u in U_vars]

            # Append the row to the augmented system
            augmented_system.append(row)

    # Convert the augmented system to a matrix
    A = Matrix(GF(q^2), augmented_system)

    return A

total_system = matrix(GF(q^2),0,d_rho^2)
for g in G:
    total_system = total_system.stack(augmented_matrix(g))

null_space = total_system.right_kernel()

print(null_space)

Vector space of degree 4 and dimension 0 over Finite Field in z2 of size 3^2
Basis matrix:
[]
dimpase commented 22 hours ago

You made a typo in conjugate_pos_char(A). Its last line should be

return matrix(GF(q**2),[[A[i][j]**q for j in range(A.nrows())] for i in range(A.nrows())])

i.e the inner loop is over j. With this change your code prints

Vector space of degree 4 and dimension 1 over Finite Field in z2 of size 3^2
Basis matrix:
[1 1 1 1]
jacksonwalters commented 13 hours ago

Ah, great! Thanks so much for spotting that. Yes, getting lots of interesting forms for different choices of n and partition, e.g. q=3, n=5, [3,1,1]:

Vector space of degree 36 and dimension 1 over Finite Field in z2 of size 3^2
Basis matrix:
[0 1 2 1 2 0 1 0 1 1 0 2 2 1 0 0 1 2 1 1 0 0 1 1 2 0 1 1 0 1 0 2 2 1 1 0]

I'll do some testing and verification to ensure these really define G-invariant sesquilinear forms tomorrow, and for which cases are we getting solutions. If for a given q, n there are solutions for every partition, it may be possible to construct a unitary DFT.