Linear Algebra

Linear algebra on Galois field arrays is supported through both native NumPy linear algebra function calls and additional linear algebra routines not included in NumPy.

In the sections below, the prime field \(\mathrm{GF}(31)\) is used.

In [1]: GF = galois.GF(31)

In [2]: print(GF)
Galois Field:
  name: GF(31)
  characteristic: 31
  degree: 1
  order: 31
  irreducible_poly: x + 28
  is_primitive_poly: True
  primitive_element: 3

NumPy routines

Most NumPy linear algebra routines are supported on vectors and matrices over finite fields. These routines are accessed using the standard numpy.linalg functions.

Expand any section for more details.

Dot product: np.dot(a, b)
In [3]: a = GF([29, 0, 2, 21]); a
Out[3]: GF([29,  0,  2, 21], order=31)

In [4]: b = GF([23, 5, 15, 12]); b
Out[4]: GF([23,  5, 15, 12], order=31)

In [5]: np.dot(a, b)
Out[5]: GF(19, order=31)
Vector dot product: np.vdot(a, b)
In [6]: a = GF([29, 0, 2, 21]); a
Out[6]: GF([29,  0,  2, 21], order=31)

In [7]: b = GF([23, 5, 15, 12]); b
Out[7]: GF([23,  5, 15, 12], order=31)

In [8]: np.vdot(a, b)
Out[8]: GF(19, order=31)
Inner product: np.inner(a, b)
In [9]: a = GF([29, 0, 2, 21]); a
Out[9]: GF([29,  0,  2, 21], order=31)

In [10]: b = GF([23, 5, 15, 12]); b
Out[10]: GF([23,  5, 15, 12], order=31)

In [11]: np.inner(a, b)
Out[11]: GF(19, order=31)
Outer product: np.outer(a, b)
In [12]: a = GF([29, 0, 2, 21]); a
Out[12]: GF([29,  0,  2, 21], order=31)

In [13]: b = GF([23, 5, 15, 12]); b
Out[13]: GF([23,  5, 15, 12], order=31)

In [14]: np.outer(a, b)
Out[14]: 
GF([[16, 21,  1,  7],
    [ 0,  0,  0,  0],
    [15, 10, 30, 24],
    [18, 12,  5,  4]], order=31)
Matrix multiplication: A @ B == np.matmul(A, B)
In [15]: A = GF([[17, 25, 18, 8], [7, 9, 21, 15], [6, 16, 6, 30]]); A
Out[15]: 
GF([[17, 25, 18,  8],
    [ 7,  9, 21, 15],
    [ 6, 16,  6, 30]], order=31)

In [16]: B = GF([[8, 18], [22, 0], [7, 8], [20, 24]]); B
Out[16]: 
GF([[ 8, 18],
    [22,  0],
    [ 7,  8],
    [20, 24]], order=31)

In [17]: A @ B
Out[17]: 
GF([[11, 22],
    [19,  3],
    [19,  8]], order=31)

In [18]: np.matmul(A, B)
Out[18]: 
GF([[11, 22],
    [19,  3],
    [19,  8]], order=31)
Matrix exponentiation: np.linalg.matrix_power(A, z)
In [19]: A = GF([[14, 1, 5], [3, 23, 6], [24, 27, 4]]); A
Out[19]: 
GF([[14,  1,  5],
    [ 3, 23,  6],
    [24, 27,  4]], order=31)

In [20]: np.linalg.matrix_power(A, 3)
Out[20]: 
GF([[ 1, 16,  4],
    [11,  9,  9],
    [ 8, 24, 29]], order=31)

In [21]: A @ A @ A
Out[21]: 
GF([[ 1, 16,  4],
    [11,  9,  9],
    [ 8, 24, 29]], order=31)
Matrix determinant: np.linalg.det(A)
In [22]: A = GF([[23, 11, 3, 3], [13, 6, 16, 4], [12, 10, 5, 3], [17, 23, 15, 28]]); A
Out[22]: 
GF([[23, 11,  3,  3],
    [13,  6, 16,  4],
    [12, 10,  5,  3],
    [17, 23, 15, 28]], order=31)

In [23]: np.linalg.det(A)
Out[23]: GF(0, order=31)
Matrix rank: np.linalg.matrix_rank(A, z)
In [24]: A = GF([[23, 11, 3, 3], [13, 6, 16, 4], [12, 10, 5, 3], [17, 23, 15, 28]]); A
Out[24]: 
GF([[23, 11,  3,  3],
    [13,  6, 16,  4],
    [12, 10,  5,  3],
    [17, 23, 15, 28]], order=31)

In [25]: np.linalg.matrix_rank(A)
Out[25]: 3

In [26]: A.row_reduce()
Out[26]: 
GF([[ 1,  0,  0, 11],
    [ 0,  1,  0, 25],
    [ 0,  0,  1, 11],
    [ 0,  0,  0,  0]], order=31)
Matrix trace: np.trace(A)
In [27]: A = GF([[23, 11, 3, 3], [13, 6, 16, 4], [12, 10, 5, 3], [17, 23, 15, 28]]); A
Out[27]: 
GF([[23, 11,  3,  3],
    [13,  6, 16,  4],
    [12, 10,  5,  3],
    [17, 23, 15, 28]], order=31)

In [28]: np.trace(A)
Out[28]: GF(0, order=31)

In [29]: A[0,0] + A[1,1] + A[2,2] + A[3,3]
Out[29]: GF(0, order=31)
Solve a system of equations: np.linalg.solve(A, b)
In [30]: A = GF([[14, 21, 14, 28], [24, 22, 23, 23], [16, 30, 26, 18], [4, 23, 18, 3]]); A
Out[30]: 
GF([[14, 21, 14, 28],
    [24, 22, 23, 23],
    [16, 30, 26, 18],
    [ 4, 23, 18,  3]], order=31)

In [31]: b = GF([15, 11, 6, 29]); b
Out[31]: GF([15, 11,  6, 29], order=31)

In [32]: x = np.linalg.solve(A, b)

In [33]: A @ x == b
Out[33]: array([ True,  True,  True,  True])
Matrix inverse: np.linalg.inv(A)
In [34]: A = GF([[14, 21, 14, 28], [24, 22, 23, 23], [16, 30, 26, 18], [4, 23, 18, 3]]); A
Out[34]: 
GF([[14, 21, 14, 28],
    [24, 22, 23, 23],
    [16, 30, 26, 18],
    [ 4, 23, 18,  3]], order=31)

In [35]: A_inv = np.linalg.inv(A); A_inv
Out[35]: 
GF([[27, 17,  9,  8],
    [20, 21, 12,  4],
    [30, 10, 23, 22],
    [13, 25,  6, 13]], order=31)

In [36]: A @ A_inv
Out[36]: 
GF([[1, 0, 0, 0],
    [0, 1, 0, 0],
    [0, 0, 1, 0],
    [0, 0, 0, 1]], order=31)

Additional routines

Row space: A.row_space()
In [37]: A = GF([[23, 11, 3, 3], [13, 6, 16, 4], [12, 10, 5, 3], [17, 23, 15, 28]]); A
Out[37]: 
GF([[23, 11,  3,  3],
    [13,  6, 16,  4],
    [12, 10,  5,  3],
    [17, 23, 15, 28]], order=31)

In [38]: A.row_space()
Out[38]: 
GF([[ 1,  0,  0, 11],
    [ 0,  1,  0, 25],
    [ 0,  0,  1, 11]], order=31)

See galois.FieldArray.row_space() for more details.

Column space: A.column_space()
In [39]: A = GF([[23, 11, 3, 3], [13, 6, 16, 4], [12, 10, 5, 3], [17, 23, 15, 28]]); A
Out[39]: 
GF([[23, 11,  3,  3],
    [13,  6, 16,  4],
    [12, 10,  5,  3],
    [17, 23, 15, 28]], order=31)

In [40]: A.column_space()
Out[40]: 
GF([[ 1,  0,  0,  0],
    [ 0,  1,  0, 11],
    [ 0,  0,  1,  5]], order=31)

See galois.FieldArray.column_space() for more details.

Left null space: A.left_null_space()
In [41]: A = GF([[23, 11, 3, 3], [13, 6, 16, 4], [12, 10, 5, 3], [17, 23, 15, 28]]); A
Out[41]: 
GF([[23, 11,  3,  3],
    [13,  6, 16,  4],
    [12, 10,  5,  3],
    [17, 23, 15, 28]], order=31)

In [42]: A.left_null_space()
Out[42]: GF([[ 0,  1, 23, 14]], order=31)

See galois.FieldArray.left_null_space() for more details.

Null space: A.null_space()
In [43]: A = GF([[23, 11, 3, 3], [13, 6, 16, 4], [12, 10, 5, 3], [17, 23, 15, 28]]); A
Out[43]: 
GF([[23, 11,  3,  3],
    [13,  6, 16,  4],
    [12, 10,  5,  3],
    [17, 23, 15, 28]], order=31)

In [44]: A.null_space()
Out[44]: GF([[ 1, 22,  1, 14]], order=31)

See galois.FieldArray.null_space() for more details.

Gaussian elimination: A.row_reduce()
In [45]: A = GF([[23, 11, 3, 3], [13, 6, 16, 4], [12, 10, 5, 3], [17, 23, 15, 28]]); A
Out[45]: 
GF([[23, 11,  3,  3],
    [13,  6, 16,  4],
    [12, 10,  5,  3],
    [17, 23, 15, 28]], order=31)

In [46]: A.row_reduce()
Out[46]: 
GF([[ 1,  0,  0, 11],
    [ 0,  1,  0, 25],
    [ 0,  0,  1, 11],
    [ 0,  0,  0,  0]], order=31)

See galois.FieldArray.row_reduce() for more details.

LU decomposition: A.lu_decompose()
In [47]: A = GF([[4, 1, 24], [7, 6, 1], [11, 20, 2]]); A
Out[47]: 
GF([[ 4,  1, 24],
    [ 7,  6,  1],
    [11, 20,  2]], order=31)

In [48]: L, U = A.lu_decompose()

In [49]: L
Out[49]: 
GF([[ 1,  0,  0],
    [25,  1,  0],
    [26, 15,  1]], order=31)

In [50]: U
Out[50]: 
GF([[ 4,  1, 24],
    [ 0, 12, 21],
    [ 0,  0, 24]], order=31)

In [51]: np.array_equal(L @ U, A)
Out[51]: True

See galois.FieldArray.lu_decompose() for more details.

PLU decomposition: A.plu_decompose()
In [52]: A = GF([[15, 4, 11], [7, 6, 1], [11, 20, 2]]); A
Out[52]: 
GF([[15,  4, 11],
    [ 7,  6,  1],
    [11, 20,  2]], order=31)

In [53]: P, L, U = A.plu_decompose()

In [54]: P
Out[54]: 
GF([[1, 0, 0],
    [0, 0, 1],
    [0, 1, 0]], order=31)

In [55]: L
Out[55]: 
GF([[ 1,  0,  0],
    [ 9,  1,  0],
    [17,  0,  1]], order=31)

In [56]: U
Out[56]: 
GF([[15,  4, 11],
    [ 0, 15, 27],
    [ 0,  0,  0]], order=31)

In [57]: np.array_equal(P @ L @ U, A)
Out[57]: True

See galois.FieldArray.plu_decompose() for more details.


Last update: Apr 03, 2022