def williamson_goethals_seidel_skew_hadamard_matrix(a, b, c, d, check=True): r""" 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 specifying the 1st row of `A` - ``b`` -- 1,-1 list specifying the 1st row of `B` - ``d`` -- 1,-1 list specifying the 1st row of `C` - ``c`` -- 1,-1 list specifying the 1st row of `D` 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) REFERENCES: .. [GS70s] \J.M. Goethals and J. J. Seidel, A skew Hadamard matrix of order 36, J. Aust. Math. Soc. 11(1970), 343-344 .. [Wall71] \J. Wallis, A skew-Hadamard matrix of order 92, Bull. Aust. Math. Soc. 5(1971), 203-204 .. [KoSt08] \C. Koukouvinos, S. Stylianou On skew-Hadamard matrices, Discrete Math. 308(2008) 2723-2731 """ n = len(a) R = matrix(ZZ, n, n, lambda i,j: 1 if i+j==n-1 else 0) A,B,C,D=map(matrix.circulant, [a,b,c,d]) if check: assert A*A.T+B*B.T+C*C.T+D*D.T==4*n*I(n) assert A+A.T==2*I(n) M = block_matrix([[ A, B*R, C*R, D*R], [-B*R, A, -D.T*R, C.T*R], [-C*R, D.T*R, A, -B.T*R], [-D*R, -C.T*R, B.T*R, A]]) if check: assert is_hadamard_matrix(M, normalized=False, skew=True) return M
def _helper_payley_matrix(n): r""" Return the marix constructed in Lemma 1.19 page 291 of [SWW72]_. This function return a `n^2` matrix `M` whose rows/columns are indexed by the element of a finite field on `n` elements `x_1,...,x_n`. The value `M_{i,j}` is equal to `\chi(x_i-x_j)`. Note that `n` must be an odd prime power. The elements `x_1,...,x_n` are ordered in such a way that the matrix is symmetric with respect to its second diagonal. The matrix is symmetric if n==4k+1, and skew-symmetric if n=4k-1. INPUT: - ``n`` -- a prime power .. SEEALSO:: :func:`rshcd_from_close_prime_powers` EXAMPLE:: sage: from sage.combinat.matrices.hadamard_matrix import _helper_payley_matrix sage: _helper_payley_matrix(5) [ 0 1 -1 -1 1] [ 1 0 1 -1 -1] [-1 1 0 1 -1] [-1 -1 1 0 1] [ 1 -1 -1 1 0] """ from sage.rings.finite_rings.constructor import FiniteField as GF K = GF(n,conway=True,prefix='x') # Order the elements of K in K_list # so that K_list[i] = -K_list[n-i-1] K_pairs = set(frozenset([x,-x]) for x in K) K_pairs.discard(frozenset([0])) K_list = [None]*n for i,(x,y) in enumerate(K_pairs): K_list[i] = x K_list[-i-1] = y K_list[n//2] = K(0) M = matrix(n,[[2*((x-y).is_square())-1 for x in K_list] for y in K_list]) M = M-I(n) assert (M*J(n)).is_zero() assert (M*M.transpose()) == n*I(n)-J(n) return M
def symmetric_conference_matrix(n, check=True): r""" 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 [BH12]_. 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 matrix - ``check`` (boolean) -- 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. Set to ``True`` by default. 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.graphs.strongly_regular_db import strongly_regular_graph as srg try: m = srg(n - 1, (n - 2) / 2, (n - 6) / 4, (n - 2) / 4) except ValueError: raise C = matrix([0] + [1] * (n - 1)).stack( matrix([1] * (n - 1)).stack(m.seidel_adjacency_matrix()).T) if check: assert (C == C.T and C**2 == (n - 1) * I(n)) return C
def rshcd_from_close_prime_powers(n): r""" 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 [GS70]_. Note that the authors of [SWW72]_ 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}` .. SEEALSO:: :func:`regular_symmetric_hadamard_matrix_with_constant_diagonal` 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] REFERENCE: .. [SWW72] A Street, W. Wallis, J. Wallis, Combinatorics: Room squares, sum-free sets, Hadamard matrices. Lecture notes in Mathematics 292 (1972). """ if n % 4: raise ValueError("n(={}) must be congruent to 0 mod 4") a, b = sorted([n - 1, n + 1], key=lambda x: -x % 4) Sa = _helper_payley_matrix(a) Sb = _helper_payley_matrix(b) U = matrix(a, [[int(i + j == a - 1) for i in range(a)] for j in range(a)]) K = (U * Sa).tensor_product(Sb) + U.tensor_product(J(b) - I(b)) - J( a).tensor_product(I(b)) F = lambda x: diagonal_matrix([-(-1)**i for i in range(x)]) G = block_diagonal_matrix([J(1), I(a).tensor_product(F(b))]) e = matrix(a * b, [1] * (a * b)) H = block_matrix(2, [-J(1), e.transpose(), e, K]) HH = G * H * G assert len(set(map(sum, HH))) == 1 assert HH**2 == n**2 * I(n**2) return HH
def _helper_payley_matrix(n, zero_position=True): r""" Return the matrix constructed in Lemma 1.19 page 291 of [SWW72]_. This function return a `n^2` matrix `M` whose rows/columns are indexed by the element of a finite field on `n` elements `x_1,...,x_n`. The value `M_{i,j}` is equal to `\chi(x_i-x_j)`. The elements `x_1,...,x_n` are ordered in such a way that the matrix (respectively, its submatrix obtained by removing first row and first column in the case ``zero_position=False``) is symmetric with respect to its second diagonal. The matrix is symmetric if `n=4k+1`, and skew-symmetric otherwise. INPUT: - ``n`` -- an odd prime power. - ``zero_position`` -- if it is true (default), place 0 of ``F_n`` in the middle, otherwise place it first. .. SEEALSO:: :func:`rshcd_from_close_prime_powers` EXAMPLES:: sage: from sage.combinat.matrices.hadamard_matrix import _helper_payley_matrix sage: _helper_payley_matrix(5) [ 0 1 -1 -1 1] [ 1 0 1 -1 -1] [-1 1 0 1 -1] [-1 -1 1 0 1] [ 1 -1 -1 1 0] TESTS:: sage: _helper_payley_matrix(11,zero_position=True) [ 0 -1 1 -1 -1 -1 1 1 1 -1 1] [ 1 0 -1 -1 1 -1 1 -1 1 1 -1] [-1 1 0 1 -1 -1 -1 -1 1 1 1] [ 1 1 -1 0 1 -1 -1 1 -1 -1 1] [ 1 -1 1 -1 0 1 -1 -1 -1 1 1] [ 1 1 1 1 -1 0 1 -1 -1 -1 -1] [-1 -1 1 1 1 -1 0 1 -1 1 -1] [-1 1 1 -1 1 1 -1 0 1 -1 -1] [-1 -1 -1 1 1 1 1 -1 0 -1 1] [ 1 -1 -1 1 -1 1 -1 1 1 0 -1] [-1 1 -1 -1 -1 1 1 1 -1 1 0] sage: _helper_payley_matrix(11,zero_position=False) [ 0 1 1 1 1 -1 1 -1 -1 -1 -1] [-1 0 -1 1 -1 -1 1 1 1 -1 1] [-1 1 0 -1 -1 1 1 -1 1 1 -1] [-1 -1 1 0 1 -1 -1 -1 1 1 1] [-1 1 1 -1 0 1 -1 1 -1 -1 1] [ 1 1 -1 1 -1 0 -1 -1 -1 1 1] [-1 -1 -1 1 1 1 0 1 -1 1 -1] [ 1 -1 1 1 -1 1 -1 0 1 -1 -1] [ 1 -1 -1 -1 1 1 1 -1 0 -1 1] [ 1 1 -1 -1 1 -1 -1 1 1 0 -1] [ 1 -1 1 -1 -1 -1 1 1 -1 1 0] """ from sage.rings.finite_rings.finite_field_constructor import FiniteField as GF K = GF(n, prefix='x') # Order the elements of K in K_list # so that K_list[i] = -K_list[n-i-1] K_pairs = set(frozenset([x, -x]) for x in K) K_pairs.discard(frozenset([0])) K_list = [None] * n if zero_position: zero_position = n // 2 shift = 0 else: shift = 1 for i, (x, y) in enumerate(K_pairs): K_list[i + shift] = x K_list[-i - 1] = y K_list[zero_position] = K(0) M = matrix(n, [[2 * ((x - y).is_square()) - 1 for x in K_list] for y in K_list]) M = M - I(n) assert (M * J(n)).is_zero() assert (M * M.transpose()) == n * I(n) - J(n) return M
def regular_symmetric_hadamard_matrix_with_constant_diagonal( n, e, existence=False): r""" 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 [HX10]_ or 10.5.1 in [BH12]_. For the case `n=324`, see :func:`RSHCD_324` and [CP16]_. INPUT: - ``n`` (integer) -- side of the matrix - ``e`` -- one of `-1` or `+1`, equal to 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] 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 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 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) 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) REFERENCE: .. [BH12] \A. Brouwer and W. Haemers, Spectra of graphs, Springer, 2012, http://homepages.cwi.nl/~aeb/math/ipm/ipm.pdf .. [HX10] \W. Haemers and Q. Xiang, Strongly regular graphs with parameters `(4m^4,2m^4+m^2,m^4+m^2,m^4+m^2)` exist for all `m>1`, European Journal of Combinatorics, Volume 31, Issue 6, August 2010, Pages 1553-1559, :doi:`10.1016/j.ejc.2009.07.009` """ if existence and (n, e) in _rshcd_cache: return _rshcd_cache[n, e] from sage.graphs.strongly_regular_db import strongly_regular_graph def true(): _rshcd_cache[n, e] = True return True M = None if abs(e) != 1: raise ValueError sqn = None if is_square(n): sqn = int(sqrt(n)) if n < 0: if existence: return False raise ValueError elif n == 4: if existence: return true() if e == 1: M = J(4) - 2 * matrix(4, [[int(i + j == 3) for i in range(4)] for j in range(4)]) else: M = -J(4) + 2 * I(4) elif n == 36: if existence: return true() if e == 1: M = strongly_regular_graph(36, 15, 6, 6).adjacency_matrix() M = J(36) - 2 * M else: M = strongly_regular_graph(36, 14, 4, 6).adjacency_matrix() M = -J(36) + 2 * M + 2 * I(36) elif n == 100: if existence: return true() if e == -1: M = strongly_regular_graph(100, 44, 18, 20).adjacency_matrix() M = 2 * M - J(100) + 2 * I(100) else: M = strongly_regular_graph(100, 45, 20, 20).adjacency_matrix() M = J(100) - 2 * M elif n == 196 and e == 1: if existence: return true() M = strongly_regular_graph(196, 91, 42, 42).adjacency_matrix() M = J(196) - 2 * M elif n == 324: if existence: return true() M = RSHCD_324(e) elif (e == 1 and n % 16 == 0 and not sqn is None and is_prime_power(sqn - 1) and is_prime_power(sqn + 1)): if existence: return true() M = -rshcd_from_close_prime_powers(sqn) elif (e == 1 and not sqn is None and sqn % 4 == 2 and True == strongly_regular_graph(sqn - 1, (sqn - 2) // 2, (sqn - 6) // 4, existence=True) and is_prime_power(ZZ(sqn + 1))): if existence: return true() M = rshcd_from_prime_power_and_conference_matrix(sqn + 1) # Recursive construction: the kronecker product of two RSHCD is a RSHCD else: from itertools import product for n1, e1 in product(divisors(n)[1:-1], [-1, 1]): e2 = e1 * e n2 = n // n1 if (regular_symmetric_hadamard_matrix_with_constant_diagonal( n1, e1, existence=True) and regular_symmetric_hadamard_matrix_with_constant_diagonal( n2, e2, existence=True)): if existence: return true() M1 = regular_symmetric_hadamard_matrix_with_constant_diagonal( n1, e1) M2 = regular_symmetric_hadamard_matrix_with_constant_diagonal( n2, e2) M = M1.tensor_product(M2) break if M is None: from sage.misc.unknown import Unknown _rshcd_cache[n, e] = Unknown if existence: return Unknown raise ValueError("I do not know how to build a {}-RSHCD".format( (n, e))) assert M * M.transpose() == n * I(n) assert set(map(sum, M)) == {ZZ(e * sqn)} return M
def rshcd_from_prime_power_and_conference_matrix(n): r""" 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 [WW72]_. In [SWW72]_ this construction (Theorem 5.15 and Corollary 5.16) is reproduced with a typo. Note that [WW72]_ refers to [Sz69]_ for the construction, provided by :func:`szekeres_difference_set_pair`, of complementary difference sets, and the latter has a typo. From a :func:`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 .. SEEALSO:: :func:`regular_symmetric_hadamard_matrix_with_constant_diagonal` 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 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 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) REFERENCE: .. [WW72] \J. Wallis and A.L. Whiteman, Some classes of Hadamard matrices with constant diagonal, Bull. Austral. Math. Soc. 7(1972), 233-249 """ from sage.graphs.strongly_regular_db import strongly_regular_graph as srg if is_prime_power(n) and 2 == (n - 1) % 4: try: M = srg(n - 2, (n - 3) // 2, (n - 7) // 4) except ValueError: return m = (n - 3) // 4 Q, X, Y = szekeres_difference_set_pair(m) B = typeI_matrix_difference_set(Q, X) A = -typeI_matrix_difference_set(Q, Y) # must be symmetric W = M.seidel_adjacency_matrix() f = J(1, 4 * m + 1) e = J(1, 2 * m + 1) JJ = J(2 * m + 1, 2 * m + 1) II = I(n - 2) Ib = I(2 * m + 1) J4m = J(4 * m + 1, 4 * m + 1) H34 = -(B + Ib).tensor_product(W) + Ib.tensor_product(J4m) + ( Ib - JJ).tensor_product(II) A_t_W = A.tensor_product(W) e_t_f = e.tensor_product(f) H = block_matrix( [[J(1, 1), f, e_t_f, -e_t_f], [f.T, J4m, e.tensor_product(W - II), e.tensor_product(W + II)], [ e_t_f.T, (e.T).tensor_product(W - II), A_t_W + JJ.tensor_product(II), H34 ], [ -e_t_f.T, (e.T).tensor_product(W + II), H34.T, -A_t_W + JJ.tensor_product(II) ]]) return H
def skew_hadamard_matrix(n, existence=False, skew_normalize=True, check=True): r""" 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 constructions from Section 14.1 of [Ha83]_ and few more exotic ones are implemented. INPUT: - ``n`` (integer) -- dimension of the matrix - ``existence`` (boolean) -- whether to build the matrix or merely query if a construction is available in Sage. When set to ``True``, the function returns: - ``True`` -- meaning that Sage knows how to build the matrix - ``Unknown`` -- meaning that Sage does not know how to build the matrix, but that the design may exist (see :mod:`sage.misc.unknown`). - ``False`` -- meaning that the matrix does not exist. - ``skew_normalize`` (boolean) -- whether to make the 1st row all-one, and adjust the 1st column accordingly. Set to ``True`` by default. - ``check`` (boolean) -- 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. Set to ``True`` by default. 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] TESTS:: sage: skew_hadamard_matrix(10,existence=True) False sage: skew_hadamard_matrix(12,existence=True) True sage: skew_hadamard_matrix(784,existence=True) True sage: skew_hadamard_matrix(10) Traceback (most recent call last): ... ValueError: A skew Hadamard matrix of order 10 does not exist sage: skew_hadamard_matrix(36) 36 x 36 dense matrix over Integer Ring... sage: skew_hadamard_matrix(36)==skew_hadamard_matrix(36,skew_normalize=False) False sage: skew_hadamard_matrix(52) 52 x 52 dense matrix over Integer Ring... sage: skew_hadamard_matrix(92) 92 x 92 dense matrix over Integer Ring... sage: skew_hadamard_matrix(816) # long time 816 x 816 dense matrix over Integer Ring... sage: skew_hadamard_matrix(100) Traceback (most recent call last): ... ValueError: A skew Hadamard matrix of order 100 is not yet implemented. sage: skew_hadamard_matrix(100,existence=True) Unknown REFERENCES: .. [Ha83] \M. Hall, Combinatorial Theory, 2nd edition, Wiley, 1983 """ def true(): _skew_had_cache[n] = True return True M = None if existence and n in _skew_had_cache: return True if not (n % 4 == 0) and (n > 2): if existence: return False raise ValueError("A skew Hadamard matrix of order %s does not exist" % n) if n == 2: if existence: return true() M = matrix([[1, 1], [-1, 1]]) elif n == 1: if existence: return true() M = matrix([1]) elif is_prime_power(n - 1) and ((n - 1) % 4 == 3): if existence: return true() M = hadamard_matrix_paleyI(n, normalize=False) elif n % 8 == 0: if skew_hadamard_matrix(n // 2, existence=True): # (Lemma 14.1.6 in [Ha83]_) if existence: return true() H = skew_hadamard_matrix(n // 2, check=False) M = block_matrix([[H, H], [-H.T, H.T]]) else: # try Williamson construction (Lemma 14.1.5 in [Ha83]_) for d in divisors(n)[2:-2]: # skip 1, 2, n/2, and n n1 = n // d if is_prime_power(d - 1) and (d % 4 == 0) and (n1 % 4 == 0)\ and skew_hadamard_matrix(n1,existence=True): if existence: return true() H = skew_hadamard_matrix(n1, check=False) - I(n1) U = matrix(ZZ, d, lambda i, j: -1 if i==j==0 else\ 1 if i==j==1 or (i>1 and j-1==d-i)\ else 0) A = block_matrix( [[matrix([0]), matrix(ZZ, 1, d - 1, [1] * (d - 1))], [ matrix(ZZ, d - 1, 1, [-1] * (d - 1)), _helper_payley_matrix(d - 1, zero_position=0) ]]) + I(d) M = A.tensor_product(I(n1)) + (U * A).tensor_product(H) break if M is None: # try Williamson-Goethals-Seidel construction if GS_skew_hadamard_smallcases(n, existence=True): if existence: return true() M = GS_skew_hadamard_smallcases(n) else: if existence: return Unknown raise ValueError( "A skew Hadamard matrix of order %s is not yet implemented." % n) if skew_normalize: dd = diagonal_matrix(M[0]) M = dd * M * dd if check: assert is_hadamard_matrix(M, normalized=False, skew=True) if skew_normalize: from sage.modules.free_module_element import vector assert M[0] == vector([1] * n) _skew_had_cache[n] = True return M