Hadamard matrices#
A Hadamard matrix is an \(n\times n\) matrix \(H\) whose entries are either \(+1\) or \(-1\) and whose rows are mutually orthogonal. For example, the matrix \(H_2\) defined by
is a Hadamard matrix. An \(n\times n\) matrix \(H\) whose entries are either \(+1\) or \(-1\) is a Hadamard matrix if and only if:
\(|det(H)|=n^{n/2}\) or
\(H*H^t = n\cdot I_n\), where \(I_n\) is the identity matrix.
In general, the tensor product of an \(m\times m\) Hadamard matrix and an \(n\times n\) Hadamard matrix is an \((mn)\times (mn)\) matrix. In particular, if there is an \(n\times n\) Hadamard matrix then there is a \((2n)\times (2n)\) Hadamard matrix (since one may tensor with \(H_2\)). This particular case is sometimes called the Sylvester construction.
The Hadamard conjecture (possibly due to Paley) states that a Hadamard matrix of order \(n\) exists if and only if \(n= 1, 2\) or \(n\) is a multiple of \(4\).
The module below implements constructions of Hadamard and skew Hadamard matrices for all known orders \(\le 1200\), plus some more greater than \(1200\). It also allows you to pull a Hadamard matrix from the database at [SloaHada].
The following code will test that a construction for all known orders \(\le 4k\)
is implemented. The assertion above can be verified by setting k=300
(note that it will take a long time to run):
sage: from sage.combinat.matrices.hadamard_matrix import (hadamard_matrix,
....: skew_hadamard_matrix, is_hadamard_matrix,
....: is_skew_hadamard_matrix)
sage: k = 20
sage: unknown_hadamard = [668, 716, 892, 1132]
sage: unknown_skew_hadamard = [356, 404, 428, 476, 596, 612, 668, 708, 712, 716,
....: 764, 772, 804, 808, 820, 836, 856, 892, 900, 916,
....: 932, 940, 952, 980, 996, 1004, 1012, 1028, 1036,
....: 1044, 1060, 1076, 1100, 1108, 1132, 1140, 1148,
....: 1156, 1180, 1192, 1196]
sage: for n in range(1, k+1):
....: if 4*n not in unknown_hadamard:
....: H = hadamard_matrix(4*n, check=False)
....: assert is_hadamard_matrix(H)
....: if 4*n not in unknown_skew_hadamard:
....: H = skew_hadamard_matrix(4*n, check=False)
....: assert is_skew_hadamard_matrix(H)
>>> from sage.all import *
>>> from sage.combinat.matrices.hadamard_matrix import (hadamard_matrix,
... skew_hadamard_matrix, is_hadamard_matrix,
... is_skew_hadamard_matrix)
>>> k = Integer(20)
>>> unknown_hadamard = [Integer(668), Integer(716), Integer(892), Integer(1132)]
>>> unknown_skew_hadamard = [Integer(356), Integer(404), Integer(428), Integer(476), Integer(596), Integer(612), Integer(668), Integer(708), Integer(712), Integer(716),
... Integer(764), Integer(772), Integer(804), Integer(808), Integer(820), Integer(836), Integer(856), Integer(892), Integer(900), Integer(916),
... Integer(932), Integer(940), Integer(952), Integer(980), Integer(996), Integer(1004), Integer(1012), Integer(1028), Integer(1036),
... Integer(1044), Integer(1060), Integer(1076), Integer(1100), Integer(1108), Integer(1132), Integer(1140), Integer(1148),
... Integer(1156), Integer(1180), Integer(1192), Integer(1196)]
>>> for n in range(Integer(1), k+Integer(1)):
... if Integer(4)*n not in unknown_hadamard:
... H = hadamard_matrix(Integer(4)*n, check=False)
... assert is_hadamard_matrix(H)
... if Integer(4)*n not in unknown_skew_hadamard:
... H = skew_hadamard_matrix(Integer(4)*n, check=False)
... assert is_skew_hadamard_matrix(H)
AUTHORS:
David Joyner (2009-05-17): initial version
Matteo Cati (2023-03-18): implemented more constructions for Hadamard and skew Hadamard matrices, to cover all known orders up to 1200.
REFERENCES:
- sage.combinat.matrices.hadamard_matrix.GS_skew_hadamard_smallcases(n, existence=False, check=True)[source]#
Data for Williamson-Goethals-Seidel construction of skew Hadamard matrices
Here we keep the data for this construction. Namely, it needs 4 circulant matrices with extra properties, as described in
sage.combinat.matrices.hadamard_matrix.williamson_goethals_seidel_skew_hadamard_matrix()
Matrices are taken from:\(n=36, 52\): [GS70s]
\(n=92\): [Wall71]
\(n=188\): [Djo2008a]
\(n=236\): [FKS2004]
\(n=276\): [Djo2023a]
Additional data is obtained from skew supplementary difference sets contained in
sage.combinat.designs.difference_family.skew_supplementary_difference_set()
, using the construction described in [Djo1992a].INPUT:
n
– integer; the order of the matrixexistence
– boolean (default:True
); ifTrue
, only check that we can do the constructioncheck
– boolean (default:False
): ifTrue
, check the result
- sage.combinat.matrices.hadamard_matrix.RSHCD_324(e)[source]#
Return a size 324x324 Regular Symmetric Hadamard Matrix with Constant Diagonal.
We build the matrix \(M\) for the case \(n=324\), \(\epsilon=1\) directly from
JankoKharaghaniTonchevGraph
and for the case \(\epsilon=-1\) from the “twist” \(M'\) of \(M\), using Lemma 11 in [HX2010]. Namely, it turns out that the matrix\[\begin{split}M'=\begin{pmatrix} M_{12} & M_{11}\\ M_{11}^\top & M_{21} \end{pmatrix}, \quad\text{where}\quad M=\begin{pmatrix} M_{11} & M_{12}\\ M_{21} & M_{22} \end{pmatrix},\end{split}\]and the \(M_{ij}\) are 162x162-blocks, also RSHCD, its diagonal blocks having zero row sums, as needed by [loc.cit.]. Interestingly, the corresponding \((324,152,70,72)\)-strongly regular graph has a vertex-transitive automorphism group of order 2592, twice the order of the (intransitive) automorphism group of the graph corresponding to \(M\). Cf. [CP2016].
INPUT:
e
– \(-1\) or \(+1\); the value of \(\epsilon\)
REFERENCE:
- sage.combinat.matrices.hadamard_matrix.amicable_hadamard_matrices(n, existence=False, check=True)[source]#
Construct amicable Hadamard matrices of order
n
using the available methods.INPUT:
n
– positive integer; the order of the amicable Hadamard matricesexistence
– boolean (default:False
); ifTrue
, only return whether amicable Hadamard matrices of order \(n\) can be constructedcheck
– boolean (default:True
); ifTrue
, check that the matrices are amicable Hadamard matrices before returning them
OUTPUT:
If
existence
is true, the function returns a boolean representing whether amicable Hadamard matrices of order \(n\) can be constructed. Ifexistence
is false, returns two amicable Hadamard matrices, or raises an error if the matrices cannot be constructed.EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import amicable_hadamard_matrices sage: amicable_hadamard_matrices(2) ( [ 1 1] [ 1 1] [-1 1], [ 1 -1] )
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import amicable_hadamard_matrices >>> amicable_hadamard_matrices(Integer(2)) ( [ 1 1] [ 1 1] [-1 1], [ 1 -1] )
If
existence
is true, the function returns a boolean:sage: amicable_hadamard_matrices(16, existence=True) False
>>> from sage.all import * >>> amicable_hadamard_matrices(Integer(16), existence=True) False
- sage.combinat.matrices.hadamard_matrix.amicable_hadamard_matrices_wallis(n, check=True)[source]#
Construct amicable Hadamard matrices of order \(n = q + 1\) where \(q\) is a prime power.
If \(q\) is a prime power \(\equiv 3 \mod 4\), then amicable Hadamard matrices of order \(q+1\) can be constructed as described in [Wal1970b].
INPUT:
n
– integer; the order of the matrices to be constructedcheck
– boolean (default:True
); ifTrue
, check that the resulting matrices are amicable Hadamard before returning them
OUTPUT:
The function returns two amicable Hadamard matrices, or raises an error if such matrices cannot be created using this construction.
EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import amicable_hadamard_matrices_wallis sage: M, N = amicable_hadamard_matrices_wallis(28)
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import amicable_hadamard_matrices_wallis >>> M, N = amicable_hadamard_matrices_wallis(Integer(28))
- sage.combinat.matrices.hadamard_matrix.are_amicable_hadamard_matrices(M, N, verbose=False)[source]#
Check if
M
andN
are amicable Hadamard matrices.Two matrices \(M\) and \(N\) of order \(n\) are called amicable if they satisfy the following conditions (see [Seb2017]):
\(M\) is a skew Hadamard matrix
\(N\) is a symmetric Hadamard matrix
\(MN^T = NM^T\)
INPUT:
M
– a square matrixN
– a square matrixverbose
– boolean (defaultFalse
); whether to be verbose when the matrices are not amicable Hadamard matrices
EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import are_amicable_hadamard_matrices sage: M = matrix([[1, 1], [-1, 1]]) sage: N = matrix([[1, 1], [1, -1]]) sage: are_amicable_hadamard_matrices(M, N) True
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import are_amicable_hadamard_matrices >>> M = matrix([[Integer(1), Integer(1)], [-Integer(1), Integer(1)]]) >>> N = matrix([[Integer(1), Integer(1)], [Integer(1), -Integer(1)]]) >>> are_amicable_hadamard_matrices(M, N) True
If
verbose
is true, the function will be verbose when returning False:sage: N = matrix([[1, 1], [1, 1]]) sage: are_amicable_hadamard_matrices(M, N, verbose=True) The second matrix is not Hadamard False
>>> from sage.all import * >>> N = matrix([[Integer(1), Integer(1)], [Integer(1), Integer(1)]]) >>> are_amicable_hadamard_matrices(M, N, verbose=True) The second matrix is not Hadamard False
- sage.combinat.matrices.hadamard_matrix.construction_four_symbol_delta_code_I(X, Y, Z, W)[source]#
Construct 4-symbol \(\delta\) code of length \(2n+1\).
The 4-symbol \(\delta\) code is constructed from sequences \(X, Y, Z, W\) of length \(n+1\), \(n+1\), \(n\), \(n\) satisfying for all \(s > 0\):
\[N_X(s) + N_Y(s) + N_Z(s) + N_W(s) = 0\]where \(N_A(s)\) is the nonperiodic correlation function:
\[N_A(s) = \sum_{i=1}^{n-s}a_ia_{i+s}\]The construction (detailed in [Tur1974]) is as follows:
\[\begin{split}\begin{aligned} T_1 &= X;Z \\ T_2 &= X;-Z \\ T_3 &= Y;W \\ T_4 &= Y;-W \end{aligned}\end{split}\]INPUT:
X
– list; the first sequence (length \(n+1\))Y
– list; the second sequence (length \(n+1\))Z
– list; the third sequence (length \(n\))W
– list; the fourth sequence (length \(n\))
- OUTPUT:
A tuple containing the 4-symbol \(\delta\) code of length \(2n+1\).
EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import construction_four_symbol_delta_code_I sage: construction_four_symbol_delta_code_I([1, 1], [1, -1], [1], [1]) ([1, 1, 1], [1, 1, -1], [1, -1, 1], [1, -1, -1])
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import construction_four_symbol_delta_code_I >>> construction_four_symbol_delta_code_I([Integer(1), Integer(1)], [Integer(1), -Integer(1)], [Integer(1)], [Integer(1)]) ([1, 1, 1], [1, 1, -1], [1, -1, 1], [1, -1, -1])
- sage.combinat.matrices.hadamard_matrix.construction_four_symbol_delta_code_II(X, Y, Z, W)[source]#
Construct 4-symbol \(\delta\) code of length \(4n+3\).
The 4-symbol \(\delta\) code is constructed from sequences \(X, Y, Z, W\) of length \(n+1\), \(n+1\), \(n\), \(n\) satisfying for all \(s > 0\):
\[N_X(s) + N_Y(s) + N_Z(s) + N_W(s) = 0\]where \(N_A(s)\) is the nonperiodic correlation function:
\[N_A(s) = \sum_{i=1}^{n-s}a_ia_{i+s}\]The construction (detailed in [Tur1974]) is as follows (writing \(A/B\) to mean \(A\) alternated with \(B\)):
\[\begin{split}\begin{aligned} T_1 &= X/Z;Y/W;1 \\ T_2 &= X/Z;Y/-W;-1 \\ T_3 &= X/Z;-Y/-W;1 \\ T_4 &= X/Z;-Y/W;-1 \end{aligned}\end{split}\]INPUT:
X
– list; the first sequence (length \(n+1\))Y
– list; the second sequence (length \(n+1\))Z
– list; the third sequence (length \(n\))W
– list; the fourth sequence (length \(n\))
- OUTPUT:
A tuple containing the four 4-symbol \(\delta\) code of length \(4n+3\).
EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import construction_four_symbol_delta_code_II sage: construction_four_symbol_delta_code_II([1, 1], [1, -1], [1], [1]) ([1, 1, 1, 1, 1, -1, 1], [1, 1, 1, 1, -1, -1, -1], [1, 1, 1, -1, -1, 1, 1], [1, 1, 1, -1, 1, 1, -1])
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import construction_four_symbol_delta_code_II >>> construction_four_symbol_delta_code_II([Integer(1), Integer(1)], [Integer(1), -Integer(1)], [Integer(1)], [Integer(1)]) ([1, 1, 1, 1, 1, -1, 1], [1, 1, 1, 1, -1, -1, -1], [1, 1, 1, -1, -1, 1, 1], [1, 1, 1, -1, 1, 1, -1])
- sage.combinat.matrices.hadamard_matrix.four_symbol_delta_code_smallcases(n, existence=False)[source]#
Return the 4-symobl \(\delta\) code of length
n
if available.The 4-symbol \(\delta\) codes are constructed using
construction_four_symbol_delta_code_I()
orconstruction_four_symbol_delta_code_II()
. The base sequences used are taken from [Tur1974].INPUT:
n
– integer; the length of the desired 4-symbol \(\delta\) codeexistence
– boolean (default:False
); ifTrue
, only check if the sequences are available
EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import four_symbol_delta_code_smallcases sage: four_symbol_delta_code_smallcases(3) ([1, -1, 1], [1, -1, -1], [1, 1, 1], [1, 1, -1]) sage: four_symbol_delta_code_smallcases(3, existence=True) True
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import four_symbol_delta_code_smallcases >>> four_symbol_delta_code_smallcases(Integer(3)) ([1, -1, 1], [1, -1, -1], [1, 1, 1], [1, 1, -1]) >>> four_symbol_delta_code_smallcases(Integer(3), existence=True) True
- sage.combinat.matrices.hadamard_matrix.hadamard_matrix(n, existence=False, check=True, construction_name=False)[source]#
This function is available as hadamard_matrix(…) and matrix.hadamard(…).
Tries to construct a Hadamard matrix using the available methods.
Currently all orders \(\le 1200\) for which a construction is known are implemented. For \(n > 1200\), only some orders are available.
INPUT:
n
– integer; dimension of the matrixexistence
– boolean (default:False
); whether to build the matrix or merely query if a construction is available in Sage. When set toTrue
, the function returns:True
– meaning that Sage knows how to build the matrixUnknown
– meaning that Sage does not know how to build the matrix, although the matrix may exist (seesage.misc.unknown
).False
– meaning that the matrix does not exist.
check
– boolean (default:True
); whether to check that output is correct before returning it. As this is expected to be useless (but we are cautious guys), you may want to disable it whenever you want speed.construction_name
– boolean (default:False
); if it isTrue
,existence
isTrue
, and a matrix exists, output the construction name. It has no effect ifexistence
is set toFalse
.
EXAMPLES:
sage: hadamard_matrix(12).det() 2985984 sage: 12^6 2985984 sage: hadamard_matrix(1) [1] sage: hadamard_matrix(2) [ 1 1] [ 1 -1] sage: hadamard_matrix(8) # random [ 1 1 1 1 1 1 1 1] [ 1 -1 1 -1 1 -1 1 -1] [ 1 1 -1 -1 1 1 -1 -1] [ 1 -1 -1 1 1 -1 -1 1] [ 1 1 1 1 -1 -1 -1 -1] [ 1 -1 1 -1 -1 1 -1 1] [ 1 1 -1 -1 -1 -1 1 1] [ 1 -1 -1 1 -1 1 1 -1] sage: hadamard_matrix(8).det() == 8^4 True
>>> from sage.all import * >>> hadamard_matrix(Integer(12)).det() 2985984 >>> Integer(12)**Integer(6) 2985984 >>> hadamard_matrix(Integer(1)) [1] >>> hadamard_matrix(Integer(2)) [ 1 1] [ 1 -1] >>> hadamard_matrix(Integer(8)) # random [ 1 1 1 1 1 1 1 1] [ 1 -1 1 -1 1 -1 1 -1] [ 1 1 -1 -1 1 1 -1 -1] [ 1 -1 -1 1 1 -1 -1 1] [ 1 1 1 1 -1 -1 -1 -1] [ 1 -1 1 -1 -1 1 -1 1] [ 1 1 -1 -1 -1 -1 1 1] [ 1 -1 -1 1 -1 1 1 -1] >>> hadamard_matrix(Integer(8)).det() == Integer(8)**Integer(4) True
We note that
hadamard_matrix()
returns a normalised Hadamard matrix (the entries in the first row and column are all +1)sage: hadamard_matrix(12) # random [ 1 1| 1 1| 1 1| 1 1| 1 1| 1 1] [ 1 -1|-1 1|-1 1|-1 1|-1 1|-1 1] [-----+-----+-----+-----+-----+-----] [ 1 -1| 1 -1| 1 1|-1 -1|-1 -1| 1 1] [ 1 1|-1 -1| 1 -1|-1 1|-1 1| 1 -1] [-----+-----+-----+-----+-----+-----] [ 1 -1| 1 1| 1 -1| 1 1|-1 -1|-1 -1] [ 1 1| 1 -1|-1 -1| 1 -1|-1 1|-1 1] [-----+-----+-----+-----+-----+-----] [ 1 -1|-1 -1| 1 1| 1 -1| 1 1|-1 -1] [ 1 1|-1 1| 1 -1|-1 -1| 1 -1|-1 1] [-----+-----+-----+-----+-----+-----] [ 1 -1|-1 -1|-1 -1| 1 1| 1 -1| 1 1] [ 1 1|-1 1|-1 1| 1 -1|-1 -1| 1 -1] [-----+-----+-----+-----+-----+-----] [ 1 -1| 1 1|-1 -1|-1 -1| 1 1| 1 -1] [ 1 1| 1 -1|-1 1|-1 1| 1 -1|-1 -1]
>>> from sage.all import * >>> hadamard_matrix(Integer(12)) # random [ 1 1| 1 1| 1 1| 1 1| 1 1| 1 1] [ 1 -1|-1 1|-1 1|-1 1|-1 1|-1 1] [-----+-----+-----+-----+-----+-----] [ 1 -1| 1 -1| 1 1|-1 -1|-1 -1| 1 1] [ 1 1|-1 -1| 1 -1|-1 1|-1 1| 1 -1] [-----+-----+-----+-----+-----+-----] [ 1 -1| 1 1| 1 -1| 1 1|-1 -1|-1 -1] [ 1 1| 1 -1|-1 -1| 1 -1|-1 1|-1 1] [-----+-----+-----+-----+-----+-----] [ 1 -1|-1 -1| 1 1| 1 -1| 1 1|-1 -1] [ 1 1|-1 1| 1 -1|-1 -1| 1 -1|-1 1] [-----+-----+-----+-----+-----+-----] [ 1 -1|-1 -1|-1 -1| 1 1| 1 -1| 1 1] [ 1 1|-1 1|-1 1| 1 -1|-1 -1| 1 -1] [-----+-----+-----+-----+-----+-----] [ 1 -1| 1 1|-1 -1|-1 -1| 1 1| 1 -1] [ 1 1| 1 -1|-1 1|-1 1| 1 -1|-1 -1]
To find how the matrix is obtained, use
construction_name
sage: matrix.hadamard(476, existence=True, construction_name=True) 'cooper-wallis small cases: 476'
>>> from sage.all import * >>> matrix.hadamard(Integer(476), existence=True, construction_name=True) 'cooper-wallis small cases: 476'
- sage.combinat.matrices.hadamard_matrix.hadamard_matrix_156()[source]#
Construct a Hadamard matrix of order 156.
The matrix is created using the construction detailed in [BH1965]. This uses four circulant matrices of size \(13\times 13\), which are composed into a \(156\times 156\) block matrix.
- sage.combinat.matrices.hadamard_matrix.hadamard_matrix_cooper_wallis_construction(x1, x2, x3, x4, A, B, C, D, check=True)[source]#
Create a Hadamard matrix using the contruction detailed in [CW1972].
Given four circulant matrices \(X_1\), X_2, X_3, X_4` of order \(n\) with entries (0, 1, -1) such that the entrywise product of two distinct matrices is always equal to \(0\) and that \(\sum_{i=1}^{4}X_iX_i^\top = nI_n\) holds, and four matrices \(A, B, C, D\) of order \(m\) with elements (1, -1) such that \(MN^\top = NM^\top\) for all distinct \(M\), \(N\) and \(AA^\top + BB^\top + CC^\top + DD^\top = 4mI_n\) holds, we construct a Hadamard matrix of order \(4nm\).
INPUT:
x1
– list or vector; the first row of the circulant matrix \(X_1\)x2
– list or vector; the first row of the circulant matrix \(X_2\)x3
– list or vector; the first row of the circulant matrix \(X_3\)x4
– list or vector; the first row of the circulant matrix \(X_4\)A
– the matrix described aboveB
– the matrix described aboveC
– the matrix described aboveD
– the matrix described abovecheck
– boolean (default:True
); ifTrue
, check that the resulting matrix is Hadamard before returing it.
EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import hadamard_matrix_cooper_wallis_construction sage: from sage.combinat.t_sequences import T_sequences_smallcases sage: seqs = T_sequences_smallcases(19) sage: hadamard_matrix_cooper_wallis_construction(seqs[0], seqs[1], seqs[2], seqs[3], matrix([1]), matrix([1]), matrix([1]), matrix([1])) 76 x 76 dense matrix over Integer Ring...
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import hadamard_matrix_cooper_wallis_construction >>> from sage.combinat.t_sequences import T_sequences_smallcases >>> seqs = T_sequences_smallcases(Integer(19)) >>> hadamard_matrix_cooper_wallis_construction(seqs[Integer(0)], seqs[Integer(1)], seqs[Integer(2)], seqs[Integer(3)], matrix([Integer(1)]), matrix([Integer(1)]), matrix([Integer(1)]), matrix([Integer(1)])) 76 x 76 dense matrix over Integer Ring...
- sage.combinat.matrices.hadamard_matrix.hadamard_matrix_cooper_wallis_smallcases(n, check=True, existence=False)[source]#
Construct Hadamard matrices using the Cooper-Wallis construction for some small values of \(n\).
This function calls the function
hadamard_matrix_cooper_wallis_construction()
with the appropriate arguments. It constructs the matrices \(X_1\), \(X_2\), \(X_3\), \(X_4\) using either T-matrices or the T-sequences fromsage.combinat.t_sequences.T_sequences_smallcases()
. The matrices \(A\), \(B\), \(C\), \(D\) are taken fromwilliamson_type_quadruples_smallcases()
.Data for T-matrices of order 67 is taken from [Saw1985].
INPUT:
n
– integer; the order of the matrix to be constructedcheck
– boolean (default:True
); ifTrue
, check that the matrix is a Hadamard matrix before returningexistence
– boolean (default:False
); ifTrue
, only check if the matrix exists.
OUTPUT:
If
existence=False
, returns the Hadamard matrix of order \(n\). It raises an error if no data is available to construct the matrix of the given order. Ifexistence=True
, returns a boolean representing whether the matrix can be constructed or not.EXAMPLES:
By default The function returns the Hadamard matrix
sage: from sage.combinat.matrices.hadamard_matrix import hadamard_matrix_cooper_wallis_smallcases sage: hadamard_matrix_cooper_wallis_smallcases(28) 28 x 28 dense matrix over Integer Ring...
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import hadamard_matrix_cooper_wallis_smallcases >>> hadamard_matrix_cooper_wallis_smallcases(Integer(28)) 28 x 28 dense matrix over Integer Ring...
If
existence
is set to True, the function returns a booleansage: hadamard_matrix_cooper_wallis_smallcases(20, existence=True) True
>>> from sage.all import * >>> hadamard_matrix_cooper_wallis_smallcases(Integer(20), existence=True) True
- sage.combinat.matrices.hadamard_matrix.hadamard_matrix_from_sds(n, existence=False, check=True)[source]#
Construction of Hadamard matrices from supplementary difference sets.
Given four SDS with parameters \(4-\{n; n_1, n_2, n_3, n_4; \lambda\}\) with \(n_1 + n_2 + n_3 + n_4 = n+\lambda\) we can construct four (-1,+1) sequences \(a_i = (a_{i,0},...,a_{i,n-1})\) where \(a_{i,j} = -1\) iff \(j \in S_i\). These will be the fist rows of four circulant matrices \(A_1, A_2, A_3, A_4\) which, when plugged into the Goethals-Seidel array, create an Hadamard matrix of order \(4n\) (see [Djo1994b]).
The supplementary difference sets are taken from
sage.combinat.designs.difference_family.supplementary_difference_set()
.INPUT:
n
– integer; the order of the matrix to be constructedcheck
– boolean (default:True
); ifTrue
, check that the matrix is a Hadamard before returningexistence
– boolean (default:False
); ifTrue
, only check if the matrix exists
OUTPUT:
If
existence=False
, returns the Hadamard matrix of order \(n\). It raises an error if no data is available to construct the matrix of the given order, or if \(n\) is not a multiple of \(4\). Ifexistence=True
, returns a boolean representing whether the matrix can be constructed or not.EXAMPLES:
By default The function returns the Hadamard matrix
sage: from sage.combinat.matrices.hadamard_matrix import hadamard_matrix_from_sds sage: hadamard_matrix_from_sds(148) 148 x 148 dense matrix over Integer Ring...
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import hadamard_matrix_from_sds >>> hadamard_matrix_from_sds(Integer(148)) 148 x 148 dense matrix over Integer Ring...
If
existence
is set to True, the function returns a booleansage: hadamard_matrix_from_sds(764, existence=True) True
>>> from sage.all import * >>> hadamard_matrix_from_sds(Integer(764), existence=True) True
- sage.combinat.matrices.hadamard_matrix.hadamard_matrix_miyamoto_construction(n, existence=False, check=True)[source]#
Construct Hadamard matrix using the Miyamoto construction.
If \(q = n/4\) is a prime power, and there exists a Hadamard matrix of order \(q-1\), then a Hadamard matrix of order \(n\) can be constructed (see [Miy1991]).
INPUT:
n
– integer; the order of the matrix to be constructedcheck
– boolean (default:True
); ifTrue
, check that the matrix is a Hadamard before returningexistence
– boolean (default:False
); ifTrue
, only check if the matrix exists
OUTPUT:
If
existence=False
, returns the Hadamard matrix of order \(n\). It raises an error if no data is available to construct the matrix of the given order, or if \(n\) does not satisfies the constraints. Ifexistence=True
, returns a boolean representing whether the matrix can be constructed or not.EXAMPLES:
By default the function returns the Hadamard matrix
sage: from sage.combinat.matrices.hadamard_matrix import hadamard_matrix_miyamoto_construction sage: hadamard_matrix_miyamoto_construction(20) 20 x 20 dense matrix over Integer Ring...
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import hadamard_matrix_miyamoto_construction >>> hadamard_matrix_miyamoto_construction(Integer(20)) 20 x 20 dense matrix over Integer Ring...
If
existence
is set to True, the function returns a booleansage: hadamard_matrix_miyamoto_construction(36, existence=True) True
>>> from sage.all import * >>> hadamard_matrix_miyamoto_construction(Integer(36), existence=True) True
- sage.combinat.matrices.hadamard_matrix.hadamard_matrix_paleyI(n, normalize=True)[source]#
Implement the Paley type I construction.
The Paley type I case corresponds to the case \(p=n-1 \cong 3 \mod{4}\) for a prime power \(p\) (see [Hora]).
INPUT:
n
– the matrix sizenormalize
– boolean (default:True
); whether to normalize the result
EXAMPLES:
We note that this method by default returns a normalised Hadamard matrix
sage: from sage.combinat.matrices.hadamard_matrix import hadamard_matrix_paleyI sage: hadamard_matrix_paleyI(4) [ 1 1 1 1] [ 1 -1 1 -1] [ 1 -1 -1 1] [ 1 1 -1 -1]
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import hadamard_matrix_paleyI >>> hadamard_matrix_paleyI(Integer(4)) [ 1 1 1 1] [ 1 -1 1 -1] [ 1 -1 -1 1] [ 1 1 -1 -1]
Otherwise, it returns a skew Hadamard matrix \(H\), i.e. \(H=S+I\), with \(S=-S^\top\)
sage: M = hadamard_matrix_paleyI(4, normalize=False); M [ 1 1 1 1] [-1 1 1 -1] [-1 -1 1 1] [-1 1 -1 1] sage: S = M - identity_matrix(4); -S == S.T True
>>> from sage.all import * >>> M = hadamard_matrix_paleyI(Integer(4), normalize=False); M [ 1 1 1 1] [-1 1 1 -1] [-1 -1 1 1] [-1 1 -1 1] >>> S = M - identity_matrix(Integer(4)); -S == S.T True
- sage.combinat.matrices.hadamard_matrix.hadamard_matrix_paleyII(n)[source]#
Implement the Paley type II construction.
The Paley type II case corresponds to the case \(p=n/2-1 \cong 1 \mod{4}\) for a prime power \(p\) (see [Hora]).
EXAMPLES:
sage: sage.combinat.matrices.hadamard_matrix.hadamard_matrix_paleyII(12).det() 2985984 sage: 12^6 2985984
>>> from sage.all import * >>> sage.combinat.matrices.hadamard_matrix.hadamard_matrix_paleyII(Integer(12)).det() 2985984 >>> Integer(12)**Integer(6) 2985984
We note that the method returns a normalised Hadamard matrix
sage: sage.combinat.matrices.hadamard_matrix.hadamard_matrix_paleyII(12) [ 1 1| 1 1| 1 1| 1 1| 1 1| 1 1] [ 1 -1|-1 1|-1 1|-1 1|-1 1|-1 1] [-----+-----+-----+-----+-----+-----] [ 1 -1| 1 -1| 1 1|-1 -1|-1 -1| 1 1] [ 1 1|-1 -1| 1 -1|-1 1|-1 1| 1 -1] [-----+-----+-----+-----+-----+-----] [ 1 -1| 1 1| 1 -1| 1 1|-1 -1|-1 -1] [ 1 1| 1 -1|-1 -1| 1 -1|-1 1|-1 1] [-----+-----+-----+-----+-----+-----] [ 1 -1|-1 -1| 1 1| 1 -1| 1 1|-1 -1] [ 1 1|-1 1| 1 -1|-1 -1| 1 -1|-1 1] [-----+-----+-----+-----+-----+-----] [ 1 -1|-1 -1|-1 -1| 1 1| 1 -1| 1 1] [ 1 1|-1 1|-1 1| 1 -1|-1 -1| 1 -1] [-----+-----+-----+-----+-----+-----] [ 1 -1| 1 1|-1 -1|-1 -1| 1 1| 1 -1] [ 1 1| 1 -1|-1 1|-1 1| 1 -1|-1 -1]
>>> from sage.all import * >>> sage.combinat.matrices.hadamard_matrix.hadamard_matrix_paleyII(Integer(12)) [ 1 1| 1 1| 1 1| 1 1| 1 1| 1 1] [ 1 -1|-1 1|-1 1|-1 1|-1 1|-1 1] [-----+-----+-----+-----+-----+-----] [ 1 -1| 1 -1| 1 1|-1 -1|-1 -1| 1 1] [ 1 1|-1 -1| 1 -1|-1 1|-1 1| 1 -1] [-----+-----+-----+-----+-----+-----] [ 1 -1| 1 1| 1 -1| 1 1|-1 -1|-1 -1] [ 1 1| 1 -1|-1 -1| 1 -1|-1 1|-1 1] [-----+-----+-----+-----+-----+-----] [ 1 -1|-1 -1| 1 1| 1 -1| 1 1|-1 -1] [ 1 1|-1 1| 1 -1|-1 -1| 1 -1|-1 1] [-----+-----+-----+-----+-----+-----] [ 1 -1|-1 -1|-1 -1| 1 1| 1 -1| 1 1] [ 1 1|-1 1|-1 1| 1 -1|-1 -1| 1 -1] [-----+-----+-----+-----+-----+-----] [ 1 -1| 1 1|-1 -1|-1 -1| 1 1| 1 -1] [ 1 1| 1 -1|-1 1|-1 1| 1 -1|-1 -1]
- sage.combinat.matrices.hadamard_matrix.hadamard_matrix_spence_construction(n, existence=False, check=True)[source]#
Create a Hadamard matrix of order \(n\) using the Spence construction.
This construction (detailed in [Spe1975]), uses supplementary difference sets implemented in
sage.combinat.designs.difference_family.supplementary_difference_set_from_rel_diff_set()
to create the desired matrix.INPUT:
n
– integer; the order of the matrix to be constructedexistence
– boolean (default:False
); ifTrue
, only check if the matrix existscheck
– bolean (default:True
); ifTrue
, check that the matrix is a Hadamard matrix before returning
OUTPUT:
If
existence=True
, returns a boolean representing whether the Hadamard matrix can be constructed. Otherwise, returns the Hadamard matrix, or raises an error if it cannot be constructed.EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import hadamard_matrix_spence_construction sage: hadamard_matrix_spence_construction(36) 36 x 36 dense matrix over Integer Ring...
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import hadamard_matrix_spence_construction >>> hadamard_matrix_spence_construction(Integer(36)) 36 x 36 dense matrix over Integer Ring...
If
existence
isTrue
, the function returns a booleansage: hadamard_matrix_spence_construction(52, existence=True) True
>>> from sage.all import * >>> hadamard_matrix_spence_construction(Integer(52), existence=True) True
- sage.combinat.matrices.hadamard_matrix.hadamard_matrix_turyn_type(a, b, c, d, e1, e2, e3, e4, check=True)[source]#
Construction of Turyn type Hadamard matrix.
Given \(n\times n\) circulant matrices \(A\), \(B\), \(C\), \(D\) with 1,-1 entries, satisfying \(AA^\top + BB^\top + CC^\top + DD^\top = 4nI\), and a set of Baumert-Hall units of order \(4t\), one can construct a Hadamard matrix of order \(4tn\) as detailed by Turyn in [Tur1974].
INPUT:
a
– 1,-1 list; the 1st row of \(A\)b
– 1,-1 list; the 1st row of \(B\)d
– 1,-1 list; the 1st row of \(C\)c
– 1,-1 list; the 1st row of \(D\)e1
– Matrix; the first Baumert-Hall unite2
– Matrix; the second Baumert-Hall unite3
– Matrix; the third Baumert-Hall unite4
– Matrix; the fourth Baumert-Hall unitcheck
– boolean (default:True
); whether to check that the output is a Hadamard matrix before returning it
EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import hadamard_matrix_turyn_type, _get_baumert_hall_units sage: A, B, C, D = _get_baumert_hall_units(28) sage: hadamard_matrix_turyn_type([1], [1], [1], [1], A, B, C, D) 28 x 28 dense matrix over Integer Ring...
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import hadamard_matrix_turyn_type, _get_baumert_hall_units >>> A, B, C, D = _get_baumert_hall_units(Integer(28)) >>> hadamard_matrix_turyn_type([Integer(1)], [Integer(1)], [Integer(1)], [Integer(1)], A, B, C, D) 28 x 28 dense matrix over Integer Ring...
- sage.combinat.matrices.hadamard_matrix.hadamard_matrix_williamson_type(a, b, c, d, check=True)[source]#
Construction of Williamson type Hadamard matrix.
Given \(n\times n\) circulant matrices \(A\), \(B\), \(C\), \(D\) with 1,-1 entries, and satisfying \(AA^\top + BB^\top + CC^\top + DD^\top = 4nI\), one can construct a Hadamard matrix of order \(4n\), cf. [Ha83].
INPUT:
a
– (1,-1) list; the 1st row of \(A\)b
– (1,-1) list; the 1st row of \(B\)d
– (1,-1) list; the 1st row of \(C\)c
– (1,-1) list; the 1st row of \(D\)check
– boolean (default:True
); whether to check that the output is a Hadamard matrix before returning it
EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import hadamard_matrix_williamson_type sage: a = [ 1, 1, 1] sage: b = [ 1, -1, -1] sage: c = [ 1, -1, -1] sage: d = [ 1, -1, -1] sage: M = hadamard_matrix_williamson_type(a,b,c,d,check=True)
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import hadamard_matrix_williamson_type >>> a = [ Integer(1), Integer(1), Integer(1)] >>> b = [ Integer(1), -Integer(1), -Integer(1)] >>> c = [ Integer(1), -Integer(1), -Integer(1)] >>> d = [ Integer(1), -Integer(1), -Integer(1)] >>> M = hadamard_matrix_williamson_type(a,b,c,d,check=True)
- sage.combinat.matrices.hadamard_matrix.hadamard_matrix_www(url_file, comments=False)[source]#
Pull file from Sloane’s database and return the corresponding Hadamard matrix as a Sage matrix.
You must input a filename of the form “had.n.xxx.txt” as described on the webpage http://neilsloane.com/hadamard/, where “xxx” could be empty or a number of some characters.
If
comments=True
then the “Automorphism…” line of the had.n.xxx.txt file is printed if it exists. Otherwise nothing is done.EXAMPLES:
sage: hadamard_matrix_www("had.4.txt") # optional - internet [ 1 1 1 1] [ 1 -1 1 -1] [ 1 1 -1 -1] [ 1 -1 -1 1] sage: hadamard_matrix_www("had.16.2.txt",comments=True) # optional - internet Automorphism group has order = 49152 = 2^14 * 3 [ 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1] [ 1 -1 1 -1 1 -1 1 -1 1 -1 1 -1 1 -1 1 -1] [ 1 1 -1 -1 1 1 -1 -1 1 1 -1 -1 1 1 -1 -1] [ 1 -1 -1 1 1 -1 -1 1 1 -1 -1 1 1 -1 -1 1] [ 1 1 1 1 -1 -1 -1 -1 1 1 1 1 -1 -1 -1 -1] [ 1 -1 1 -1 -1 1 -1 1 1 -1 1 -1 -1 1 -1 1] [ 1 1 -1 -1 -1 -1 1 1 1 1 -1 -1 -1 -1 1 1] [ 1 -1 -1 1 -1 1 1 -1 1 -1 -1 1 -1 1 1 -1] [ 1 1 1 1 1 1 1 1 -1 -1 -1 -1 -1 -1 -1 -1] [ 1 1 1 1 -1 -1 -1 -1 -1 -1 -1 -1 1 1 1 1] [ 1 1 -1 -1 1 -1 1 -1 -1 -1 1 1 -1 1 -1 1] [ 1 1 -1 -1 -1 1 -1 1 -1 -1 1 1 1 -1 1 -1] [ 1 -1 1 -1 1 -1 -1 1 -1 1 -1 1 -1 1 1 -1] [ 1 -1 1 -1 -1 1 1 -1 -1 1 -1 1 1 -1 -1 1] [ 1 -1 -1 1 1 1 -1 -1 -1 1 1 -1 -1 -1 1 1] [ 1 -1 -1 1 -1 -1 1 1 -1 1 1 -1 1 1 -1 -1]
>>> from sage.all import * >>> hadamard_matrix_www("had.4.txt") # optional - internet [ 1 1 1 1] [ 1 -1 1 -1] [ 1 1 -1 -1] [ 1 -1 -1 1] >>> hadamard_matrix_www("had.16.2.txt",comments=True) # optional - internet Automorphism group has order = 49152 = 2^14 * 3 [ 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1] [ 1 -1 1 -1 1 -1 1 -1 1 -1 1 -1 1 -1 1 -1] [ 1 1 -1 -1 1 1 -1 -1 1 1 -1 -1 1 1 -1 -1] [ 1 -1 -1 1 1 -1 -1 1 1 -1 -1 1 1 -1 -1 1] [ 1 1 1 1 -1 -1 -1 -1 1 1 1 1 -1 -1 -1 -1] [ 1 -1 1 -1 -1 1 -1 1 1 -1 1 -1 -1 1 -1 1] [ 1 1 -1 -1 -1 -1 1 1 1 1 -1 -1 -1 -1 1 1] [ 1 -1 -1 1 -1 1 1 -1 1 -1 -1 1 -1 1 1 -1] [ 1 1 1 1 1 1 1 1 -1 -1 -1 -1 -1 -1 -1 -1] [ 1 1 1 1 -1 -1 -1 -1 -1 -1 -1 -1 1 1 1 1] [ 1 1 -1 -1 1 -1 1 -1 -1 -1 1 1 -1 1 -1 1] [ 1 1 -1 -1 -1 1 -1 1 -1 -1 1 1 1 -1 1 -1] [ 1 -1 1 -1 1 -1 -1 1 -1 1 -1 1 -1 1 1 -1] [ 1 -1 1 -1 -1 1 1 -1 -1 1 -1 1 1 -1 -1 1] [ 1 -1 -1 1 1 1 -1 -1 -1 1 1 -1 -1 -1 1 1] [ 1 -1 -1 1 -1 -1 1 1 -1 1 1 -1 1 1 -1 -1]
- sage.combinat.matrices.hadamard_matrix.is_hadamard_matrix(M, normalized=False, skew=False, verbose=False)[source]#
Test if
M
is a Hadamard matrix.INPUT:
M
– a matrixnormalized
– boolean (default:False
); whether to test ifM
is a normalized Hadamard matrix, i.e. has its first row/column filled with +1skew
– boolean (default:False
); whether to test ifM
is a skew Hadamard matrix, i.e. \(M=S+I\) for \(-S=S^\top\), and \(I\) the identity matrixverbose
– boolean (default:False
); whether to be verbose when the matrix is not Hadamard
EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import is_hadamard_matrix sage: h = matrix.hadamard(12) sage: is_hadamard_matrix(h) True sage: from sage.combinat.matrices.hadamard_matrix import skew_hadamard_matrix sage: h = skew_hadamard_matrix(12) sage: is_hadamard_matrix(h, skew=True) True sage: h = matrix.hadamard(12) sage: h[0,0] = 2 sage: is_hadamard_matrix(h,verbose=True) The matrix does not only contain +1 and -1 entries, e.g. 2 False sage: h = matrix.hadamard(12) sage: for i in range(12): ....: h[i,2] = -h[i,2] sage: is_hadamard_matrix(h,verbose=True,normalized=True) The matrix is not normalized False
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import is_hadamard_matrix >>> h = matrix.hadamard(Integer(12)) >>> is_hadamard_matrix(h) True >>> from sage.combinat.matrices.hadamard_matrix import skew_hadamard_matrix >>> h = skew_hadamard_matrix(Integer(12)) >>> is_hadamard_matrix(h, skew=True) True >>> h = matrix.hadamard(Integer(12)) >>> h[Integer(0),Integer(0)] = Integer(2) >>> is_hadamard_matrix(h,verbose=True) The matrix does not only contain +1 and -1 entries, e.g. 2 False >>> h = matrix.hadamard(Integer(12)) >>> for i in range(Integer(12)): ... h[i,Integer(2)] = -h[i,Integer(2)] >>> is_hadamard_matrix(h,verbose=True,normalized=True) The matrix is not normalized False
- sage.combinat.matrices.hadamard_matrix.is_skew_hadamard_matrix(M, normalized=False, verbose=False)[source]#
Test if
M
is a skew Hadamard matrix.this is a wrapper around the function
is_hadamard_matrix()
INPUT:
M
– a matrixnormalized
– boolean (default:False
); whether to test ifM
is a skew-normalized Hadamard matrix, i.e. has its first row filled with +1verbose
– boolean (default:False
); whether to be verbose when the matrix is not skew Hadamard
EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import is_skew_hadamard_matrix, skew_hadamard_matrix sage: h = matrix.hadamard(12) sage: is_skew_hadamard_matrix(h, verbose=True) The matrix is not skew False sage: h = skew_hadamard_matrix(12) sage: is_skew_hadamard_matrix(h) True sage: from sage.combinat.matrices.hadamard_matrix import normalise_hadamard sage: h = normalise_hadamard(skew_hadamard_matrix(12), skew=True) sage: is_skew_hadamard_matrix(h, verbose=True, normalized=True) True
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import is_skew_hadamard_matrix, skew_hadamard_matrix >>> h = matrix.hadamard(Integer(12)) >>> is_skew_hadamard_matrix(h, verbose=True) The matrix is not skew False >>> h = skew_hadamard_matrix(Integer(12)) >>> is_skew_hadamard_matrix(h) True >>> from sage.combinat.matrices.hadamard_matrix import normalise_hadamard >>> h = normalise_hadamard(skew_hadamard_matrix(Integer(12)), skew=True) >>> is_skew_hadamard_matrix(h, verbose=True, normalized=True) True
- sage.combinat.matrices.hadamard_matrix.normalise_hadamard(H, skew=False)[source]#
Return the normalised Hadamard matrix corresponding to
H
.The normalised Hadamard matrix corresponding to a Hadamard matrix \(H\) is a matrix whose every entry in the first row and column is +1.
If
skew
is True, the matrix returned will be skew-normal: a skew Hadamard matrix with first row of all \(+1\).EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import normalise_hadamard, is_hadamard_matrix, skew_hadamard_matrix sage: H = normalise_hadamard(hadamard_matrix(4)) sage: H == hadamard_matrix(4) True sage: H = normalise_hadamard(skew_hadamard_matrix(20, skew_normalize=False), skew=True) sage: is_hadamard_matrix(H, skew=True, normalized=True) True
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import normalise_hadamard, is_hadamard_matrix, skew_hadamard_matrix >>> H = normalise_hadamard(hadamard_matrix(Integer(4))) >>> H == hadamard_matrix(Integer(4)) True >>> H = normalise_hadamard(skew_hadamard_matrix(Integer(20), skew_normalize=False), skew=True) >>> is_hadamard_matrix(H, skew=True, normalized=True) True
If
skew
is True but the Hadamard matrix is not skew, the matrix returned will not be normalized:sage: H = normalise_hadamard(hadamard_matrix(92), skew=True) sage: is_hadamard_matrix(H, normalized=True) False
>>> from sage.all import * >>> H = normalise_hadamard(hadamard_matrix(Integer(92)), skew=True) >>> is_hadamard_matrix(H, normalized=True) False
- sage.combinat.matrices.hadamard_matrix.regular_symmetric_hadamard_matrix_with_constant_diagonal(n, e, existence=False)[source]#
Return a Regular Symmetric Hadamard Matrix with Constant Diagonal.
A Hadamard matrix is said to be regular if its rows all sum to the same value.
For \(\epsilon\in\{-1,+1\}\), we say that \(M\) is a \((n,\epsilon)-RSHCD\) if \(M\) is a regular symmetric Hadamard matrix with constant diagonal \(\delta\in\{-1,+1\}\) and row sums all equal to \(\delta \epsilon \sqrt(n)\). For more information, see [HX2010] or 10.5.1 in [BH2012]. For the case \(n=324\), see
RSHCD_324()
and [CP2016].INPUT:
n
– integer; side of the matrixe
– \(-1\) or \(+1\); the value of \(\epsilon\)
EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import regular_symmetric_hadamard_matrix_with_constant_diagonal sage: regular_symmetric_hadamard_matrix_with_constant_diagonal(4,1) [ 1 1 1 -1] [ 1 1 -1 1] [ 1 -1 1 1] [-1 1 1 1] sage: regular_symmetric_hadamard_matrix_with_constant_diagonal(4,-1) [ 1 -1 -1 -1] [-1 1 -1 -1] [-1 -1 1 -1] [-1 -1 -1 1]
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import regular_symmetric_hadamard_matrix_with_constant_diagonal >>> regular_symmetric_hadamard_matrix_with_constant_diagonal(Integer(4),Integer(1)) [ 1 1 1 -1] [ 1 1 -1 1] [ 1 -1 1 1] [-1 1 1 1] >>> regular_symmetric_hadamard_matrix_with_constant_diagonal(Integer(4),-Integer(1)) [ 1 -1 -1 -1] [-1 1 -1 -1] [-1 -1 1 -1] [-1 -1 -1 1]
Other hardcoded values:
sage: for n,e in [(36,1),(36,-1),(100,1),(100,-1),(196, 1)]: # long time ....: print(repr(regular_symmetric_hadamard_matrix_with_constant_diagonal(n,e))) 36 x 36 dense matrix over Integer Ring 36 x 36 dense matrix over Integer Ring 100 x 100 dense matrix over Integer Ring 100 x 100 dense matrix over Integer Ring 196 x 196 dense matrix over Integer Ring sage: for n,e in [(324,1),(324,-1)]: # not tested - long time, tested in RSHCD_324 ....: print(repr(regular_symmetric_hadamard_matrix_with_constant_diagonal(n,e))) 324 x 324 dense matrix over Integer Ring 324 x 324 dense matrix over Integer Ring
>>> from sage.all import * >>> for n,e in [(Integer(36),Integer(1)),(Integer(36),-Integer(1)),(Integer(100),Integer(1)),(Integer(100),-Integer(1)),(Integer(196), Integer(1))]: # long time ... print(repr(regular_symmetric_hadamard_matrix_with_constant_diagonal(n,e))) 36 x 36 dense matrix over Integer Ring 36 x 36 dense matrix over Integer Ring 100 x 100 dense matrix over Integer Ring 100 x 100 dense matrix over Integer Ring 196 x 196 dense matrix over Integer Ring >>> for n,e in [(Integer(324),Integer(1)),(Integer(324),-Integer(1))]: # not tested - long time, tested in RSHCD_324 ... print(repr(regular_symmetric_hadamard_matrix_with_constant_diagonal(n,e))) 324 x 324 dense matrix over Integer Ring 324 x 324 dense matrix over Integer Ring
From two close prime powers:
sage: regular_symmetric_hadamard_matrix_with_constant_diagonal(64,-1) 64 x 64 dense matrix over Integer Ring (use the '.str()' method to see the entries)
>>> from sage.all import * >>> regular_symmetric_hadamard_matrix_with_constant_diagonal(Integer(64),-Integer(1)) 64 x 64 dense matrix over Integer Ring (use the '.str()' method to see the entries)
From a prime power and a conference matrix:
sage: regular_symmetric_hadamard_matrix_with_constant_diagonal(676,1) # long time 676 x 676 dense matrix over Integer Ring (use the '.str()' method to see the entries)
>>> from sage.all import * >>> regular_symmetric_hadamard_matrix_with_constant_diagonal(Integer(676),Integer(1)) # long time 676 x 676 dense matrix over Integer Ring (use the '.str()' method to see the entries)
Recursive construction:
sage: regular_symmetric_hadamard_matrix_with_constant_diagonal(144,-1) 144 x 144 dense matrix over Integer Ring (use the '.str()' method to see the entries)
>>> from sage.all import * >>> regular_symmetric_hadamard_matrix_with_constant_diagonal(Integer(144),-Integer(1)) 144 x 144 dense matrix over Integer Ring (use the '.str()' method to see the entries)
REFERENCE:
- sage.combinat.matrices.hadamard_matrix.rshcd_from_close_prime_powers(n)[source]#
Return a \((n^2,1)\)-RSHCD when \(n-1\) and \(n+1\) are odd prime powers and \(n=0\pmod{4}\).
The construction implemented here appears in Theorem 4.3 from [GS1970].
Note that the authors of [SWW1972] claim in Corollary 5.12 (page 342) to have proved the same result without the \(n=0\pmod{4}\) restriction with a very similar construction. So far, however, I (Nathann Cohen) have not been able to make it work.
INPUT:
n
– an integer congruent to \(0\pmod{4}\)
EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import rshcd_from_close_prime_powers sage: rshcd_from_close_prime_powers(4) [-1 -1 1 -1 1 -1 -1 1 -1 1 -1 -1 1 -1 1 -1] [-1 -1 -1 1 1 -1 -1 1 -1 -1 1 -1 -1 1 -1 1] [ 1 -1 -1 1 1 1 1 -1 -1 -1 -1 -1 -1 -1 1 -1] [-1 1 1 -1 1 1 -1 -1 -1 -1 -1 1 -1 -1 -1 1] [ 1 1 1 1 -1 -1 -1 -1 -1 -1 1 -1 1 -1 -1 -1] [-1 -1 1 1 -1 -1 1 -1 -1 1 -1 1 -1 1 -1 -1] [-1 -1 1 -1 -1 1 -1 -1 1 -1 1 -1 -1 1 1 -1] [ 1 1 -1 -1 -1 -1 -1 -1 -1 1 -1 -1 -1 1 1 1] [-1 -1 -1 -1 -1 -1 1 -1 -1 -1 1 1 1 -1 1 1] [ 1 -1 -1 -1 -1 1 -1 1 -1 -1 -1 1 1 1 -1 -1] [-1 1 -1 -1 1 -1 1 -1 1 -1 -1 -1 1 1 -1 -1] [-1 -1 -1 1 -1 1 -1 -1 1 1 -1 -1 1 -1 -1 1] [ 1 -1 -1 -1 1 -1 -1 -1 1 1 1 1 -1 -1 -1 -1] [-1 1 -1 -1 -1 1 1 1 -1 1 1 -1 -1 -1 -1 -1] [ 1 -1 1 -1 -1 -1 1 1 1 -1 -1 -1 -1 -1 -1 1] [-1 1 -1 1 -1 -1 -1 1 1 -1 -1 1 -1 -1 1 -1]
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import rshcd_from_close_prime_powers >>> rshcd_from_close_prime_powers(Integer(4)) [-1 -1 1 -1 1 -1 -1 1 -1 1 -1 -1 1 -1 1 -1] [-1 -1 -1 1 1 -1 -1 1 -1 -1 1 -1 -1 1 -1 1] [ 1 -1 -1 1 1 1 1 -1 -1 -1 -1 -1 -1 -1 1 -1] [-1 1 1 -1 1 1 -1 -1 -1 -1 -1 1 -1 -1 -1 1] [ 1 1 1 1 -1 -1 -1 -1 -1 -1 1 -1 1 -1 -1 -1] [-1 -1 1 1 -1 -1 1 -1 -1 1 -1 1 -1 1 -1 -1] [-1 -1 1 -1 -1 1 -1 -1 1 -1 1 -1 -1 1 1 -1] [ 1 1 -1 -1 -1 -1 -1 -1 -1 1 -1 -1 -1 1 1 1] [-1 -1 -1 -1 -1 -1 1 -1 -1 -1 1 1 1 -1 1 1] [ 1 -1 -1 -1 -1 1 -1 1 -1 -1 -1 1 1 1 -1 -1] [-1 1 -1 -1 1 -1 1 -1 1 -1 -1 -1 1 1 -1 -1] [-1 -1 -1 1 -1 1 -1 -1 1 1 -1 -1 1 -1 -1 1] [ 1 -1 -1 -1 1 -1 -1 -1 1 1 1 1 -1 -1 -1 -1] [-1 1 -1 -1 -1 1 1 1 -1 1 1 -1 -1 -1 -1 -1] [ 1 -1 1 -1 -1 -1 1 1 1 -1 -1 -1 -1 -1 -1 1] [-1 1 -1 1 -1 -1 -1 1 1 -1 -1 1 -1 -1 1 -1]
REFERENCE:
- sage.combinat.matrices.hadamard_matrix.rshcd_from_prime_power_and_conference_matrix(n)[source]#
Return a \(((n-1)^2,1)\)-RSHCD if \(n\) is prime power, and symmetric \((n-1)\)-conference matrix exists
The construction implemented here is Theorem 16 (and Corollary 17) from [WW1972].
In [SWW1972] this construction (Theorem 5.15 and Corollary 5.16) is reproduced with a typo. Note that [WW1972] refers to [Sz1969] for the construction, provided by
szekeres_difference_set_pair()
, of complementary difference sets, and the latter has a typo.From a
symmetric_conference_matrix()
, we only need the Seidel adjacency matrix of the underlying strongly regular conference (i.e. Paley type) graph, which we construct directly.INPUT:
n
– an integer
EXAMPLES:
A 36x36 example
sage: from sage.combinat.matrices.hadamard_matrix import rshcd_from_prime_power_and_conference_matrix sage: from sage.combinat.matrices.hadamard_matrix import is_hadamard_matrix sage: H = rshcd_from_prime_power_and_conference_matrix(7); H 36 x 36 dense matrix over Integer Ring (use the '.str()' method to see the entries) sage: H == H.T and is_hadamard_matrix(H) and H.diagonal() == [1]*36 and list(sum(H)) == [6]*36 True
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import rshcd_from_prime_power_and_conference_matrix >>> from sage.combinat.matrices.hadamard_matrix import is_hadamard_matrix >>> H = rshcd_from_prime_power_and_conference_matrix(Integer(7)); H 36 x 36 dense matrix over Integer Ring (use the '.str()' method to see the entries) >>> H == H.T and is_hadamard_matrix(H) and H.diagonal() == [Integer(1)]*Integer(36) and list(sum(H)) == [Integer(6)]*Integer(36) True
Bigger examples, only provided by this construction
sage: H = rshcd_from_prime_power_and_conference_matrix(27) # long time sage: H == H.T and is_hadamard_matrix(H) # long time True sage: H.diagonal() == [1]*676 and list(sum(H)) == [26]*676 # long time True
>>> from sage.all import * >>> H = rshcd_from_prime_power_and_conference_matrix(Integer(27)) # long time >>> H == H.T and is_hadamard_matrix(H) # long time True >>> H.diagonal() == [Integer(1)]*Integer(676) and list(sum(H)) == [Integer(26)]*Integer(676) # long time True
In this example the conference matrix is not Paley, as 45 is not a prime power
sage: H = rshcd_from_prime_power_and_conference_matrix(47) # not tested (long time)
>>> from sage.all import * >>> H = rshcd_from_prime_power_and_conference_matrix(Integer(47)) # not tested (long time)
REFERENCE:
- sage.combinat.matrices.hadamard_matrix.skew_hadamard_matrix(n, existence=False, skew_normalize=True, check=True, construction_name=False)[source]#
Tries to construct a skew Hadamard matrix.
A Hadamard matrix \(H\) is called skew if \(H=S-I\), for \(I\) the identity matrix and \(-S=S^\top\). Currently all orders \(\le 1200\) for which a construction is known are implemented. For \(n > 1200\), only some orders are available.
INPUT:
n
– integer; dimension of the matrixexistence
– boolean (default:False
); whether to build the matrix or merely query if a construction is available in Sage. When set toTrue
, the function returns:True
– meaning that Sage knows how to build the matrixUnknown
– meaning that Sage does not know how to build the matrix, but that the design may exist (seesage.misc.unknown
).False
– meaning that the matrix does not exist.
skew_normalize
– boolean (default:True
); whether to make the 1st row all-one, and adjust the 1st column accordinglycheck
– boolean (default:True
); whether to check that output is correct before returning it. As this is expected to be useless (but we are cautious guys), you may want to disable it whenever you want speedconstruction_name
– boolean (default:False
); if it isTrue
,existence
isTrue
, and a matrix exists, output the construction name. It has no effect ifexistence
is set toFalse
.
EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import skew_hadamard_matrix sage: skew_hadamard_matrix(12).det() 2985984 sage: 12^6 2985984 sage: skew_hadamard_matrix(1) [1] sage: skew_hadamard_matrix(2) [ 1 1] [-1 1] sage: skew_hadamard_matrix(196, existence=True, construction_name=True) 'Williamson-Goethals-Seidel: 196'
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import skew_hadamard_matrix >>> skew_hadamard_matrix(Integer(12)).det() 2985984 >>> Integer(12)**Integer(6) 2985984 >>> skew_hadamard_matrix(Integer(1)) [1] >>> skew_hadamard_matrix(Integer(2)) [ 1 1] [-1 1] >>> skew_hadamard_matrix(Integer(196), existence=True, construction_name=True) 'Williamson-Goethals-Seidel: 196'
REFERENCES:
- sage.combinat.matrices.hadamard_matrix.skew_hadamard_matrix_from_complementary_difference_sets(n, existence=False, check=True)[source]#
Construct a skew Hadamard matrix of order \(n=4(m+1)\) from complementary difference sets.
If \(A, B\) are complementary difference sets over a group of order \(2m+1\), then they can be used to construct a skew Hadamard matrix, as described in [BS1969].
The complementary difference sets are constructed using the function
sage.combinat.designs.difference_family.complementary_difference_sets()
.INPUT:
n
– positive integer; the order of the matrix to be constructedexistence
– boolean (default:False
); ifTrue
, only return whether the skew Hadamard matrix can be constructedcheck
– boolean (default:True
); ifTrue
, check that the result is a skew Hadamard matrix before returning it
OUTPUT:
If
existence=False
, returns the skew Hadamard matrix of order \(n\). It raises an error if \(n\) does not satisfy the required conditions. Ifexistence=True
, returns a boolean representing whether the matrix can be constructed or not.EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import skew_hadamard_matrix_from_complementary_difference_sets sage: skew_hadamard_matrix_from_complementary_difference_sets(20) 20 x 20 dense matrix over Integer Ring... sage: skew_hadamard_matrix_from_complementary_difference_sets(52, existence=True) True
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import skew_hadamard_matrix_from_complementary_difference_sets >>> skew_hadamard_matrix_from_complementary_difference_sets(Integer(20)) 20 x 20 dense matrix over Integer Ring... >>> skew_hadamard_matrix_from_complementary_difference_sets(Integer(52), existence=True) True
- sage.combinat.matrices.hadamard_matrix.skew_hadamard_matrix_from_good_matrices(a, b, c, d, check=True)[source]#
Construct skew Hadamard matrix from good matrices.
Given good matrices \(A\), \(B\), \(C\), \(D\) (\(A\) circulant, \(B, C, D\) back-circulant) they can be used to construct a skew Hadamard matrix using the following block matrix (as described in [Sze1988]):
\[\begin{split}\left(\begin{array}{rrrr} A & B & C & D \\ -B & A & D & -C \\ -C & -D & A & B \\ -D & C & -B & A \end{array}\right)\end{split}\]INPUT:
a
– (1,-1) list; the 1st row of \(A\)b
– (1,-1) list; the 1st row of \(B\)d
– (1,-1) list; the 1st row of \(C\)c
– (1,-1) list; the 1st row of \(D\)check
– boolean (default:True
); ifTrue
, check that the matrix is a skew Hadamard matrix before returning it
EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import skew_hadamard_matrix_from_good_matrices sage: a, b, c, d = ([1, 1, 1, -1, -1], [1, -1, 1, 1, -1], [1, -1, -1, -1, -1], [1, -1, -1, -1, -1]) sage: skew_hadamard_matrix_from_good_matrices(a, b, c, d) 20 x 20 dense matrix over Integer Ring...
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import skew_hadamard_matrix_from_good_matrices >>> a, b, c, d = ([Integer(1), Integer(1), Integer(1), -Integer(1), -Integer(1)], [Integer(1), -Integer(1), Integer(1), Integer(1), -Integer(1)], [Integer(1), -Integer(1), -Integer(1), -Integer(1), -Integer(1)], [Integer(1), -Integer(1), -Integer(1), -Integer(1), -Integer(1)]) >>> skew_hadamard_matrix_from_good_matrices(a, b, c, d) 20 x 20 dense matrix over Integer Ring...
- sage.combinat.matrices.hadamard_matrix.skew_hadamard_matrix_from_good_matrices_smallcases(n, existence=False, check=True)[source]#
Construct skew Hadamard matrices from good matrices for some small values of \(n=4m\), with \(m\) odd.
The function stores good matrices of odd orders \(\le 31\), taken from [Sze1988]. These are used to create skew Hadamard matrices of order \(4m\), \(1 \le m \le 31\) (\(m\) odd), using the function
skew_hadamard_matrix_from_good_matrices()
.ALGORITHM:
Given four sequences (stored in
E_sequences
) of length \(m\), they can be used to construct four \(E-sequences\) of length \(n=2m+1\), as follows:\[\begin{split}\begin{aligned} a &= 1, a_0, a_1, ..., a_m, -a_m, -a_{m-1}, ..., -a_0 \\ b &= 1, b_0, b_1, ..., b_m, b_m, b_{m-1}, ..., b_0 \\ c &= 1, c_0, c_1, ..., c_m, c_m, c_{m-1}, ..., c_0 \\ d &= 1, d_0, d_1, ..., d_m, d_m, d_{m-1}, ..., d_0 \\ \end{aligned}\end{split}\]These E-sequences will be the first rows of the four good matrices needed to construct a skew Hadamard matrix of order \(4n\).
INPUT:
n
– integer; the order of the skew Hadamard matrix to be constructedexistence
– boolean (default:False
); IfTrue
, only return whether the Hadamard matrix can be constructedcheck
– boolean (default:True
): ifTrue
, check that the matrix is a Hadamard matrix before returning it
OUTPUT:
If
existence=False
, returns the skew Hadamard matrix of order \(n\). It raises an error if no data is available to construct the matrix of the given order. Ifexistence=True
, returns a boolean representing whether the matrix can be constructed or not.EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import skew_hadamard_matrix_from_good_matrices_smallcases sage: skew_hadamard_matrix_from_good_matrices_smallcases(20) 20 x 20 dense matrix over Integer Ring... sage: skew_hadamard_matrix_from_good_matrices_smallcases(20, existence=True) True
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import skew_hadamard_matrix_from_good_matrices_smallcases >>> skew_hadamard_matrix_from_good_matrices_smallcases(Integer(20)) 20 x 20 dense matrix over Integer Ring... >>> skew_hadamard_matrix_from_good_matrices_smallcases(Integer(20), existence=True) True
- sage.combinat.matrices.hadamard_matrix.skew_hadamard_matrix_from_orthogonal_design(n, existence=False, check=True)[source]#
Construct skew Hadamard matrices of order \(mn(n - 1)\) if suitable orthogonal designs exist.
In [Seb1978] is proved that if amicable Hadamard matrices of order \(n\) and an orthogonal design of type \((1, m, mn - m - 1)\) in order \(mn\) exist, then a skew Hadamard matrix of order \(mn(n - 1)\) can be constructed. The paper uses amicable orthogonal designs instead of amicable Hadamard matrices, but the two are equivalent (see [Seb2017]).
Amicable Hadamard matrices are constructed using
amicable_hadamard_matrices()
, and the orthogonal designs are constructed using the Goethals-Seidel array, with data taken from [Seb2017].INPUT:
n
– positive integer; the order of the matrix to be constructedexistence
– boolean (default:False
); ifTrue
, only return whether the skew Hadamard matrix can be constructedcheck
– boolean (default:True
); ifTrue
, check that the result is a skew Hadamard matrix before returning it
OUTPUT:
If
existence=False
, returns the skew Hadamard matrix of order \(n\). It raises an error if a construction for order \(n\) is not yet implemented, or if \(n\) does not satisfy the constraint. Ifexistence=True
, returns a boolean representing whether the matrix can be constructed or not.EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import skew_hadamard_matrix_from_orthogonal_design sage: skew_hadamard_matrix_from_orthogonal_design(756) 756 x 756 dense matrix over Integer Ring...
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import skew_hadamard_matrix_from_orthogonal_design >>> skew_hadamard_matrix_from_orthogonal_design(Integer(756)) 756 x 756 dense matrix over Integer Ring...
If
existence
is True, the function returns a boolean:sage: skew_hadamard_matrix_from_orthogonal_design(200, existence=True) False
>>> from sage.all import * >>> skew_hadamard_matrix_from_orthogonal_design(Integer(200), existence=True) False
- sage.combinat.matrices.hadamard_matrix.skew_hadamard_matrix_spence_1975(n, existence=False, check=True)[source]#
Construct a skew Hadamard matrix of order \(n = 4(1 + q + q^2)\) using the Spence construction.
If \(n = 4(1 + q + q^2)\) where \(q\) is a prime power such that either \(1 + q + q^2\) is a prime congruent to \(3, 5, 7 \mod 8\) or \(3 + 2q + 2q^2\) is a prime power, then a skew Hadamard matrix of order \(n\) can be constructed using the Goethals Seidel array. The four matrices \(A, B, C, D\) plugged into the GS-array are created using complementary difference sets of order \(1 + q + q^2\) (which exist if \(q\) satisfies the conditions above), and a cyclic planar difference set with parameters \((1 + q^2 + q^4, 1 + q^2, 1)\). These are constructed by the functions
sage.combinat.designs.difference_family.complementary_difference_sets()
andsage.combinat.designs.difference_family.difference_family()
.For more details, see [Spe1975b].
INPUT:
n
– positive integer; the order of the matrix to be constructedexistence
– boolean (default:False
); ifTrue
, only return whether the Hadamard matrix can be constructedcheck
– boolean (default:True
); check that the result is a skew Hadamard matrix before returning it
OUTPUT:
If
existence=False
, returns the skew Hadamard matrix of order \(n\). It raises an error if \(n\) does not satisfy the required conditions. Ifexistence=True
, returns a boolean representing whether the matrix can be constructed or not.EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import skew_hadamard_matrix_spence_1975 sage: skew_hadamard_matrix_spence_1975(52) 52 x 52 dense matrix over Integer Ring...
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import skew_hadamard_matrix_spence_1975 >>> skew_hadamard_matrix_spence_1975(Integer(52)) 52 x 52 dense matrix over Integer Ring...
If
existence
is True, the function returns a boolean:sage: skew_hadamard_matrix_spence_1975(52, existence=True) True
>>> from sage.all import * >>> skew_hadamard_matrix_spence_1975(Integer(52), existence=True) True
- sage.combinat.matrices.hadamard_matrix.skew_hadamard_matrix_spence_construction(n, check=True)[source]#
Construct skew Hadamard matrix of order \(n\) using Spence constrution.
This function will construct skew Hadamard matrix of order \(n=2(q+1)\) where \(q\) is a prime power with \(q = 5\) (mod 8). The construction is taken from [Spe1977], and the relative difference sets are constructed using
sage.combinat.designs.difference_family.relative_difference_set_from_homomorphism()
.INPUT:
n
– positive integercheck
– boolean (default:True
); ifTrue
, check that the resulting matrix is Hadamard before returning it
OUTPUT:
If
n
satisfies the requirements described above, the function returns a \(n\times n\) Hadamard matrix. Otherwise, an exception is raised.EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import skew_hadamard_matrix_spence_construction sage: skew_hadamard_matrix_spence_construction(28) 28 x 28 dense matrix over Integer Ring...
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import skew_hadamard_matrix_spence_construction >>> skew_hadamard_matrix_spence_construction(Integer(28)) 28 x 28 dense matrix over Integer Ring...
- sage.combinat.matrices.hadamard_matrix.skew_hadamard_matrix_whiteman_construction(n, existence=False, check=True)[source]#
Construct a skew Hadamard matrix of order \(n=2(q+1)\) where \(q=p^t\) is a prime power with \(p \cong 5 \mod 8\) and \(t \cong 2 \mod 4\).
Assuming \(n\) satisfies the conditions above, it is possible to construct two supplementary difference sets \(A, B\) (see [Whi1971]), and these can be used to construct a skew Hadamard matrix, as described in [BS1969].
INPUT:
n
– positive integer; the order of the matrix to be constructedexistence
– boolean (default:False
); IfTrue
, only return whether the Hadamard matrix can be constructedcheck
– boolean (default:True
); ifTrue
, check that the result is a skew Hadamard matrix before returning it
OUTPUT:
If
existence=False
, returns the skew Hadamard matrix of order \(n\). It raises an error if \(n\) does not satisfy the required conditions. Ifexistence=True
, returns a boolean representing whether the matrix can be constructed or not.EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import skew_hadamard_matrix_whiteman_construction sage: skew_hadamard_matrix_whiteman_construction(52) 52 x 52 dense matrix over Integer Ring... sage: skew_hadamard_matrix_whiteman_construction(52, existence=True) True
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import skew_hadamard_matrix_whiteman_construction >>> skew_hadamard_matrix_whiteman_construction(Integer(52)) 52 x 52 dense matrix over Integer Ring... >>> skew_hadamard_matrix_whiteman_construction(Integer(52), existence=True) True
Note
A more general version of this construction is
skew_hadamard_matrix_from_complementary_difference_sets()
.
- sage.combinat.matrices.hadamard_matrix.symmetric_conference_matrix(n, check=True)[source]#
Tries to construct a symmetric conference matrix
A conference matrix is an \(n\times n\) matrix \(C\) with 0s on the main diagonal and 1s and -1s elsewhere, satisfying \(CC^\top=(n-1)I\). If \(C=C^\top\) then \(n \cong 2 \mod 4\) and \(C\) is Seidel adjacency matrix of a graph, whose descendent graphs are strongly regular graphs with parameters \((n-1,(n-2)/2,(n-6)/4,(n-2)/4)\), see Sec.10.4 of [BH2012]. Thus we build \(C\) from the Seidel adjacency matrix of the latter by adding row and column of 1s.
INPUT:
n
– integer; dimension of the matrixcheck
– boolean (default:True
); whether to check that output is correct before returning it. As this is expected to be useless (but we are cautious guys), you may want to disable it whenever you want speed
EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import symmetric_conference_matrix sage: C = symmetric_conference_matrix(10); C [ 0 1 1 1 1 1 1 1 1 1] [ 1 0 -1 -1 1 -1 1 1 1 -1] [ 1 -1 0 -1 1 1 -1 -1 1 1] [ 1 -1 -1 0 -1 1 1 1 -1 1] [ 1 1 1 -1 0 -1 -1 1 -1 1] [ 1 -1 1 1 -1 0 -1 1 1 -1] [ 1 1 -1 1 -1 -1 0 -1 1 1] [ 1 1 -1 1 1 1 -1 0 -1 -1] [ 1 1 1 -1 -1 1 1 -1 0 -1] [ 1 -1 1 1 1 -1 1 -1 -1 0] sage: C^2 == 9*identity_matrix(10) and C == C.T True
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import symmetric_conference_matrix >>> C = symmetric_conference_matrix(Integer(10)); C [ 0 1 1 1 1 1 1 1 1 1] [ 1 0 -1 -1 1 -1 1 1 1 -1] [ 1 -1 0 -1 1 1 -1 -1 1 1] [ 1 -1 -1 0 -1 1 1 1 -1 1] [ 1 1 1 -1 0 -1 -1 1 -1 1] [ 1 -1 1 1 -1 0 -1 1 1 -1] [ 1 1 -1 1 -1 -1 0 -1 1 1] [ 1 1 -1 1 1 1 -1 0 -1 -1] [ 1 1 1 -1 -1 1 1 -1 0 -1] [ 1 -1 1 1 1 -1 1 -1 -1 0] >>> C**Integer(2) == Integer(9)*identity_matrix(Integer(10)) and C == C.T True
- sage.combinat.matrices.hadamard_matrix.symmetric_conference_matrix_paley(n)[source]#
Construct a symmetric conference matrix of order n.
A conference matrix is an \(n\times n\) matrix \(C\) with 0s on the main diagonal and 1s and -1s elsewhere, satisfying \(CC^\top=(n-1)I\). This construction assumes that \(q = n-1\) is a prime power, with \(q \cong 1 \mod 4\). See [Hora] or [Lon2013].
These matrices are used in
hadamard_matrix_paleyII()
.INPUT:
n
– integer; the order of the symmetric conference matrix to construct
EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import symmetric_conference_matrix_paley sage: symmetric_conference_matrix_paley(6) [ 0 1 1 1 1 1] [ 1 0 1 -1 -1 1] [ 1 1 0 1 -1 -1] [ 1 -1 1 0 1 -1] [ 1 -1 -1 1 0 1] [ 1 1 -1 -1 1 0]
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import symmetric_conference_matrix_paley >>> symmetric_conference_matrix_paley(Integer(6)) [ 0 1 1 1 1 1] [ 1 0 1 -1 -1 1] [ 1 1 0 1 -1 -1] [ 1 -1 1 0 1 -1] [ 1 -1 -1 1 0 1] [ 1 1 -1 -1 1 0]
- sage.combinat.matrices.hadamard_matrix.szekeres_difference_set_pair(m, check=True)[source]#
Construct Szekeres \((2m+1,m,1)\)-cyclic difference family
Let \(4m+3\) be a prime power. Theorem 3 in [Sz1969] contains a construction of a pair of complementary difference sets \(A\), \(B\) in the subgroup \(G\) of the quadratic residues in \(F_{4m+3}^*\). Namely \(|A|=|B|=m\), \(a\in A\) whenever \(a-1\in G\), \(b\in B\) whenever \(b+1 \in G\). See also Theorem 2.6 in [SWW1972] (there the formula for \(B\) is correct, as opposed to (4.2) in [Sz1969], where the sign before \(1\) is wrong.
In modern terminology, for \(m>1\) the sets \(A\) and \(B\) form a
difference family
with parameters \((2m+1,m,1)\). I.e. each non-identity \(g \in G\) can be expressed uniquely as \(xy^{-1}\) for \(x,y \in A\) or \(x,y \in B\). Other, specific to this construction, properties of \(A\) and \(B\) are: for \(a\) in \(A\) one has \(a^{-1}\) not in \(A\), whereas for \(b\) in \(B\) one has \(b^{-1}\) in \(B\).INPUT:
m
– integer; dimension of the matrixcheck
– boolean (default:True
); whether to check \(A\) and \(B\) for correctness
EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import szekeres_difference_set_pair sage: G,A,B=szekeres_difference_set_pair(6) sage: G,A,B=szekeres_difference_set_pair(7)
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import szekeres_difference_set_pair >>> G,A,B=szekeres_difference_set_pair(Integer(6)) >>> G,A,B=szekeres_difference_set_pair(Integer(7))
REFERENCE:
- sage.combinat.matrices.hadamard_matrix.turyn_type_hadamard_matrix_smallcases(n, existence=False, check=True)[source]#
Construct a Hadamard matrix of order \(n\) from available 4-symbol \(\delta\) codes and Williamson quadruples.
The function looks for Baumert-Hall units and Williamson type matrices from
four_symbol_delta_code_smallcases()
andwilliamson_type_quadruples_smallcases()
and use them to construct a Hadamard matrix with the Turyn construction defined inhadamard_matrix_turyn_type()
.INPUT:
n
– integer; the order of the matrix to be constructedexistence
– boolean (default:False
): ifTrue
, only check if the matrix existscheck
– boolean (default:True
): ifTrue
, check that the matrix is a Hadamard matrix before returning
EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import turyn_type_hadamard_matrix_smallcases sage: turyn_type_hadamard_matrix_smallcases(28, existence=True) True sage: turyn_type_hadamard_matrix_smallcases(28) 28 x 28 dense matrix over Integer Ring...
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import turyn_type_hadamard_matrix_smallcases >>> turyn_type_hadamard_matrix_smallcases(Integer(28), existence=True) True >>> turyn_type_hadamard_matrix_smallcases(Integer(28)) 28 x 28 dense matrix over Integer Ring...
- sage.combinat.matrices.hadamard_matrix.typeI_matrix_difference_set(G, A)[source]#
(1,-1)-incidence type I matrix of a difference set \(A\) in \(G\)
Let \(A\) be a difference set in a group \(G\) of order \(n\). Return \(n\times n\) matrix \(M\) with \(M_{ij}=1\) if \(A_i A_j^{-1} \in A\), and \(M_{ij}=-1\) otherwise.
EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import szekeres_difference_set_pair sage: from sage.combinat.matrices.hadamard_matrix import typeI_matrix_difference_set sage: G,A,B=szekeres_difference_set_pair(2) sage: typeI_matrix_difference_set(G,A) [-1 1 -1 -1 1] [-1 -1 -1 1 1] [ 1 1 -1 -1 -1] [ 1 -1 1 -1 -1] [-1 -1 1 1 -1]
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import szekeres_difference_set_pair >>> from sage.combinat.matrices.hadamard_matrix import typeI_matrix_difference_set >>> G,A,B=szekeres_difference_set_pair(Integer(2)) >>> typeI_matrix_difference_set(G,A) [-1 1 -1 -1 1] [-1 -1 -1 1 1] [ 1 1 -1 -1 -1] [ 1 -1 1 -1 -1] [-1 -1 1 1 -1]
- sage.combinat.matrices.hadamard_matrix.williamson_goethals_seidel_skew_hadamard_matrix(a, b, c, d, check=True)[source]#
Williamson-Goethals-Seidel construction of a skew Hadamard matrix
Given \(n\times n\) (anti)circulant matrices \(A\), \(B\), \(C\), \(D\) with 1,-1 entries, and satisfying \(A+A^\top = 2I\), \(AA^\top + BB^\top + CC^\top + DD^\top = 4nI\), one can construct a skew Hadamard matrix of order \(4n\), cf. [GS70s].
INPUT:
a
– 1,-1 list; the 1st row of \(A\)b
– 1,-1 list; the 1st row of \(B\)d
– 1,-1 list; the 1st row of \(C\)c
– 1,-1 list; the 1st row of \(D\)check
– boolean (default:True
); ifTrue
, check that the resulting matrix is skew Hadamard before returning it
EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import williamson_goethals_seidel_skew_hadamard_matrix as WGS sage: a = [ 1, 1, 1, -1, 1, -1, 1, -1, -1] sage: b = [ 1, -1, 1, 1, -1, -1, 1, 1, -1] sage: c = [-1, -1]+[1]*6+[-1] sage: d = [ 1, 1, 1, -1, 1, 1, -1, 1, 1] sage: M = WGS(a,b,c,d,check=True)
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import williamson_goethals_seidel_skew_hadamard_matrix as WGS >>> a = [ Integer(1), Integer(1), Integer(1), -Integer(1), Integer(1), -Integer(1), Integer(1), -Integer(1), -Integer(1)] >>> b = [ Integer(1), -Integer(1), Integer(1), Integer(1), -Integer(1), -Integer(1), Integer(1), Integer(1), -Integer(1)] >>> c = [-Integer(1), -Integer(1)]+[Integer(1)]*Integer(6)+[-Integer(1)] >>> d = [ Integer(1), Integer(1), Integer(1), -Integer(1), Integer(1), Integer(1), -Integer(1), Integer(1), Integer(1)] >>> M = WGS(a,b,c,d,check=True)
REFERENCES:
- sage.combinat.matrices.hadamard_matrix.williamson_hadamard_matrix_smallcases(n, existence=False, check=True)[source]#
Construct Williamson type Hadamard matrices for some small values of n.
This function uses the data contained in
sage.combinat.matrices.hadamard_matrix.williamson_type_quadruples_smallcases()
to create Hadamard matrices of the Williamson type, using the construction fromsage.combinat.matrices.hadamard_matrix.hadamard_matrix_williamson_type()
.INPUT:
n
– integer; the order of the matrixexistence
– boolean (dafault:False
); ifTrue
, only check that we can do the constructioncheck
– boolean (default:True
); ifTrue
check the result
- sage.combinat.matrices.hadamard_matrix.williamson_type_quadruples_smallcases(n, existence=False)[source]#
Quadruples of matrices that can be used to construct Williamson type Hadamard matrices.
This function contains for some values of n, four \(n\times n\) matrices used in the Williamson construction of Hadamard matrices. Namely, the function returns the first row of 4 \(n\times n\) circulant matrices with the properties described in
sage.combinat.matrices.hadamard_matrix.hadamard_matrix_williamson_type()
. The matrices for \(n = 3, 5, ..., 29, 37, 43\) are given in [Ha83]. The matrices for \(n = 31, 33, 39, 41, 45, 49, 51, 55, 57, 61, 63\) are given in [Lon2013].INPUT:
n
– integer; the order of the matrices to be returnedexistence
– boolean (dafault:False
); ifTrue
, only check that we have the quadruple
OUTPUT:
If
existence
is false, returns a tuple containing four vectors, each being the first line of one of the four matrices. It raises an error if no such matrices are available. Ifexistence
is true, returns a boolean representing whether the matrices are available or not.EXAMPLES:
sage: from sage.combinat.matrices.hadamard_matrix import williamson_type_quadruples_smallcases sage: williamson_type_quadruples_smallcases(29) ((1, 1, 1, -1, -1, -1, 1, 1, -1, -1, 1, -1, 1, -1, -1, -1, -1, 1, -1, 1, -1, -1, 1, 1, -1, -1, -1, 1, 1), (1, -1, 1, -1, -1, -1, 1, 1, -1, -1, 1, -1, 1, 1, 1, 1, 1, 1, -1, 1, -1, -1, 1, 1, -1, -1, -1, 1, -1), (1, 1, 1, 1, -1, 1, 1, -1, 1, -1, -1, -1, 1, 1, 1, 1, 1, 1, -1, -1, -1, 1, -1, 1, 1, -1, 1, 1, 1), (1, 1, -1, -1, 1, -1, -1, 1, -1, 1, 1, 1, -1, 1, 1, 1, 1, -1, 1, 1, 1, -1, 1, -1, -1, 1, -1, -1, 1)) sage: williamson_type_quadruples_smallcases(43, existence=True) True
>>> from sage.all import * >>> from sage.combinat.matrices.hadamard_matrix import williamson_type_quadruples_smallcases >>> williamson_type_quadruples_smallcases(Integer(29)) ((1, 1, 1, -1, -1, -1, 1, 1, -1, -1, 1, -1, 1, -1, -1, -1, -1, 1, -1, 1, -1, -1, 1, 1, -1, -1, -1, 1, 1), (1, -1, 1, -1, -1, -1, 1, 1, -1, -1, 1, -1, 1, 1, 1, 1, 1, 1, -1, 1, -1, -1, 1, 1, -1, -1, -1, 1, -1), (1, 1, 1, 1, -1, 1, 1, -1, 1, -1, -1, -1, 1, 1, 1, 1, 1, 1, -1, -1, -1, 1, -1, 1, 1, -1, 1, 1, 1), (1, 1, -1, -1, 1, -1, -1, 1, -1, 1, 1, 1, -1, 1, 1, 1, 1, -1, 1, 1, 1, -1, 1, -1, -1, 1, -1, -1, 1)) >>> williamson_type_quadruples_smallcases(Integer(43), existence=True) True