在这篇笔记介绍了计算出给定有限群的所有不可约表示特征标(复数域上)的算法。并介绍了Dixon算法,约化可约表示。

计算特征标

对于有限群,考虑他的 conjugation class 一共 个 conjugation class。在群代数 中,定义 conjugation class 的平均值是

这里 表示集合 中元素的个数。于是

这里 , ,证明见附录。

在不可约表示中, 被表示成矩阵 , 对应的矩阵是

根据 conjugation class 的定义,对 ,有

是不可约表示,由 schur’s lemma, 只能等于 , 是单位矩阵。对 () 式取 trace, 得到 , 所以 , 是表示空间的维数。

() 式对应的表示矩阵等式是

定义矩阵形式 ,那么上式可看作矩阵的特征方程

特征向量是 ,他是 这些矩阵的共同特征向量。

根据有限群表示的性质:不等价不可约表示的个数等于 conjugation class 的个数,且不同特征标之间满足正交完备关系,而 个矩阵都是 大小的,如果他们具有唯一的一组 个共同特征向量,那么从这组特征向量就可以直接计算出 个不可约表示的特征标了。根据 可以确定 的相位和归一化系数。这就是Burnside算法。

计算程序如下

1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
import numpy as np

class FiniteGroupRepresentation:
def __init__(self, multi_table, identity=1):
"""
Initialize with a multiplication table of the group.
@param multi_table: 2D list or numpy array representing the group's multiplication table.
@param identity: The identity element of the group (default is 1).
"""
self.multi_table = multi_table
self.g_order = len(multi_table)
self.identity = identity # Assuming the identity element is represented by 1

self.__conjugacy_classes = None # List of sets representing conjugacy classes
self.__conj_class_multi_coeff = None # 3D numpy array for conjugacy class multiplication coefficients

def find_conjugacy_classes(self)->list[set[int]]:
"""
Find and return the conjugacy classes of the group.
@return: A list of sets, each set representing a conjugacy class.
"""
visited = [False] * (self.g_order + 1)
self.__conjugacy_classes = []

for g in range(1, self.g_order + 1):
if not visited[g]:
conj_class = set()
for x in range(1, self.g_order + 1):
# Compute xgx^{-1}
x_inv = self.find_inverse(x)
conj_element = self.multi_table[self.multi_table[x - 1][g - 1] - 1][x_inv - 1]
conj_class.add(conj_element)
for elem in conj_class:
visited[elem] = True
self.__conjugacy_classes.append(conj_class)

return self.__conjugacy_classes

def find_inverse(self, element):
"""
Find the inverse of a given element in the group.
@param element: The element to find the inverse for.
@return: The inverse element.
"""
for i in range(1, self.g_order + 1):
if self.multi_table[element - 1][i - 1] == self.identity:
return i
raise ValueError(f"No inverse found for element {element}")

def compute_conj_class_multi_coeff(self):
"""
Compute the conjugacy class multiplication coefficients.
@return: A 3D numpy array where entry [i][j][k] represents the coefficient for the product of
conjugacy classes i and j resulting in class k.
"""
if self.__conjugacy_classes is None:
self.__conjugacy_classes = self.find_conjugacy_classes()
conj_classes = self.__conjugacy_classes
num_conj_classes = len(conj_classes)
self.__conj_class_multi_coeff = np.zeros((num_conj_classes, num_conj_classes, num_conj_classes), dtype=int)
for i, class_i in enumerate(conj_classes):
for j, class_j in enumerate(conj_classes):
for g in class_i:
for h in class_j:
gh = self.multi_table[g - 1][h - 1]
for k, class_k in enumerate(conj_classes):
if gh in class_k:
self.__conj_class_multi_coeff[i][j][k] += 1
break
for k in range(num_conj_classes):
for i in range(num_conj_classes):
for j in range(num_conj_classes):
if self.__conj_class_multi_coeff[i][j][k] % len(conj_classes[k]) != 0:
raise ValueError("Conjugacy class matrix entry not divisible by class size")
self.__conj_class_multi_coeff[i][j][k] //= len(conj_classes[k])
return self.__conj_class_multi_coeff


def __find_common_eigenvectors(self, matrices, tol=1e-12, index = 0, current_sol=None) -> dict[tuple[float,...], list[np.ndarray]]:
"""
Find common eigenvectors for a list of matrices.
@param matrices: List of numpy arrays (matrices).
@param tol: Tolerance for numerical comparisons.
@param index: Current index in the recursion.
@param current_sol: Current solution dictionary mapping eigenvalue tuples to lists of eigenvectors.
@return: A dictionary mapping tuples of eigenvalues to lists of common eigenvectors.
"""
ndigits = int(-np.log10(tol))
if index >= len(matrices):
sol = dict()
for eigval, eigvecs in current_sol.items():
if len(eigvecs) != 1:
raise ValueError("Degenerate eigenvalue found; cannot proceed.")
sol[eigval] = eigvecs[0]
return sol
if index == 0 or current_sol is None:
sol = dict()
eigval, eigvec = np.linalg.eig(matrices[0])
for i, val in enumerate(eigval):
key = (round(val, ndigits),)
if key not in sol:
sol[key] = []
sol[key].append(eigvec[:, i])
return self.__find_common_eigenvectors(matrices, tol, 1, sol)

sol = dict()
for eigval in current_sol:
eigvecs = current_sol[eigval]

# solve eigen problem in the subspace spanned by degenerate eigenvectors
project_dim = len(eigvecs)

# Orthonormalize eigvecs using Gram-Schmidt process
orthonorm_eigvecs = []
for v in eigvecs:
w = v.copy()
for u in orthonorm_eigvecs:
w -= np.dot(u.conj().T, w) * u
norm = np.linalg.norm(w)
if norm > tol:
orthonorm_eigvecs.append(w / norm)
eigvecs = orthonorm_eigvecs

projected_mat = np.zeros((project_dim, project_dim), dtype=np.complex128)
for i in range(project_dim):
for j in range(project_dim):
projected_mat[i][j] = np.dot(eigvecs[i].conj().T, np.dot(matrices[index], eigvecs[j]))
eigval_new, eigvec_new = np.linalg.eig(projected_mat)
for i, val in enumerate(eigval_new):
newkey = eigval + (np.round(val, ndigits),)
if newkey not in sol:
sol[newkey] = []
sol[newkey].append(sum(eigvecs[j] * eigvec_new[j][i] for j in range(project_dim)))
return self.__find_common_eigenvectors(matrices, tol, index + 1, sol)


def compute_character_table(self):
"""
Compute the character table of the group.
@return: A dictionary mapping tuples of eigenvalues to normalized character vectors.
[tuple of eigenvalues] -> [character vector as numpy array]
"""
if self.__conj_class_multi_coeff is None:
self.compute_conj_class_multi_coeff()
# find common eigenvectors of the conjugacy class multiplication coefficient matrices
common_eigvecs = self.__find_common_eigenvectors(self.__conj_class_multi_coeff)
# compute character table
for i, eigval in enumerate(common_eigvecs):
for j in range(len(common_eigvecs[eigval])):
common_eigvecs[eigval][j] /= len(self.__conjugacy_classes[j])
common_eigvecs[eigval] /= common_eigvecs[eigval][0] # set character of identity to 1
# normalize character vectors
norm = sum([len(self.__conjugacy_classes[k]) * abs(common_eigvecs[eigval][k])**2 for k in range(len(self.__conjugacy_classes))]) / self.g_order
common_eigvecs[eigval] /= np.sqrt(norm)
for i, eigval in enumerate(common_eigvecs):
print(f"Representation {i+1}, Eigenvalues: {eigval},\n\tCharacter: {common_eigvecs[eigval]}")
return common_eigvecs


if __name__ == "__main__":
import time
time_start = time.time()
multi_table = [[1,2,3,4,5,6],
[2,1,4,3,6,5],
[3,5,1,6,2,4],
[4,6,2,5,1,3],
[5,3,6,1,4,2],
[6,4,5,2,3,1]]
group = FiniteGroupRepresentation(multi_table)
conj_classes = group.find_conjugacy_classes()
for i, cls in enumerate(conj_classes):
print(f"Conjugacy Class {i+1}: {cls}")

conj_matrix = group.compute_conj_class_multi_coeff()
print("Conjugacy class multi coeff matrix:")
print(conj_matrix)

group.compute_character_table()

time_end = time.time()
print('time cost', time_end - time_start, 's')

作为例子,输入3阶置换群 的乘法表,程序输出结果如下

1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
Conjugacy Class 1: {1}
Conjugacy Class 2: {2, 3, 6}
Conjugacy Class 3: {4, 5}
Conjugacy class multi coeff matrix:
[[[1 0 0]
[0 1 0]
[0 0 1]]

[[0 1 0]
[3 0 3]
[0 2 0]]

[[0 0 1]
[0 2 0]
[2 0 1]]]
Representation 1, Eigenvalues: (1.0, (3+0j), (2+0j)),
Character: [1.+0.j 1.+0.j 1.+0.j]
Representation 2, Eigenvalues: (1.0, (-0+0j), (-1+0j)),
Character: [ 2.00000000e+00+0.j 2.01951476e-16+0.j -1.00000000e+00+0.j]
Representation 3, Eigenvalues: (1.0, (-3+0j), (2+0j)),
Character: [ 1.-0.j -1.+0.j 1.-0.j]
time cost 0.0007150173187255859 s

得到 的特征标表

表示的编号 (1) (2,3,6) (4,5)
1 1 1 1
2 2 0 -1
3 1 -1 1

从特征标约化正则表示

有了特征标之后,就可以构造出 idempotent (投影算子) ,约化正则表示。

由特征标的正交性质,可以看出 是完备的,

于是正则表示被约化为 , 不过除非 是 Abelian group,否则这并不能直接将 约化为一组不可约表示。

Dixon算法约化表示

Dixon给出了一个可以约化可约表示的算法。

对于一个unitray表示 。对于 ,定义

其中 是一个 的矩阵,除了第 行第 列元素为 外,其余元素都是 。那么 的 Hermitian basis。

对每一个 ,计算出矩阵

容易验证 和所有的 对易。如果 是不可约表示, 就只能是一个常数乘以单位矩阵。否则,一定会存在某个 ,使得 不正比于单位阵。这是因为所有的 对易的 Hermitian 矩阵 都可以写成 () 式的形式,,而 是 Hermitian basis。

找出这个不正比于 后,将其对角化 ,那么在 的变换下,所有的 也会成为分块矩阵,这样就将表示空间约化了。对于约化后的子空间,继续递归使用该算法,就可以将整个空间约化成若干个不可约表示。

程序实现如下

1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
import numpy as np

class ReduceGroupRepresentation:
def __init__(self) -> None:
pass


def reduce_representation(self, group_rep: list[np.ndarray], tol=1e-12) -> list[list[np.ndarray]]:
"""
Reduce a given group representation into its irreducible components.
@param group_rep: List of numpy arrays representing the group representation matrices.
@param tol: Tolerance for numerical comparisons.
@return: List of lists, where each sublist contains numpy arrays representing an irreducible representation.
"""
self.irrep = []
self.conj_table = []
self.__reduce_representation_impl(group_rep, tol)
return self.irrep


def __add_to_irr_rep(self, rep: list[np.ndarray]):
"""
Add a new irreducible representation if it is not already present.
@param rep: List of numpy arrays representing the irreducible representation matrices.
"""
conj_table = []
for Dg in rep:
conj_table.append(np.trace(Dg))
for i, table in enumerate(self.conj_table):
if np.allclose(table, conj_table):
return
self.conj_table.append(conj_table)
self.irrep.append(rep)


def __reduce_representation_impl(self, group_rep: list[np.ndarray], tol=1e-12):
"""
Reduce a given group representation into its irreducible components.
@param group_rep: List of numpy arrays representing the group representation matrices.
@param tol: Tolerance for numerical comparisons.
@return: List of lists, where each sublist contains numpy arrays representing an irreducible representation.

Note: This implementation uses the Dixon method for reducing representations.
"""

# Dixon method implementation
dim = group_rep[0].shape[0]
assert all(mat.shape == (dim, dim) for mat in group_rep), "All representation matrices must have the same dimensions."

is_irr = self.is_irr_rep(group_rep, tol)
print(f"Is irreducible: {is_irr}")
if is_irr:
self.__add_to_irr_rep(group_rep)
return

res_reps = self.__decompose_into_block(group_rep, tol)

for rep in res_reps:
self.__reduce_representation_impl(rep, tol)


def __get_hermitian_basis(self, r, s, dim):
"""
Get a Hermitian basis matrix with 1 at (r,s) and (s,r) positions.
@param r: Row index.
@param s: Column index.
@param dim: Dimension of the square matrix.
@return: Hermitian basis matrix.
"""
hermitian_basis = np.zeros((dim, dim), dtype=complex)
if r == s:
hermitian_basis[r, s] = 1
elif r > s:
hermitian_basis[r, s] = 1
hermitian_basis[s, r] = 1
else:
hermitian_basis[r, s] = 1j
hermitian_basis[s, r] = -1j
return hermitian_basis


def __decompose_into_block(self, group_rep: list[np.ndarray], tol=1e-12) -> list[list[np.ndarray]]:
"""
Decompose a group representation into block diagonal form.
@param group_rep: List of numpy arrays representing the group representation matrices.
@param tol: Tolerance for numerical comparisons.
@return: List of lists, where each sublist contains numpy arrays representing a block of the representation.
"""
# If not irreducible, proceed to find invariant subspaces
# find basis transformation to block diagonal form
dim = group_rep[0].shape[0]
g_order = len(group_rep)
mat_sum = []
found_hermtian_basis = False
for r in range(dim):
if found_hermtian_basis:
break
for s in range(dim):
mat_sum = np.zeros((dim, dim), dtype=complex)
hermitian_basis = self.__get_hermitian_basis(r, s, dim)
for g in range(g_order):
mat_sum += group_rep[g].conj().T @ hermitian_basis @ group_rep[g]
# Check if mat_sum is a scalar multiple of the identity matrix
if not np.allclose(mat_sum, np.trace(mat_sum) / dim * np.eye(dim), atol=tol):
found_hermtian_basis = True
break
# Find eigenvalues and eigenvectors of Hermitian matrix mat_sum
eigvals, eigvecs = np.linalg.eigh(mat_sum)
print(eigvals)
rep_vectors = []
for rep_mat in group_rep:
rep_vectors.append(eigvecs.conj().T @ rep_mat @ eigvecs)

# import matplotlib.pyplot as plt
# for i in range(len(rep_vectors)):
# plt.subplot(1, len(rep_vectors), i+1)
# plt.imshow(np.abs(rep_vectors[i]))
# plt.colorbar()
# plt.title(f'group elem {i+1}')
# plt.show()

# split into blocks
block_sizes = []
last_idx = 0
for i in range(1, dim):
if eigvals[i] - eigvals[i-1] > tol:
block_sizes.append(i - last_idx)
last_idx = i
if last_idx < dim:
block_sizes.append(dim - last_idx)
print("Block sizes:", block_sizes)
res_reps = []
start_idx = 0
for size in block_sizes:
sub_rep = []
for rep_mat in rep_vectors:
sub_rep.append(rep_mat[start_idx:start_idx+size, start_idx:start_idx+size])
res_reps.append(sub_rep)
start_idx += size
return res_reps


def is_irr_rep(self, group_rep: list[np.ndarray], tol=1e-12) -> bool:
"""
Check if a given group representation is irreducible.
@param group_rep: List of numpy arrays representing the group representation matrices.
@param tol: Tolerance for numerical comparisons.
@return: True if the representation is irreducible, False otherwise.
Note: This implementation uses the character orthogonality relation to determine irreducibility.
"""
g_order = len(group_rep)
dim = group_rep[0].shape[0]
assert all(mat.shape == (dim, dim) for mat in group_rep), "All representation matrices must have the same dimensions."

charac_sum = 0
for g in range(g_order):
charac_sum += np.abs(np.trace(group_rep[g]))**2
charac_sum /= g_order
print(f"charac sum {charac_sum}")
return np.abs(charac_sum - 1) < tol


if __name__ == "__main__":
# multiplication table for S3
multi_table = [[1,2,3,4,5,6],
[2,1,4,3,6,5],
[3,5,1,6,2,4],
[4,6,2,5,1,3],
[5,3,6,1,4,2],
[6,4,5,2,3,1]]
from finite_group_rep import FiniteGroupRepresentation
group = FiniteGroupRepresentation(multi_table)
conj_class = group.find_conjugacy_classes()
print("Conjugacy classes:", conj_class)

### construct regular representation
reg_rep = []
size = len(multi_table)
for i in range(size):
mat = np.zeros((size, size), dtype=int)
for j in range(size):
mat[multi_table[i][j]-1][j] = 1
reg_rep.append(mat)
print("Regular representation matrices:")
for mat in reg_rep:
print(mat)
import time
start_time = time.time()

reducer = ReduceGroupRepresentation()
irreps = reducer.reduce_representation(reg_rep)
print(f"Found {len(irreps)} irreducible representations:")
for i, irrep in enumerate(irreps):
print(f"Irreducible Representation {i+1}: dim = {irrep[0].shape[0]}")
for mat in irrep:
print(mat)

end_time = time.time()
print(f"Time taken: {end_time - start_time} seconds")

附录

证明

, 记 是所有乘积等于 组成的集合。于是

所属的 conjugation class 中有另一个不同于 的元素 ,那么 之间存在一一对应 ,所以 ,可以将 所属 conjugation class 的下标 标记为 。所以

注意到

命题得证。

验证 是一组正交的 idempotent

验证是正交的 idempotent:

其中用到了群表示的正交性质