Generic structures for linear codes over the Hamming metric#

Linear Codes#

Let \(F = \GF{q}\) be a finite field. A rank \(k\) linear subspace of the vector space \(F^n\) is called an \([n, k]\)-linear code, \(n\) being the length of the code and \(k\) its dimension. Elements of a code \(C\) are called codewords.

A linear map from \(F^k\) to an \([n,k]\) code \(C\) is called an “encoding”, and it can be represented as a \(k \times n\) matrix, called a generator matrix. Alternatively, \(C\) can be represented by its orthogonal complement in \(F^n\), i.e. the \((n-k)\)-dimensional vector space \(C^\perp\) such that the inner product of any element from \(C\) and any element from \(C^\perp\) is zero. \(C^\perp\) is called the dual code of \(C\), and any generator matrix for \(C^\perp\) is called a parity check matrix for \(C\).

We commonly endow \(F^n\) with the Hamming metric, i.e. the weight of a vector is the number of non-zero elements in it. The central operation of a linear code is then “decoding”: given a linear code \(C \subset F^n\) and a “received word” \(r \in F^n\) , retrieve the codeword \(c \in C\) such that the Hamming distance between \(r\) and \(c\) is minimal.

Families or Generic codes#

Linear codes are either studied as generic vector spaces without any known structure, or as particular sub-families with special properties.

The class sage.coding.linear_code.LinearCode is used to represent the former.

For the latter, these will be represented by specialised classes; for instance, the family of Hamming codes are represented by the class sage.coding.hamming_code.HammingCode. Type codes.<tab> for a list of all code families known to Sage. Such code family classes should inherit from the abstract base class sage.coding.linear_code.AbstractLinearCode.

AbstractLinearCode#

This is a base class designed to contain methods, features and parameters shared by every linear code. For instance, generic algorithms for computing the minimum distance, the covering radius, etc. Many of these algorithms are slow, e.g. exponential in the code length. For specific subfamilies, better algorithms or even closed formulas might be known, in which case the respective method should be overridden.

AbstractLinearCode is an abstract class for linear codes, so any linear code class should inherit from this class. Also AbstractLinearCode should never itself be instantiated.

See sage.coding.linear_code.AbstractLinearCode for details and examples.

LinearCode#

This class is used to represent arbitrary and unstructured linear codes. It mostly rely directly on generic methods provided by AbstractLinearCode, which means that basic operations on the code (e.g. computation of the minimum distance) will use slow algorithms.

A LinearCode is instantiated by providing a generator matrix:

sage: M = matrix(GF(2), [[1, 0, 0, 1, 0],\
                         [0, 1, 0, 1, 1],\
                         [0, 0, 1, 1, 1]])
sage: C = codes.LinearCode(M)
sage: C
[5, 3] linear code over GF(2)
sage: C.generator_matrix()
[1 0 0 1 0]
[0 1 0 1 1]
[0 0 1 1 1]

sage: MS = MatrixSpace(GF(2),4,7)
sage: G = MS([[1,1,1,0,0,0,0], [1,0,0,1,1,0,0], [0,1,0,1,0,1,0], [1,1,0,1,0,0,1]])
sage: C = LinearCode(G)
sage: C.basis()
[
(1, 1, 1, 0, 0, 0, 0),
(1, 0, 0, 1, 1, 0, 0),
(0, 1, 0, 1, 0, 1, 0),
(1, 1, 0, 1, 0, 0, 1)
]
sage: c = C.basis()[1]
sage: c in C
True
sage: c.nonzero_positions()
[0, 3, 4]
sage: c.support()
[0, 3, 4]
sage: c.parent()
Vector space of dimension 7 over Finite Field of size 2
>>> from sage.all import *
>>> M = matrix(GF(Integer(2)), [[Integer(1), Integer(0), Integer(0), Integer(1), Integer(0)],\
                         [0, 1, 0, 1, 1],\
                         [0, 0, 1, 1, 1]])
>>> C = codes.LinearCode(M)
>>> C
[5, 3] linear code over GF(2)
>>> C.generator_matrix()
[1 0 0 1 0]
[0 1 0 1 1]
[0 0 1 1 1]

>>> MS = MatrixSpace(GF(Integer(2)),Integer(4),Integer(7))
>>> G = MS([[Integer(1),Integer(1),Integer(1),Integer(0),Integer(0),Integer(0),Integer(0)], [Integer(1),Integer(0),Integer(0),Integer(1),Integer(1),Integer(0),Integer(0)], [Integer(0),Integer(1),Integer(0),Integer(1),Integer(0),Integer(1),Integer(0)], [Integer(1),Integer(1),Integer(0),Integer(1),Integer(0),Integer(0),Integer(1)]])
>>> C = LinearCode(G)
>>> C.basis()
[
(1, 1, 1, 0, 0, 0, 0),
(1, 0, 0, 1, 1, 0, 0),
(0, 1, 0, 1, 0, 1, 0),
(1, 1, 0, 1, 0, 0, 1)
]
>>> c = C.basis()[Integer(1)]
>>> c in C
True
>>> c.nonzero_positions()
[0, 3, 4]
>>> c.support()
[0, 3, 4]
>>> c.parent()
Vector space of dimension 7 over Finite Field of size 2

Further references#

If you want to get started on Sage’s linear codes library, see https://doc.sagemath.org/html/en/thematic_tutorials/coding_theory.html

If you want to learn more on the design of this library, see https://doc.sagemath.org/html/en/thematic_tutorials/structures_in_coding_theory.html

REFERENCES:

AUTHORS:

  • David Joyner (2005-11-22, 2006-12-03): initial version

  • William Stein (2006-01-23): Inclusion in Sage

  • David Joyner (2006-01-30, 2006-04): small fixes

  • David Joyner (2006-07): added documentation, group-theoretical methods, ToricCode

  • David Joyner (2006-08): hopeful latex fixes to documentation, added list and __iter__ methods to LinearCode and examples, added hamming_weight function, fixed random method to return a vector, TrivialCode, fixed subtle bug in dual_code, added galois_closure method, fixed mysterious bug in permutation_automorphism_group (GAP was over-using “G” somehow?)

  • David Joyner (2006-08): hopeful latex fixes to documentation, added CyclicCode, best_known_linear_code, bounds_minimum_distance, assmus_mattson_designs (implementing Assmus-Mattson Theorem).

  • David Joyner (2006-09): modified decode syntax, fixed bug in is_galois_closed, added LinearCode_from_vectorspace, extended_code, zeta_function

  • Nick Alexander (2006-12-10): factor GUAVA code to guava.py

  • David Joyner (2007-05): added methods punctured, shortened, divisor, characteristic_polynomial, binomial_moment, support for LinearCode. Completely rewritten zeta_function (old version is now zeta_function2) and a new function, LinearCodeFromVectorSpace.

  • David Joyner (2007-11): added zeta_polynomial, weight_enumerator, chinen_polynomial; improved best_known_code; made some pythonic revisions; added is_equivalent (for binary codes)

  • David Joyner (2008-01): fixed bug in decode reported by Harald Schilly, (with Mike Hansen) added some doctests.

  • David Joyner (2008-02): translated standard_form, dual_code to Python.

  • David Joyner (2008-03): translated punctured, shortened, extended_code, random (and renamed random to random_element), deleted zeta_function2, zeta_function3, added wrapper automorphism_group_binary_code to Robert Miller’s code), added direct_sum_code, is_subcode, is_self_dual, is_self_orthogonal, redundancy_matrix, did some alphabetical reorganizing to make the file more readable. Fixed a bug in permutation_automorphism_group which caused it to crash.

  • David Joyner (2008-03): fixed bugs in spectrum and zeta_polynomial, which misbehaved over non-prime base rings.

  • David Joyner (2008-10): use CJ Tjhal’s MinimumWeight if char = 2 or 3 for min_dist; add is_permutation_equivalent and improve permutation_automorphism_group using an interface with Robert Miller’s code; added interface with Leon’s code for the spectrum method.

  • David Joyner (2009-02): added native decoding methods (see module_decoder.py)

  • David Joyner (2009-05): removed dependence on Guava, allowing it to be an option. Fixed errors in some docstrings.

  • Kwankyu Lee (2010-01): added methods generator_matrix_systematic, information_set, and magma interface for linear codes.

  • Niles Johnson (2010-08): Issue #3893: random_element() should pass on *args and **kwds.

  • Thomas Feulner (2012-11): Issue #13723: deprecation of hamming_weight()

  • Thomas Feulner (2013-10): added methods to compute a canonical representative and the automorphism group

class sage.coding.linear_code.AbstractLinearCode(base_field, length, default_encoder_name, default_decoder_name)[source]#

Bases: AbstractLinearCodeNoMetric

Abstract base class for linear codes.

This class contains all methods that can be used on Linear Codes and on Linear Codes families. So, every Linear Code-related class should inherit from this abstract class.

To implement a linear code, you need to:

  • inherit from AbstractLinearCode

  • call AbstractLinearCode __init__ method in the subclass constructor. Example: super().__init__(base_field, length, "EncoderName", "DecoderName"). By doing that, your subclass will have its length parameter initialized and will be properly set as a member of the category framework. You need of course to complete the constructor by adding any additional parameter needed to describe properly the code defined in the subclass.

  • Add the following two lines on the class level:

    _registered_encoders = {}
    _registered_decoders = {}
    
  • fill the dictionary of its encoders in sage.coding.__init__.py file. Example: I want to link the encoder MyEncoderClass to MyNewCodeClass under the name MyEncoderName. All I need to do is to write this line in the __init__.py file: MyNewCodeClass._registered_encoders["NameOfMyEncoder"] = MyEncoderClass and all instances of MyNewCodeClass will be able to use instances of MyEncoderClass.

  • fill the dictionary of its decoders in sage.coding.__init__ file. Example: I want to link the encoder MyDecoderClass to MyNewCodeClass under the name MyDecoderName. All I need to do is to write this line in the __init__.py file: MyNewCodeClass._registered_decoders["NameOfMyDecoder"] = MyDecoderClass and all instances of MyNewCodeClass will be able to use instances of MyDecoderClass.

As the class AbstractLinearCode is not designed to be instantiated, it does not have any representation methods. You should implement _repr_ and _latex_ methods in the subclass.

Note

AbstractLinearCode has a generic implementation of the method __eq__ which uses the generator matrix and is quite slow. In subclasses you are encouraged to override __eq__ and __hash__.

Warning

The default encoder should always have \(F^{k}\) as message space, with \(k\) the dimension of the code and \(F\) is the base ring of the code.

A lot of methods of the abstract class rely on the knowledge of a generator matrix. It is thus strongly recommended to set an encoder with a generator matrix implemented as a default encoder.

assmus_mattson_designs(t, mode=None)[source]#

Assmus and Mattson Theorem (section 8.4, page 303 of [HP2003]): Let \(A_0, A_1, ..., A_n\) be the weights of the codewords in a binary linear \([n , k, d]\) code \(C\), and let \(A_0^*, A_1^*, ..., A_n^*\) be the weights of the codewords in its dual \([n, n-k, d^*]\) code \(C^*\). Fix a \(t\), \(0<t<d\), and let

\[s = |\{ i\ |\ A_i^* \not= 0, 0< i \leq n-t\}|.\]

Assume \(s\leq d-t\).

  1. If \(A_i\not= 0\) and \(d\leq i\leq n\) then \(C_i = \{ c \in C\ |\ wt(c) = i\}\) holds a simple t-design.

  2. If \(A_i^*\not= 0\) and \(d*\leq i\leq n-t\) then \(C_i^* = \{ c \in C^*\ |\ wt(c) = i\}\) holds a simple t-design.

A block design is a pair \((X,B)\), where \(X\) is a non-empty finite set of \(v>0\) elements called points, and \(B\) is a non-empty finite multiset of size b whose elements are called blocks, such that each block is a non-empty finite multiset of \(k\) points. \(A\) design without repeated blocks is called a simple block design. If every subset of points of size \(t\) is contained in exactly \(\lambda\) blocks the block design is called a \(t-(v,k,\lambda)\) design (or simply a \(t\)-design when the parameters are not specified). When \(\lambda=1\) then the block design is called a \(S(t,k,v)\) Steiner system.

In the Assmus and Mattson Theorem (1), \(X\) is the set \(\{1,2,...,n\}\) of coordinate locations and \(B = \{supp(c)\ |\ c \in C_i\}\) is the set of supports of the codewords of \(C\) of weight \(i\). Therefore, the parameters of the \(t\)-design for \(C_i\) are

t =       given
v =       n
k =       i   (k not to be confused with dim(C))
b =       Ai
lambda = b*binomial(k,t)/binomial(v,t) (by Theorem 8.1.6,
                                           p 294, in [HP2003]_)

Setting the mode="verbose" option prints out the values of the parameters.

The first example below means that the binary [24,12,8]-code C has the property that the (support of the) codewords of weight 8 (resp., 12, 16) form a 5-design. Similarly for its dual code \(C^*\) (of course \(C=C^*\) in this case, so this info is extraneous). The test fails to produce 6-designs (ie, the hypotheses of the theorem fail to hold, not that the 6-designs definitely don’t exist). The command assmus_mattson_designs(C,5,mode="verbose") returns the same value but prints out more detailed information.

The second example below illustrates the blocks of the 5-(24, 8, 1) design (i.e., the S(5,8,24) Steiner system).

EXAMPLES:

sage: C = codes.GolayCode(GF(2))             #  example 1
sage: C.assmus_mattson_designs(5)
['weights from C: ', [8, 12, 16, 24],
 'designs from C: ',  [[5, (24, 8, 1)], [5, (24, 12, 48)],
                       [5, (24, 16, 78)], [5, (24, 24, 1)]],
 'weights from C*: ', [8, 12, 16],
 'designs from C*: ', [[5, (24, 8, 1)], [5, (24, 12, 48)], [5, (24, 16, 78)]]]
sage: C.assmus_mattson_designs(6)
0
sage: X = range(24)                          #  example 2
sage: blocks = [c.support()  # long time
....:           for c in C if c.hamming_weight()==8]; len(blocks)
759
>>> from sage.all import *
>>> C = codes.GolayCode(GF(Integer(2)))             #  example 1
>>> C.assmus_mattson_designs(Integer(5))
['weights from C: ', [8, 12, 16, 24],
 'designs from C: ',  [[5, (24, 8, 1)], [5, (24, 12, 48)],
                       [5, (24, 16, 78)], [5, (24, 24, 1)]],
 'weights from C*: ', [8, 12, 16],
 'designs from C*: ', [[5, (24, 8, 1)], [5, (24, 12, 48)], [5, (24, 16, 78)]]]
>>> C.assmus_mattson_designs(Integer(6))
0
>>> X = range(Integer(24))                          #  example 2
>>> blocks = [c.support()  # long time
...           for c in C if c.hamming_weight()==Integer(8)]; len(blocks)
759
automorphism_group_gens(equivalence='semilinear')[source]#

Return generators of the automorphism group of self.

INPUT:

  • equivalence (optional) – which defines the acting group, either

    • "permutational"

    • "linear"

    • "semilinear"

OUTPUT:

  • generators of the automorphism group of self

  • the order of the automorphism group of self

EXAMPLES:

Note, this result can depend on the PRNG state in libgap in a way that depends on which packages are loaded, so we must re-seed GAP to ensure a consistent result for this example:

sage: # needs sage.libs.gap
sage: libgap.set_seed(0)
0
sage: C = codes.HammingCode(GF(4, 'z'), 3)
sage: C.automorphism_group_gens()
([((1, 1, 1, z, z + 1, 1, 1, 1, 1, z + 1, z, z, z + 1, z + 1,
    z + 1, 1, z + 1, z, z, 1, z);
   (1,13,14,20)(2,21,8,18,7,16,19,15)(3,10,5,12,17,9,6,4),
   Ring endomorphism of Finite Field in z of size 2^2
     Defn: z |--> z + 1),
  ((z, 1, z, z, z, z + 1, z, z, z, z, z, z, z + 1, z, z, z,
    z, z + 1, z, z, z);
   (1,11,5,12,3,19)(2,8)(6,18,13)(7,17,15)(9,10,14,16,20,21),
   Ring endomorphism of Finite Field in z of size 2^2
     Defn: z |--> z + 1),
  ((z, z, z, z, z, z, z, z, z, z, z, z, z, z, z, z, z, z, z, z, z);
   (),
   Ring endomorphism of Finite Field in z of size 2^2
     Defn: z |--> z)],
 362880)
sage: C.automorphism_group_gens(equivalence="linear")
([((z, 1, z + 1, z + 1, 1, z + 1, z, 1, z + 1, z + 1, 1, z, 1, z + 1,
    z, 1, z, 1, z + 1, 1, 1);
   (1,12,11,10,6,8,9,20,13,21,5,14,3,16,17,19,7,4,2,15,18),
   Ring endomorphism of Finite Field in z of size 2^2
     Defn: z |--> z),
  ((z + 1, z + 1, z + 1, z, 1, 1, z, z, 1, z + 1, z, 1, 1, z, 1, z + 1,
    z, z + 1, z + 1, 1, z);
   (1,3,18,2,17,6,19)(4,15,13,20,7,14,16)(5,11,8,21,12,9,10),
   Ring endomorphism of Finite Field in z of size 2^2
     Defn: z |--> z),
  ((z + 1, z + 1, z + 1, z + 1, z + 1, z + 1, z + 1, z + 1, z + 1, z + 1,
    z + 1, z + 1, z + 1, z + 1, z + 1, z + 1, z + 1, z + 1, z + 1, z + 1, z + 1);
   (),
   Ring endomorphism of Finite Field in z of size 2^2
     Defn: z |--> z)],
 181440)
sage: C.automorphism_group_gens(equivalence="permutational")
([((1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1);
   (1,11)(3,10)(4,9)(5,7)(12,21)(14,20)(15,19)(16,17),
   Ring endomorphism of Finite Field in z of size 2^2
     Defn: z |--> z),
  ((1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1);
   (2,18)(3,19)(4,10)(5,16)(8,13)(9,14)(11,21)(15,20),
   Ring endomorphism of Finite Field in z of size 2^2
     Defn: z |--> z),
  ((1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1);
   (1,19)(3,17)(4,21)(5,20)(7,14)(9,12)(10,16)(11,15),
   Ring endomorphism of Finite Field in z of size 2^2
     Defn: z |--> z),
  ((1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1);
   (2,13)(3,14)(4,20)(5,11)(8,18)(9,19)(10,15)(16,21),
   Ring endomorphism of Finite Field in z of size 2^2
     Defn: z |--> z)],
64)
>>> from sage.all import *
>>> # needs sage.libs.gap
>>> libgap.set_seed(Integer(0))
0
>>> C = codes.HammingCode(GF(Integer(4), 'z'), Integer(3))
>>> C.automorphism_group_gens()
([((1, 1, 1, z, z + 1, 1, 1, 1, 1, z + 1, z, z, z + 1, z + 1,
    z + 1, 1, z + 1, z, z, 1, z);
   (1,13,14,20)(2,21,8,18,7,16,19,15)(3,10,5,12,17,9,6,4),
   Ring endomorphism of Finite Field in z of size 2^2
     Defn: z |--> z + 1),
  ((z, 1, z, z, z, z + 1, z, z, z, z, z, z, z + 1, z, z, z,
    z, z + 1, z, z, z);
   (1,11,5,12,3,19)(2,8)(6,18,13)(7,17,15)(9,10,14,16,20,21),
   Ring endomorphism of Finite Field in z of size 2^2
     Defn: z |--> z + 1),
  ((z, z, z, z, z, z, z, z, z, z, z, z, z, z, z, z, z, z, z, z, z);
   (),
   Ring endomorphism of Finite Field in z of size 2^2
     Defn: z |--> z)],
 362880)
>>> C.automorphism_group_gens(equivalence="linear")
([((z, 1, z + 1, z + 1, 1, z + 1, z, 1, z + 1, z + 1, 1, z, 1, z + 1,
    z, 1, z, 1, z + 1, 1, 1);
   (1,12,11,10,6,8,9,20,13,21,5,14,3,16,17,19,7,4,2,15,18),
   Ring endomorphism of Finite Field in z of size 2^2
     Defn: z |--> z),
  ((z + 1, z + 1, z + 1, z, 1, 1, z, z, 1, z + 1, z, 1, 1, z, 1, z + 1,
    z, z + 1, z + 1, 1, z);
   (1,3,18,2,17,6,19)(4,15,13,20,7,14,16)(5,11,8,21,12,9,10),
   Ring endomorphism of Finite Field in z of size 2^2
     Defn: z |--> z),
  ((z + 1, z + 1, z + 1, z + 1, z + 1, z + 1, z + 1, z + 1, z + 1, z + 1,
    z + 1, z + 1, z + 1, z + 1, z + 1, z + 1, z + 1, z + 1, z + 1, z + 1, z + 1);
   (),
   Ring endomorphism of Finite Field in z of size 2^2
     Defn: z |--> z)],
 181440)
>>> C.automorphism_group_gens(equivalence="permutational")
([((1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1);
   (1,11)(3,10)(4,9)(5,7)(12,21)(14,20)(15,19)(16,17),
   Ring endomorphism of Finite Field in z of size 2^2
     Defn: z |--> z),
  ((1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1);
   (2,18)(3,19)(4,10)(5,16)(8,13)(9,14)(11,21)(15,20),
   Ring endomorphism of Finite Field in z of size 2^2
     Defn: z |--> z),
  ((1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1);
   (1,19)(3,17)(4,21)(5,20)(7,14)(9,12)(10,16)(11,15),
   Ring endomorphism of Finite Field in z of size 2^2
     Defn: z |--> z),
  ((1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1);
   (2,13)(3,14)(4,20)(5,11)(8,18)(9,19)(10,15)(16,21),
   Ring endomorphism of Finite Field in z of size 2^2
     Defn: z |--> z)],
64)
binomial_moment(i)[source]#

Return the i-th binomial moment of the \([n,k,d]_q\)-code \(C\):

\[B_i(C) = \sum_{S, |S|=i} \frac{q^{k_S}-1}{q-1}\]

where \(k_S\) is the dimension of the shortened code \(C_{J-S}\), \(J=[1,2,...,n]\). (The normalized binomial moment is \(b_i(C) = \binom{n}{d+i})^{-1}B_{d+i}(C)\).) In other words, \(C_{J-S}\) is isomorphic to the subcode of C of codewords supported on S.

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: C.binomial_moment(2)                                                  # needs sage.libs.gap
0
sage: C.binomial_moment(4)          # long time                             # needs sage.libs.gap
35
>>> from sage.all import *
>>> C = codes.HammingCode(GF(Integer(2)), Integer(3))
>>> C.binomial_moment(Integer(2))                                                  # needs sage.libs.gap
0
>>> C.binomial_moment(Integer(4))          # long time                             # needs sage.libs.gap
35

Warning

This is slow.

REFERENCE:

canonical_representative(equivalence='semilinear')[source]#

Compute a canonical orbit representative under the action of the semimonomial transformation group.

See sage.coding.codecan.autgroup_can_label for more details, for example if you would like to compute a canonical form under some more restrictive notion of equivalence, i.e. if you would like to restrict the permutation group to a Young subgroup.

INPUT:

  • equivalence (optional) – which defines the acting group, either

    • "permutational"

    • "linear"

    • "semilinear"

OUTPUT:

  • a canonical representative of self

  • a semimonomial transformation mapping self onto its representative

EXAMPLES:

sage: F.<z> = GF(4)
sage: C = codes.HammingCode(F, 3)
sage: CanRep, transp = C.canonical_representative()                         # needs sage.libs.gap
>>> from sage.all import *
>>> F = GF(Integer(4), names=('z',)); (z,) = F._first_ngens(1)
>>> C = codes.HammingCode(F, Integer(3))
>>> CanRep, transp = C.canonical_representative()                         # needs sage.libs.gap

Check that the transporter element is correct:

sage: LinearCode(transp*C.generator_matrix()) == CanRep                     # needs sage.libs.gap
True
>>> from sage.all import *
>>> LinearCode(transp*C.generator_matrix()) == CanRep                     # needs sage.libs.gap
True

Check if an equivalent code has the same canonical representative:

sage: f = F.hom([z**2])
sage: C_iso = LinearCode(C.generator_matrix().apply_map(f))
sage: CanRep_iso, _ = C_iso.canonical_representative()                      # needs sage.libs.gap
sage: CanRep_iso == CanRep                                                  # needs sage.libs.gap
True
>>> from sage.all import *
>>> f = F.hom([z**Integer(2)])
>>> C_iso = LinearCode(C.generator_matrix().apply_map(f))
>>> CanRep_iso, _ = C_iso.canonical_representative()                      # needs sage.libs.gap
>>> CanRep_iso == CanRep                                                  # needs sage.libs.gap
True

Since applying the Frobenius automorphism could be extended to an automorphism of \(C\), the following must also yield True:

sage: CanRep1, _ = C.canonical_representative("linear")                     # needs sage.libs.gap
sage: CanRep2, _ = C_iso.canonical_representative("linear")                 # needs sage.libs.gap
sage: CanRep2 == CanRep1                                                    # needs sage.libs.gap
True
>>> from sage.all import *
>>> CanRep1, _ = C.canonical_representative("linear")                     # needs sage.libs.gap
>>> CanRep2, _ = C_iso.canonical_representative("linear")                 # needs sage.libs.gap
>>> CanRep2 == CanRep1                                                    # needs sage.libs.gap
True
characteristic()[source]#

Return the characteristic of the base ring of self.

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: C.characteristic()
2
>>> from sage.all import *
>>> C = codes.HammingCode(GF(Integer(2)), Integer(3))
>>> C.characteristic()
2
characteristic_polynomial()[source]#

Return the characteristic polynomial of a linear code, as defined in [Lin1999].

EXAMPLES:

sage: C = codes.GolayCode(GF(2))
sage: C.characteristic_polynomial()
-4/3*x^3 + 64*x^2 - 2816/3*x + 4096
>>> from sage.all import *
>>> C = codes.GolayCode(GF(Integer(2)))
>>> C.characteristic_polynomial()
-4/3*x^3 + 64*x^2 - 2816/3*x + 4096
chinen_polynomial()[source]#

Return the Chinen zeta polynomial of the code.

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: C.chinen_polynomial()       # long time
1/5*(2*sqrt(2)*t^3 + 2*sqrt(2)*t^2 + 2*t^2 + sqrt(2)*t + 2*t + 1)/(sqrt(2) + 1)
sage: C = codes.GolayCode(GF(3), False)
sage: C.chinen_polynomial()       # long time
1/7*(3*sqrt(3)*t^3 + 3*sqrt(3)*t^2 + 3*t^2 + sqrt(3)*t + 3*t + 1)/(sqrt(3) + 1)
>>> from sage.all import *
>>> C = codes.HammingCode(GF(Integer(2)), Integer(3))
>>> C.chinen_polynomial()       # long time
1/5*(2*sqrt(2)*t^3 + 2*sqrt(2)*t^2 + 2*t^2 + sqrt(2)*t + 2*t + 1)/(sqrt(2) + 1)
>>> C = codes.GolayCode(GF(Integer(3)), False)
>>> C.chinen_polynomial()       # long time
1/7*(3*sqrt(3)*t^3 + 3*sqrt(3)*t^2 + 3*t^2 + sqrt(3)*t + 3*t + 1)/(sqrt(3) + 1)

This last output agrees with the corresponding example given in Chinen’s paper below.

REFERENCES:

  • Chinen, K. “An abundance of invariant polynomials satisfying the Riemann hypothesis”, April 2007 preprint.

construction_x(other, aux)[source]#

Construction X applied to self=C_1, other=C_2 and aux=C_a.

other must be a subcode of self.

If \(C_1\) is a \([n, k_1, d_1]\) linear code and \(C_2\) is a \([n, k_2, d_2]\) linear code, then \(k_1 > k_2\) and \(d_1 < d_2\). \(C_a\) must be a \([n_a, k_a, d_a]\) linear code, such that \(k_a + k_2 = k_1\) and \(d_a + d_1 \leq d_2\).

The method will then return a \([n+n_a, k_1, d_a+d_1]\) linear code.

EXAMPLES:

sage: C = codes.BCHCode(GF(2),15,7)
sage: C
[15, 5] BCH Code over GF(2) with designed distance 7
sage: D = codes.BCHCode(GF(2),15,5)
sage: D
[15, 7] BCH Code over GF(2) with designed distance 5
sage: C.is_subcode(D)
True

sage: # needs sage.libs.gap
sage: C.minimum_distance()
7
sage: D.minimum_distance()
5
sage: aux = codes.HammingCode(GF(2),2)
sage: aux = aux.dual_code()
sage: aux.minimum_distance()
2
sage: Cx = D.construction_x(C,aux)
sage: Cx
[18, 7] linear code over GF(2)
sage: Cx.minimum_distance()
7
>>> from sage.all import *
>>> C = codes.BCHCode(GF(Integer(2)),Integer(15),Integer(7))
>>> C
[15, 5] BCH Code over GF(2) with designed distance 7
>>> D = codes.BCHCode(GF(Integer(2)),Integer(15),Integer(5))
>>> D
[15, 7] BCH Code over GF(2) with designed distance 5
>>> C.is_subcode(D)
True

>>> # needs sage.libs.gap
>>> C.minimum_distance()
7
>>> D.minimum_distance()
5
>>> aux = codes.HammingCode(GF(Integer(2)),Integer(2))
>>> aux = aux.dual_code()
>>> aux.minimum_distance()
2
>>> Cx = D.construction_x(C,aux)
>>> Cx
[18, 7] linear code over GF(2)
>>> Cx.minimum_distance()
7
cosetGraph()[source]#

Return the coset graph of this linear code.

The coset graph of a linear code \(C\) is the graph whose vertices are the cosets of \(C\), considered as a subgroup of the additive group of the ambient vector space, and two cosets are adjacent if they have representatives that differ in exactly one coordinate.

EXAMPLES:

sage: # needs sage.graphs
sage: C = codes.GolayCode(GF(3))
sage: G = C.cosetGraph()
sage: G.is_distance_regular()
True
sage: C = codes.KasamiCode(8,2)
sage: G = C.cosetGraph()
sage: G.is_distance_regular()
True
>>> from sage.all import *
>>> # needs sage.graphs
>>> C = codes.GolayCode(GF(Integer(3)))
>>> G = C.cosetGraph()
>>> G.is_distance_regular()
True
>>> C = codes.KasamiCode(Integer(8),Integer(2))
>>> G = C.cosetGraph()
>>> G.is_distance_regular()
True

ALGORITHM:

Instead of working with cosets we compute a (direct sum) complement of \(C\). Let \(P\) be the projection of the cosets to the newly found subspace. Then two vectors are adjacent if they differ by \(\lambda P(e_i)\) for some \(i\).

covering_radius()[source]#

Return the minimal integer \(r\) such that any element in the ambient space of self has distance at most \(r\) to a codeword of self.

This method requires the optional GAP package Guava.

If the covering radius of a code equals its minimum distance, then the code is called perfect.

Note

This method is currently not implemented on codes over base fields of cardinality greater than 256 due to limitations in the underlying algorithm of GAP.

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 5)
sage: C.covering_radius()                           # optional - gap_package_guava
...1

sage: C = codes.random_linear_code(GF(263), 5, 1)
sage: C.covering_radius()                           # optional - gap_package_guava
Traceback (most recent call last):
...
NotImplementedError: the GAP algorithm that Sage is using
is limited to computing with fields of size at most 256
>>> from sage.all import *
>>> C = codes.HammingCode(GF(Integer(2)), Integer(5))
>>> C.covering_radius()                           # optional - gap_package_guava
...1

>>> C = codes.random_linear_code(GF(Integer(263)), Integer(5), Integer(1))
>>> C.covering_radius()                           # optional - gap_package_guava
Traceback (most recent call last):
...
NotImplementedError: the GAP algorithm that Sage is using
is limited to computing with fields of size at most 256
direct_sum(other)[source]#

Return the direct sum of the codes self and other.

This returns the code given by the direct sum of the codes self and other, which must be linear codes defined over the same base ring.

EXAMPLES:

sage: C1 = codes.HammingCode(GF(2), 3)
sage: C2 = C1.direct_sum(C1); C2
[14, 8] linear code over GF(2)
sage: C3 = C1.direct_sum(C2); C3
[21, 12] linear code over GF(2)
>>> from sage.all import *
>>> C1 = codes.HammingCode(GF(Integer(2)), Integer(3))
>>> C2 = C1.direct_sum(C1); C2
[14, 8] linear code over GF(2)
>>> C3 = C1.direct_sum(C2); C3
[21, 12] linear code over GF(2)
divisor()[source]#

Return the greatest common divisor of the weights of the nonzero codewords.

EXAMPLES:

sage: C = codes.GolayCode(GF(2))
sage: C.divisor()   # Type II self-dual
4
sage: C = codes.QuadraticResidueCodeEvenPair(17, GF(2))[0]
sage: C.divisor()
2
>>> from sage.all import *
>>> C = codes.GolayCode(GF(Integer(2)))
>>> C.divisor()   # Type II self-dual
4
>>> C = codes.QuadraticResidueCodeEvenPair(Integer(17), GF(Integer(2)))[Integer(0)]
>>> C.divisor()
2
extended_code()[source]#

Return self as an extended code.

See documentation of sage.coding.extended_code.ExtendedCode for details.

EXAMPLES:

sage: C = codes.HammingCode(GF(4,'a'), 3)
sage: C
[21, 18] Hamming Code over GF(4)
sage: Cx = C.extended_code()
sage: Cx
Extension of [21, 18] Hamming Code over GF(4)
>>> from sage.all import *
>>> C = codes.HammingCode(GF(Integer(4),'a'), Integer(3))
>>> C
[21, 18] Hamming Code over GF(4)
>>> Cx = C.extended_code()
>>> Cx
Extension of [21, 18] Hamming Code over GF(4)
galois_closure(F0)[source]#

If self is a linear code defined over \(F\) and \(F_0\) is a subfield with Galois group \(G = Gal(F/F_0)\) then this returns the \(G\)-module \(C^-\) containing \(C\).

EXAMPLES:

sage: C = codes.HammingCode(GF(4,'a'), 3)
sage: Cc = C.galois_closure(GF(2))
sage: C; Cc
[21, 18] Hamming Code over GF(4)
[21, 20] linear code over GF(4)
sage: c = C.basis()[2]
sage: V = VectorSpace(GF(4,'a'),21)
sage: c2 = V([x^2 for x in c.list()])
sage: c2 in C
False
sage: c2 in Cc
True
>>> from sage.all import *
>>> C = codes.HammingCode(GF(Integer(4),'a'), Integer(3))
>>> Cc = C.galois_closure(GF(Integer(2)))
>>> C; Cc
[21, 18] Hamming Code over GF(4)
[21, 20] linear code over GF(4)
>>> c = C.basis()[Integer(2)]
>>> V = VectorSpace(GF(Integer(4),'a'),Integer(21))
>>> c2 = V([x**Integer(2) for x in c.list()])
>>> c2 in C
False
>>> c2 in Cc
True
genus()[source]#

Return the “Duursma genus” of the code, \(\gamma_C = n+1-k-d\).

EXAMPLES:

sage: C1 = codes.HammingCode(GF(2), 3); C1
[7, 4] Hamming Code over GF(2)
sage: C1.genus()
1
sage: C2 = codes.HammingCode(GF(4,"a"), 2); C2
[5, 3] Hamming Code over GF(4)
sage: C2.genus()
0
>>> from sage.all import *
>>> C1 = codes.HammingCode(GF(Integer(2)), Integer(3)); C1
[7, 4] Hamming Code over GF(2)
>>> C1.genus()
1
>>> C2 = codes.HammingCode(GF(Integer(4),"a"), Integer(2)); C2
[5, 3] Hamming Code over GF(4)
>>> C2.genus()
0

Since all Hamming codes have minimum distance 3, these computations agree with the definition, \(n+1-k-d\).

is_galois_closed()[source]#

Checks if self is equal to its Galois closure.

EXAMPLES:

sage: C = codes.HammingCode(GF(4,"a"), 3)
sage: C.is_galois_closed()
False
>>> from sage.all import *
>>> C = codes.HammingCode(GF(Integer(4),"a"), Integer(3))
>>> C.is_galois_closed()
False
is_permutation_equivalent(other, algorithm=None)[source]#

Return True if self and other are permutation equivalent codes and False otherwise.

The algorithm="verbose" option also returns a permutation (if True) sending self to other.

Uses Robert Miller’s double coset partition refinement work.

EXAMPLES:

sage: P.<x> = PolynomialRing(GF(2),"x")
sage: g = x^3 + x + 1
sage: C1 = codes.CyclicCode(length=7, generator_pol=g); C1
[7, 4] Cyclic Code over GF(2)
sage: C2 = codes.HammingCode(GF(2), 3); C2
[7, 4] Hamming Code over GF(2)
sage: C1.is_permutation_equivalent(C2)
True
sage: C1.is_permutation_equivalent(C2, algorithm="verbose")                 # needs sage.groups
(True, (3,4)(5,7,6))
sage: C1 = codes.random_linear_code(GF(2), 10, 5)
sage: C2 = codes.random_linear_code(GF(3), 10, 5)
sage: C1.is_permutation_equivalent(C2)
False
>>> from sage.all import *
>>> P = PolynomialRing(GF(Integer(2)),"x", names=('x',)); (x,) = P._first_ngens(1)
>>> g = x**Integer(3) + x + Integer(1)
>>> C1 = codes.CyclicCode(length=Integer(7), generator_pol=g); C1
[7, 4] Cyclic Code over GF(2)
>>> C2 = codes.HammingCode(GF(Integer(2)), Integer(3)); C2
[7, 4] Hamming Code over GF(2)
>>> C1.is_permutation_equivalent(C2)
True
>>> C1.is_permutation_equivalent(C2, algorithm="verbose")                 # needs sage.groups
(True, (3,4)(5,7,6))
>>> C1 = codes.random_linear_code(GF(Integer(2)), Integer(10), Integer(5))
>>> C2 = codes.random_linear_code(GF(Integer(3)), Integer(10), Integer(5))
>>> C1.is_permutation_equivalent(C2)
False
is_projective()[source]#

Test whether the code is projective.

A linear code \(C\) over a field is called projective when its dual \(Cd\) has minimum weight \(\geq 3\), i.e. when no two coordinate positions of \(C\) are linearly independent (cf. definition 3 from [BS2011] or 9.8.1 from [BH2012]).

EXAMPLES:

sage: C = codes.GolayCode(GF(2), False)
sage: C.is_projective()
True
sage: C.dual_code().minimum_distance()                                      # needs sage.libs.gap
8
>>> from sage.all import *
>>> C = codes.GolayCode(GF(Integer(2)), False)
>>> C.is_projective()
True
>>> C.dual_code().minimum_distance()                                      # needs sage.libs.gap
8

A non-projective code:

sage: C = codes.LinearCode(matrix(GF(2), [[1,0,1],[1,1,1]]))
sage: C.is_projective()
False
>>> from sage.all import *
>>> C = codes.LinearCode(matrix(GF(Integer(2)), [[Integer(1),Integer(0),Integer(1)],[Integer(1),Integer(1),Integer(1)]]))
>>> C.is_projective()
False
juxtapose(other)[source]#

Juxtaposition of self and other

The two codes must have equal dimension.

EXAMPLES:

sage: C1 = codes.HammingCode(GF(2), 3)
sage: C2 = C1.juxtapose(C1)
sage: C2
[14, 4] linear code over GF(2)
>>> from sage.all import *
>>> C1 = codes.HammingCode(GF(Integer(2)), Integer(3))
>>> C2 = C1.juxtapose(C1)
>>> C2
[14, 4] linear code over GF(2)
minimum_distance(algorithm=None)[source]#

Return the minimum distance of self.

Note

When using GAP, this raises a NotImplementedError if the base field of the code has size greater than 256 due to limitations in GAP.

INPUT:

  • algorithm – (default: None) the name of the algorithm to use to perform minimum distance computation. algorithm can be:

    • None, to use GAP methods (but not Guava)

    • "Guava", to use the optional GAP package Guava

OUTPUT:

  • Integer, minimum distance of this code

EXAMPLES:

sage: MS = MatrixSpace(GF(3),4,7)
sage: G = MS([[1,1,1,0,0,0,0], [1,0,0,1,1,0,0], [0,1,0,1,0,1,0], [1,1,0,1,0,0,1]])
sage: C = LinearCode(G)
sage: C.minimum_distance()                                                  # needs sage.libs.gap
3
>>> from sage.all import *
>>> MS = MatrixSpace(GF(Integer(3)),Integer(4),Integer(7))
>>> G = MS([[Integer(1),Integer(1),Integer(1),Integer(0),Integer(0),Integer(0),Integer(0)], [Integer(1),Integer(0),Integer(0),Integer(1),Integer(1),Integer(0),Integer(0)], [Integer(0),Integer(1),Integer(0),Integer(1),Integer(0),Integer(1),Integer(0)], [Integer(1),Integer(1),Integer(0),Integer(1),Integer(0),Integer(0),Integer(1)]])
>>> C = LinearCode(G)
>>> C.minimum_distance()                                                  # needs sage.libs.gap
3

If algorithm is provided, then the minimum distance will be recomputed even if there is a stored value from a previous run.:

sage: C.minimum_distance(algorithm="gap")                                   # needs sage.libs.gap
3
sage: libgap.SetAllInfoLevels(0)         # to suppress extra info messages  # needs sage.libs.gap
sage: C.minimum_distance(algorithm="guava")         # optional - gap_package_guava
...3
>>> from sage.all import *
>>> C.minimum_distance(algorithm="gap")                                   # needs sage.libs.gap
3
>>> libgap.SetAllInfoLevels(Integer(0))         # to suppress extra info messages  # needs sage.libs.gap
>>> C.minimum_distance(algorithm="guava")         # optional - gap_package_guava
...3
module_composition_factors(gp)[source]#

Print the GAP record of the Meataxe composition factors module.

This is displayed in Meataxe notation.

This uses GAP but not Guava.

EXAMPLES:

sage: MS = MatrixSpace(GF(2),4,8)
sage: G  = MS([[1,0,0,0,1,1,1,0], [0,1,1,1,0,0,0,0],
....:          [0,0,0,0,0,0,0,1], [0,0,0,0,0,1,0,0]])
sage: C  = LinearCode(G)
sage: gp = C.permutation_automorphism_group()                               # needs sage.libs.gap
sage: C.module_composition_factors(gp)                                      # needs sage.libs.gap
[ rec(
      IsIrreducible := true,
      IsOverFiniteField := true,
...) ]
>>> from sage.all import *
>>> MS = MatrixSpace(GF(Integer(2)),Integer(4),Integer(8))
>>> G  = MS([[Integer(1),Integer(0),Integer(0),Integer(0),Integer(1),Integer(1),Integer(1),Integer(0)], [Integer(0),Integer(1),Integer(1),Integer(1),Integer(0),Integer(0),Integer(0),Integer(0)],
...          [Integer(0),Integer(0),Integer(0),Integer(0),Integer(0),Integer(0),Integer(0),Integer(1)], [Integer(0),Integer(0),Integer(0),Integer(0),Integer(0),Integer(1),Integer(0),Integer(0)]])
>>> C  = LinearCode(G)
>>> gp = C.permutation_automorphism_group()                               # needs sage.libs.gap
>>> C.module_composition_factors(gp)                                      # needs sage.libs.gap
[ rec(
      IsIrreducible := true,
      IsOverFiniteField := true,
...) ]
permutation_automorphism_group(algorithm='partition')[source]#

If \(C\) is an \([n,k,d]\) code over \(F\), this function computes the subgroup \(Aut(C) \subset S_n\) of all permutation automorphisms of \(C\). The binary case always uses the (default) partition refinement algorithm of Robert Miller.

Note that if the base ring of \(C\) is \(GF(2)\) then this is the full automorphism group. Otherwise, you could use automorphism_group_gens() to compute generators of the full automorphism group.

INPUT:

  • algorithm – If "gap" then GAP’s MatrixAutomorphism function (written by Thomas Breuer) is used. The implementation combines an idea of mine with an improvement suggested by Cary Huffman. If "gap+verbose" then code-theoretic data is printed out at several stages of the computation. If "partition" then the (default) partition refinement algorithm of Robert Miller is used. Finally, if "codecan" then the partition refinement algorithm of Thomas Feulner is used, which also computes a canonical representative of self (call canonical_representative() to access it).

OUTPUT:

  • Permutation automorphism group

EXAMPLES:

sage: MS = MatrixSpace(GF(2),4,8)
sage: G  = MS([[1,0,0,0,1,1,1,0], [0,1,1,1,0,0,0,0],
....:          [0,0,0,0,0,0,0,1], [0,0,0,0,0,1,0,0]])
sage: C  = LinearCode(G); C
[8, 4] linear code over GF(2)

sage: # needs sage.groups
sage: G = C.permutation_automorphism_group()
sage: G.order()
144
sage: GG = C.permutation_automorphism_group("codecan")
sage: GG == G
True
>>> from sage.all import *
>>> MS = MatrixSpace(GF(Integer(2)),Integer(4),Integer(8))
>>> G  = MS([[Integer(1),Integer(0),Integer(0),Integer(0),Integer(1),Integer(1),Integer(1),Integer(0)], [Integer(0),Integer(1),Integer(1),Integer(1),Integer(0),Integer(0),Integer(0),Integer(0)],
...          [Integer(0),Integer(0),Integer(0),Integer(0),Integer(0),Integer(0),Integer(0),Integer(1)], [Integer(0),Integer(0),Integer(0),Integer(0),Integer(0),Integer(1),Integer(0),Integer(0)]])
>>> C  = LinearCode(G); C
[8, 4] linear code over GF(2)

>>> # needs sage.groups
>>> G = C.permutation_automorphism_group()
>>> G.order()
144
>>> GG = C.permutation_automorphism_group("codecan")
>>> GG == G
True

A less easy example involves showing that the permutation automorphism group of the extended ternary Golay code is the Mathieu group \(M_{11}\).

sage: # needs sage.groups
sage: C = codes.GolayCode(GF(3))
sage: M11 = MathieuGroup(11)
sage: M11.order()
7920
sage: G = C.permutation_automorphism_group()                # long time (6s on sage.math, 2011)
sage: G.is_isomorphic(M11)                                  # long time
True
sage: GG = C.permutation_automorphism_group("codecan")      # long time
sage: GG == G                                               # long time
True
>>> from sage.all import *
>>> # needs sage.groups
>>> C = codes.GolayCode(GF(Integer(3)))
>>> M11 = MathieuGroup(Integer(11))
>>> M11.order()
7920
>>> G = C.permutation_automorphism_group()                # long time (6s on sage.math, 2011)
>>> G.is_isomorphic(M11)                                  # long time
True
>>> GG = C.permutation_automorphism_group("codecan")      # long time
>>> GG == G                                               # long time
True

Other examples:

sage: # needs sage.groups
sage: C = codes.GolayCode(GF(2))
sage: G = C.permutation_automorphism_group()
sage: G.order()
244823040
sage: C = codes.HammingCode(GF(2), 5)
sage: G = C.permutation_automorphism_group()
sage: G.order()
9999360
sage: C = codes.HammingCode(GF(3), 2); C
[4, 2] Hamming Code over GF(3)
sage: C.permutation_automorphism_group(algorithm="partition")
Permutation Group with generators [(1,3,4)]
sage: C = codes.HammingCode(GF(4,"z"), 2); C
[5, 3] Hamming Code over GF(4)
sage: G = C.permutation_automorphism_group(algorithm="partition"); G
Permutation Group with generators [(1,3)(4,5), (1,4)(3,5)]
sage: GG = C.permutation_automorphism_group(algorithm="codecan")    # long time
sage: GG == G                                                       # long time
True
sage: C.permutation_automorphism_group(algorithm="gap")     # optional - gap_package_guava
Permutation Group with generators [(1,3)(4,5), (1,4)(3,5)]
sage: C = codes.GolayCode(GF(3), True)
sage: C.permutation_automorphism_group(algorithm="gap")     # optional - gap_package_guava
Permutation Group with generators
 [(5,7)(6,11)(8,9)(10,12), (4,6,11)(5,8,12)(7,10,9), (3,4)(6,8)(9,11)(10,12),
  (2,3)(6,11)(8,12)(9,10), (1,2)(5,10)(7,12)(8,9)]
>>> from sage.all import *
>>> # needs sage.groups
>>> C = codes.GolayCode(GF(Integer(2)))
>>> G = C.permutation_automorphism_group()
>>> G.order()
244823040
>>> C = codes.HammingCode(GF(Integer(2)), Integer(5))
>>> G = C.permutation_automorphism_group()
>>> G.order()
9999360
>>> C = codes.HammingCode(GF(Integer(3)), Integer(2)); C
[4, 2] Hamming Code over GF(3)
>>> C.permutation_automorphism_group(algorithm="partition")
Permutation Group with generators [(1,3,4)]
>>> C = codes.HammingCode(GF(Integer(4),"z"), Integer(2)); C
[5, 3] Hamming Code over GF(4)
>>> G = C.permutation_automorphism_group(algorithm="partition"); G
Permutation Group with generators [(1,3)(4,5), (1,4)(3,5)]
>>> GG = C.permutation_automorphism_group(algorithm="codecan")    # long time
>>> GG == G                                                       # long time
True
>>> C.permutation_automorphism_group(algorithm="gap")     # optional - gap_package_guava
Permutation Group with generators [(1,3)(4,5), (1,4)(3,5)]
>>> C = codes.GolayCode(GF(Integer(3)), True)
>>> C.permutation_automorphism_group(algorithm="gap")     # optional - gap_package_guava
Permutation Group with generators
 [(5,7)(6,11)(8,9)(10,12), (4,6,11)(5,8,12)(7,10,9), (3,4)(6,8)(9,11)(10,12),
  (2,3)(6,11)(8,12)(9,10), (1,2)(5,10)(7,12)(8,9)]

However, the option algorithm="gap+verbose", will print out:

Minimum distance: 5 Weight distribution: [1, 0, 0, 0, 0, 132, 132,
0, 330, 110, 0, 24]

Using the 132 codewords of weight 5 Supergroup size: 39916800

in addition to the output of C.permutation_automorphism_group(algorithm="gap").

product_code(other)[source]#

Combines self with other to give the tensor product code.

If self is a \([n_1, k_1, d_1]\)-code and other is a \([n_2, k_2, d_2]\)-code, the product is a \([n_1n_2, k_1k_2, d_1d_2]\)-code.

Note that the two codes have to be over the same field.

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: C
[7, 4] Hamming Code over GF(2)
sage: D = codes.ReedMullerCode(GF(2), 2, 2)
sage: D
Binary Reed-Muller Code of order 2 and number of variables 2
sage: A = C.product_code(D)
sage: A
[28, 16] linear code over GF(2)
sage: A.length() == C.length()*D.length()
True
sage: A.dimension() == C.dimension()*D.dimension()
True
sage: A.minimum_distance() == C.minimum_distance()*D.minimum_distance()     # needs sage.libs.gap
True
>>> from sage.all import *
>>> C = codes.HammingCode(GF(Integer(2)), Integer(3))
>>> C
[7, 4] Hamming Code over GF(2)
>>> D = codes.ReedMullerCode(GF(Integer(2)), Integer(2), Integer(2))
>>> D
Binary Reed-Muller Code of order 2 and number of variables 2
>>> A = C.product_code(D)
>>> A
[28, 16] linear code over GF(2)
>>> A.length() == C.length()*D.length()
True
>>> A.dimension() == C.dimension()*D.dimension()
True
>>> A.minimum_distance() == C.minimum_distance()*D.minimum_distance()     # needs sage.libs.gap
True
punctured(L)[source]#

Return a sage.coding.punctured_code object from L.

INPUT:

  • L – List of positions to puncture

OUTPUT:

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: C.punctured([1,2])
Puncturing of [7, 4] Hamming Code over GF(2) on position(s) [1, 2]
>>> from sage.all import *
>>> C = codes.HammingCode(GF(Integer(2)), Integer(3))
>>> C.punctured([Integer(1),Integer(2)])
Puncturing of [7, 4] Hamming Code over GF(2) on position(s) [1, 2]
relative_distance()[source]#

Return the ratio of the minimum distance to the code length.

EXAMPLES:

sage: C = codes.HammingCode(GF(2),3)
sage: C.relative_distance()
3/7
>>> from sage.all import *
>>> C = codes.HammingCode(GF(Integer(2)),Integer(3))
>>> C.relative_distance()
3/7
shortened(L)[source]#

Return the code shortened at the positions L, where \(L \subset \{1,2,...,n\}\).

Consider the subcode \(C(L)\) consisting of all codewords \(c\in C\) which satisfy \(c_i=0\) for all \(i\in L\). The punctured code \(C(L)^L\) is called the shortened code on \(L\) and is denoted \(C_L\). The code constructed is actually only isomorphic to the shortened code defined in this way.

By Theorem 1.5.7 in [HP2003], \(C_L\) is \(((C^\perp)^L)^\perp\). This is used in the construction below.

INPUT:

  • L – Subset of \(\{1,...,n\}\), where \(n\) is the length of this code

OUTPUT:

  • Linear code, the shortened code described above

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: C.shortened([1,2])
[5, 2] linear code over GF(2)
>>> from sage.all import *
>>> C = codes.HammingCode(GF(Integer(2)), Integer(3))
>>> C.shortened([Integer(1),Integer(2)])
[5, 2] linear code over GF(2)
spectrum(algorithm=None)[source]#

Return the weight distribution, or spectrum, of self as a list.

The weight distribution a code of length \(n\) is the sequence \(A_0, A_1,..., A_n\) where \(A_i\) is the number of codewords of weight \(i\).

INPUT:

  • algorithm – (default: None) If set to "gap", call GAP. If set to "leon", call the option GAP package GUAVA and call a function therein by Jeffrey Leon (see warning below). If set to "binary", use an algorithm optimized for binary codes. The default is to use "binary" for binary codes and "gap" otherwise.

OUTPUT:

  • A list of non-negative integers: the weight distribution.

Warning

Specifying algorithm="leon" sometimes prints a traceback related to a stack smashing error in the C library. The result appears to be computed correctly, however. It appears to run much faster than the GAP algorithm in small examples and much slower than the GAP algorithm in larger examples.

EXAMPLES:

sage: MS = MatrixSpace(GF(2),4,7)
sage: G = MS([[1,1,1,0,0,0,0],[1,0,0,1,1,0,0],[0,1,0,1,0,1,0],[1,1,0,1,0,0,1]])
sage: C = LinearCode(G)
sage: C.weight_distribution()
[1, 0, 0, 7, 7, 0, 0, 1]
sage: F.<z> = GF(2^2,"z")
sage: C = codes.HammingCode(F, 2); C
[5, 3] Hamming Code over GF(4)
sage: C.weight_distribution()                                               # needs sage.libs.gap
[1, 0, 0, 30, 15, 18]
sage: C = codes.HammingCode(GF(2), 3); C
[7, 4] Hamming Code over GF(2)
sage: C.weight_distribution(algorithm="leon")   # optional - gap_package_guava
[1, 0, 0, 7, 7, 0, 0, 1]
sage: C.weight_distribution(algorithm="gap")                                # needs sage.libs.gap
[1, 0, 0, 7, 7, 0, 0, 1]
sage: C.weight_distribution(algorithm="binary")
[1, 0, 0, 7, 7, 0, 0, 1]

sage: # optional - gap_package_guava
sage: C = codes.HammingCode(GF(3), 3); C
[13, 10] Hamming Code over GF(3)
sage: C.weight_distribution() == C.weight_distribution(algorithm="leon")
True
sage: C = codes.HammingCode(GF(5), 2); C
[6, 4] Hamming Code over GF(5)
sage: C.weight_distribution() == C.weight_distribution(algorithm="leon")
True
sage: C = codes.HammingCode(GF(7), 2); C
[8, 6] Hamming Code over GF(7)
sage: C.weight_distribution() == C.weight_distribution(algorithm="leon")
True
>>> from sage.all import *
>>> MS = MatrixSpace(GF(Integer(2)),Integer(4),Integer(7))
>>> G = MS([[Integer(1),Integer(1),Integer(1),Integer(0),Integer(0),Integer(0),Integer(0)],[Integer(1),Integer(0),Integer(0),Integer(1),Integer(1),Integer(0),Integer(0)],[Integer(0),Integer(1),Integer(0),Integer(1),Integer(0),Integer(1),Integer(0)],[Integer(1),Integer(1),Integer(0),Integer(1),Integer(0),Integer(0),Integer(1)]])
>>> C = LinearCode(G)
>>> C.weight_distribution()
[1, 0, 0, 7, 7, 0, 0, 1]
>>> F = GF(Integer(2)**Integer(2),"z", names=('z',)); (z,) = F._first_ngens(1)
>>> C = codes.HammingCode(F, Integer(2)); C
[5, 3] Hamming Code over GF(4)
>>> C.weight_distribution()                                               # needs sage.libs.gap
[1, 0, 0, 30, 15, 18]
>>> C = codes.HammingCode(GF(Integer(2)), Integer(3)); C
[7, 4] Hamming Code over GF(2)
>>> C.weight_distribution(algorithm="leon")   # optional - gap_package_guava
[1, 0, 0, 7, 7, 0, 0, 1]
>>> C.weight_distribution(algorithm="gap")                                # needs sage.libs.gap
[1, 0, 0, 7, 7, 0, 0, 1]
>>> C.weight_distribution(algorithm="binary")
[1, 0, 0, 7, 7, 0, 0, 1]

>>> # optional - gap_package_guava
>>> C = codes.HammingCode(GF(Integer(3)), Integer(3)); C
[13, 10] Hamming Code over GF(3)
>>> C.weight_distribution() == C.weight_distribution(algorithm="leon")
True
>>> C = codes.HammingCode(GF(Integer(5)), Integer(2)); C
[6, 4] Hamming Code over GF(5)
>>> C.weight_distribution() == C.weight_distribution(algorithm="leon")
True
>>> C = codes.HammingCode(GF(Integer(7)), Integer(2)); C
[8, 6] Hamming Code over GF(7)
>>> C.weight_distribution() == C.weight_distribution(algorithm="leon")
True
support()[source]#

Return the set of indices \(j\) where \(A_j\) is nonzero, where \(A_j\) is the number of codewords in self of Hamming weight \(j\).

OUTPUT:

  • List of integers

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: C.weight_distribution()
[1, 0, 0, 7, 7, 0, 0, 1]
sage: C.support()
[0, 3, 4, 7]
>>> from sage.all import *
>>> C = codes.HammingCode(GF(Integer(2)), Integer(3))
>>> C.weight_distribution()
[1, 0, 0, 7, 7, 0, 0, 1]
>>> C.support()
[0, 3, 4, 7]
u_u_plus_v_code(other)[source]#

Return the \((u|u+v)\)-construction with self=u and other=v.

This returns the code obtained through \((u|u+v)\)-construction with self as \(u\) and other as \(v\). Note that \(u\) and \(v\) must have equal lengths. For \(u\) a \([n, k_1, d_1]\)-code and \(v\) a \([n, k_2, d_2]\)-code this returns a \([2n, k_1+k_2, d]\)-code, where \(d=\min(2d_1,d_2)\).

EXAMPLES:

sage: C1 = codes.HammingCode(GF(2), 3)
sage: C2 = codes.HammingCode(GF(2), 3)
sage: D = C1.u_u_plus_v_code(C2)
sage: D
[14, 8] linear code over GF(2)
>>> from sage.all import *
>>> C1 = codes.HammingCode(GF(Integer(2)), Integer(3))
>>> C2 = codes.HammingCode(GF(Integer(2)), Integer(3))
>>> D = C1.u_u_plus_v_code(C2)
>>> D
[14, 8] linear code over GF(2)
weight_distribution(algorithm=None)[source]#

Return the weight distribution, or spectrum, of self as a list.

The weight distribution a code of length \(n\) is the sequence \(A_0, A_1,..., A_n\) where \(A_i\) is the number of codewords of weight \(i\).

INPUT:

  • algorithm – (default: None) If set to "gap", call GAP. If set to "leon", call the option GAP package GUAVA and call a function therein by Jeffrey Leon (see warning below). If set to "binary", use an algorithm optimized for binary codes. The default is to use "binary" for binary codes and "gap" otherwise.

OUTPUT:

  • A list of non-negative integers: the weight distribution.

Warning

Specifying algorithm="leon" sometimes prints a traceback related to a stack smashing error in the C library. The result appears to be computed correctly, however. It appears to run much faster than the GAP algorithm in small examples and much slower than the GAP algorithm in larger examples.

EXAMPLES:

sage: MS = MatrixSpace(GF(2),4,7)
sage: G = MS([[1,1,1,0,0,0,0],[1,0,0,1,1,0,0],[0,1,0,1,0,1,0],[1,1,0,1,0,0,1]])
sage: C = LinearCode(G)
sage: C.weight_distribution()
[1, 0, 0, 7, 7, 0, 0, 1]
sage: F.<z> = GF(2^2,"z")
sage: C = codes.HammingCode(F, 2); C
[5, 3] Hamming Code over GF(4)
sage: C.weight_distribution()                                               # needs sage.libs.gap
[1, 0, 0, 30, 15, 18]
sage: C = codes.HammingCode(GF(2), 3); C
[7, 4] Hamming Code over GF(2)
sage: C.weight_distribution(algorithm="leon")   # optional - gap_package_guava
[1, 0, 0, 7, 7, 0, 0, 1]
sage: C.weight_distribution(algorithm="gap")                                # needs sage.libs.gap
[1, 0, 0, 7, 7, 0, 0, 1]
sage: C.weight_distribution(algorithm="binary")
[1, 0, 0, 7, 7, 0, 0, 1]

sage: # optional - gap_package_guava
sage: C = codes.HammingCode(GF(3), 3); C
[13, 10] Hamming Code over GF(3)
sage: C.weight_distribution() == C.weight_distribution(algorithm="leon")
True
sage: C = codes.HammingCode(GF(5), 2); C
[6, 4] Hamming Code over GF(5)
sage: C.weight_distribution() == C.weight_distribution(algorithm="leon")
True
sage: C = codes.HammingCode(GF(7), 2); C
[8, 6] Hamming Code over GF(7)
sage: C.weight_distribution() == C.weight_distribution(algorithm="leon")
True
>>> from sage.all import *
>>> MS = MatrixSpace(GF(Integer(2)),Integer(4),Integer(7))
>>> G = MS([[Integer(1),Integer(1),Integer(1),Integer(0),Integer(0),Integer(0),Integer(0)],[Integer(1),Integer(0),Integer(0),Integer(1),Integer(1),Integer(0),Integer(0)],[Integer(0),Integer(1),Integer(0),Integer(1),Integer(0),Integer(1),Integer(0)],[Integer(1),Integer(1),Integer(0),Integer(1),Integer(0),Integer(0),Integer(1)]])
>>> C = LinearCode(G)
>>> C.weight_distribution()
[1, 0, 0, 7, 7, 0, 0, 1]
>>> F = GF(Integer(2)**Integer(2),"z", names=('z',)); (z,) = F._first_ngens(1)
>>> C = codes.HammingCode(F, Integer(2)); C
[5, 3] Hamming Code over GF(4)
>>> C.weight_distribution()                                               # needs sage.libs.gap
[1, 0, 0, 30, 15, 18]
>>> C = codes.HammingCode(GF(Integer(2)), Integer(3)); C
[7, 4] Hamming Code over GF(2)
>>> C.weight_distribution(algorithm="leon")   # optional - gap_package_guava
[1, 0, 0, 7, 7, 0, 0, 1]
>>> C.weight_distribution(algorithm="gap")                                # needs sage.libs.gap
[1, 0, 0, 7, 7, 0, 0, 1]
>>> C.weight_distribution(algorithm="binary")
[1, 0, 0, 7, 7, 0, 0, 1]

>>> # optional - gap_package_guava
>>> C = codes.HammingCode(GF(Integer(3)), Integer(3)); C
[13, 10] Hamming Code over GF(3)
>>> C.weight_distribution() == C.weight_distribution(algorithm="leon")
True
>>> C = codes.HammingCode(GF(Integer(5)), Integer(2)); C
[6, 4] Hamming Code over GF(5)
>>> C.weight_distribution() == C.weight_distribution(algorithm="leon")
True
>>> C = codes.HammingCode(GF(Integer(7)), Integer(2)); C
[8, 6] Hamming Code over GF(7)
>>> C.weight_distribution() == C.weight_distribution(algorithm="leon")
True
weight_enumerator(names=None, bivariate=True)[source]#

Return the weight enumerator polynomial of self.

This is the bivariate, homogeneous polynomial in \(x\) and \(y\) whose coefficient to \(x^i y^{n-i}\) is the number of codewords of self of Hamming weight \(i\). Here, \(n\) is the length of self.

INPUT:

  • names – (default: "xy") The names of the variables in the homogeneous polynomial. Can be given as a single string of length 2, or a single string with a comma, or as a tuple or list of two strings.

  • bivariate – (default: True) Whether to return a bivariate, homogeneous polynomial or just a univariate polynomial. If set to False, then names will be interpreted as a single variable name and default to "x".

OUTPUT:

  • The weight enumerator polynomial over \(\ZZ\).

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: C.weight_enumerator()
x^7 + 7*x^4*y^3 + 7*x^3*y^4 + y^7
sage: C.weight_enumerator(names="st")
s^7 + 7*s^4*t^3 + 7*s^3*t^4 + t^7
sage: C.weight_enumerator(names="var1, var2")
var1^7 + 7*var1^4*var2^3 + 7*var1^3*var2^4 + var2^7
sage: C.weight_enumerator(names=('var1', 'var2'))
var1^7 + 7*var1^4*var2^3 + 7*var1^3*var2^4 + var2^7
sage: C.weight_enumerator(bivariate=False)
x^7 + 7*x^4 + 7*x^3 + 1
>>> from sage.all import *
>>> C = codes.HammingCode(GF(Integer(2)), Integer(3))
>>> C.weight_enumerator()
x^7 + 7*x^4*y^3 + 7*x^3*y^4 + y^7
>>> C.weight_enumerator(names="st")
s^7 + 7*s^4*t^3 + 7*s^3*t^4 + t^7
>>> C.weight_enumerator(names="var1, var2")
var1^7 + 7*var1^4*var2^3 + 7*var1^3*var2^4 + var2^7
>>> C.weight_enumerator(names=('var1', 'var2'))
var1^7 + 7*var1^4*var2^3 + 7*var1^3*var2^4 + var2^7
>>> C.weight_enumerator(bivariate=False)
x^7 + 7*x^4 + 7*x^3 + 1

An example of a code with a non-symmetrical weight enumerator:

sage: C = codes.GolayCode(GF(3), extended=False)
sage: C.weight_enumerator()
24*x^11 + 110*x^9*y^2 + 330*x^8*y^3 + 132*x^6*y^5 + 132*x^5*y^6 + y^11
>>> from sage.all import *
>>> C = codes.GolayCode(GF(Integer(3)), extended=False)
>>> C.weight_enumerator()
24*x^11 + 110*x^9*y^2 + 330*x^8*y^3 + 132*x^6*y^5 + 132*x^5*y^6 + y^11
zeta_function(name='T')[source]#

Return the Duursma zeta function of the code.

INPUT:

  • name – String, variable name (default: "T")

OUTPUT:

Element of \(\QQ(T)\)

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: C.zeta_function()                                                     # needs sage.libs.gap
(1/5*T^2 + 1/5*T + 1/10)/(T^2 - 3/2*T + 1/2)
>>> from sage.all import *
>>> C = codes.HammingCode(GF(Integer(2)), Integer(3))
>>> C.zeta_function()                                                     # needs sage.libs.gap
(1/5*T^2 + 1/5*T + 1/10)/(T^2 - 3/2*T + 1/2)
zeta_polynomial(name='T')[source]#

Return the Duursma zeta polynomial of this code.

Assumes that the minimum distances of this code and its dual are greater than 1. Prints a warning to stdout otherwise.

INPUT:

  • name – String, variable name (default: "T")

OUTPUT:

  • Polynomial over \(\QQ\)

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: C.zeta_polynomial()                                                   # needs sage.libs.gap
2/5*T^2 + 2/5*T + 1/5

sage: C = codes.databases.best_linear_code_in_guava(6, 3, GF(2))    # optional - gap_package_guava
sage: C.minimum_distance()                                          # optional - gap_package_guava
3
sage: C.zeta_polynomial()                                           # optional - gap_package_guava
2/5*T^2 + 2/5*T + 1/5

sage: C = codes.HammingCode(GF(2), 4)
sage: C.zeta_polynomial()                                                   # needs sage.libs.gap
16/429*T^6 + 16/143*T^5 + 80/429*T^4 + 32/143*T^3 + 30/143*T^2 + 2/13*T + 1/13

sage: F.<z> = GF(4,"z")
sage: MS = MatrixSpace(F, 3, 6)
sage: G = MS([[1,0,0,1,z,z],[0,1,0,z,1,z],[0,0,1,z,z,1]])
sage: C = LinearCode(G)  # the "hexacode"
sage: C.zeta_polynomial()                                                   # needs sage.libs.gap
1
>>> from sage.all import *
>>> C = codes.HammingCode(GF(Integer(2)), Integer(3))
>>> C.zeta_polynomial()                                                   # needs sage.libs.gap
2/5*T^2 + 2/5*T + 1/5

>>> C = codes.databases.best_linear_code_in_guava(Integer(6), Integer(3), GF(Integer(2)))    # optional - gap_package_guava
>>> C.minimum_distance()                                          # optional - gap_package_guava
3
>>> C.zeta_polynomial()                                           # optional - gap_package_guava
2/5*T^2 + 2/5*T + 1/5

>>> C = codes.HammingCode(GF(Integer(2)), Integer(4))
>>> C.zeta_polynomial()                                                   # needs sage.libs.gap
16/429*T^6 + 16/143*T^5 + 80/429*T^4 + 32/143*T^3 + 30/143*T^2 + 2/13*T + 1/13

>>> F = GF(Integer(4),"z", names=('z',)); (z,) = F._first_ngens(1)
>>> MS = MatrixSpace(F, Integer(3), Integer(6))
>>> G = MS([[Integer(1),Integer(0),Integer(0),Integer(1),z,z],[Integer(0),Integer(1),Integer(0),z,Integer(1),z],[Integer(0),Integer(0),Integer(1),z,z,Integer(1)]])
>>> C = LinearCode(G)  # the "hexacode"
>>> C.zeta_polynomial()                                                   # needs sage.libs.gap
1

REFERENCES:

class sage.coding.linear_code.LinearCode(generator, d=None)[source]#

Bases: AbstractLinearCode

Linear codes over a finite field or finite ring, represented using a generator matrix.

This class should be used for arbitrary and unstructured linear codes. This means that basic operations on the code, such as the computation of the minimum distance, will use generic, slow algorithms.

If you are looking for constructing a code from a more specific family, see if the family has been implemented by investigating codes.<tab>. These more specific classes use properties particular to that family to allow faster algorithms, and could also have family-specific methods.

See Wikipedia article Linear_code for more information on unstructured linear codes.

INPUT:

  • generator – a generator matrix over a finite field (G can be defined over a finite ring but the matrices over that ring must have certain attributes, such as rank); or a code over a finite field

  • d – (default: None) the minimum distance of the code

Note

The veracity of the minimum distance d, if provided, is not checked.

EXAMPLES:

sage: MS = MatrixSpace(GF(2),4,7)
sage: G  = MS([[1,1,1,0,0,0,0], [1,0,0,1,1,0,0], [0,1,0,1,0,1,0], [1,1,0,1,0,0,1]])
sage: C  = LinearCode(G)
sage: C
[7, 4] linear code over GF(2)
sage: C.base_ring()
Finite Field of size 2
sage: C.dimension()
4
sage: C.length()
7
sage: C.minimum_distance()                                                      # needs sage.libs.gap
3
sage: C.spectrum()
[1, 0, 0, 7, 7, 0, 0, 1]
sage: C.weight_distribution()
[1, 0, 0, 7, 7, 0, 0, 1]
>>> from sage.all import *
>>> MS = MatrixSpace(GF(Integer(2)),Integer(4),Integer(7))
>>> G  = MS([[Integer(1),Integer(1),Integer(1),Integer(0),Integer(0),Integer(0),Integer(0)], [Integer(1),Integer(0),Integer(0),Integer(1),Integer(1),Integer(0),Integer(0)], [Integer(0),Integer(1),Integer(0),Integer(1),Integer(0),Integer(1),Integer(0)], [Integer(1),Integer(1),Integer(0),Integer(1),Integer(0),Integer(0),Integer(1)]])
>>> C  = LinearCode(G)
>>> C
[7, 4] linear code over GF(2)
>>> C.base_ring()
Finite Field of size 2
>>> C.dimension()
4
>>> C.length()
7
>>> C.minimum_distance()                                                      # needs sage.libs.gap
3
>>> C.spectrum()
[1, 0, 0, 7, 7, 0, 0, 1]
>>> C.weight_distribution()
[1, 0, 0, 7, 7, 0, 0, 1]

The minimum distance of the code, if known, can be provided as an optional parameter.:

sage: C  = LinearCode(G, d=3)
sage: C.minimum_distance()                                                      # needs sage.libs.gap
3
>>> from sage.all import *
>>> C  = LinearCode(G, d=Integer(3))
>>> C.minimum_distance()                                                      # needs sage.libs.gap
3

Another example.:

sage: MS = MatrixSpace(GF(5),4,7)
sage: G  = MS([[1,1,1,0,0,0,0], [1,0,0,1,1,0,0], [0,1,0,1,0,1,0], [1,1,0,1,0,0,1]])
sage: C  = LinearCode(G); C
[7, 4] linear code over GF(5)
>>> from sage.all import *
>>> MS = MatrixSpace(GF(Integer(5)),Integer(4),Integer(7))
>>> G  = MS([[Integer(1),Integer(1),Integer(1),Integer(0),Integer(0),Integer(0),Integer(0)], [Integer(1),Integer(0),Integer(0),Integer(1),Integer(1),Integer(0),Integer(0)], [Integer(0),Integer(1),Integer(0),Integer(1),Integer(0),Integer(1),Integer(0)], [Integer(1),Integer(1),Integer(0),Integer(1),Integer(0),Integer(0),Integer(1)]])
>>> C  = LinearCode(G); C
[7, 4] linear code over GF(5)

Providing a code as the parameter in order to “forget” its structure (see Issue #20198):

sage: C = codes.GeneralizedReedSolomonCode(GF(23).list(), 12)
sage: LinearCode(C)
[23, 12] linear code over GF(23)
>>> from sage.all import *
>>> C = codes.GeneralizedReedSolomonCode(GF(Integer(23)).list(), Integer(12))
>>> LinearCode(C)
[23, 12] linear code over GF(23)

Another example:

sage: C = codes.HammingCode(GF(7), 3)
sage: C
[57, 54] Hamming Code over GF(7)
sage: LinearCode(C)
[57, 54] linear code over GF(7)
>>> from sage.all import *
>>> C = codes.HammingCode(GF(Integer(7)), Integer(3))
>>> C
[57, 54] Hamming Code over GF(7)
>>> LinearCode(C)
[57, 54] linear code over GF(7)

AUTHORS:

  • David Joyner (11-2005)

  • Charles Prior (03-2016): Issue #20198, LinearCode from a code

generator_matrix(encoder_name=None, **kwargs)[source]#

Return a generator matrix of self.

INPUT:

  • encoder_name – (default: None) name of the encoder which will be used to compute the generator matrix. self._generator_matrix will be returned if default value is kept.

  • kwargs – all additional arguments are forwarded to the construction of the encoder that is used.

EXAMPLES:

sage: G = matrix(GF(3),2,[1,-1,1,-1,1,1])
sage: code = LinearCode(G)
sage: code.generator_matrix()
[1 2 1]
[2 1 1]
>>> from sage.all import *
>>> G = matrix(GF(Integer(3)),Integer(2),[Integer(1),-Integer(1),Integer(1),-Integer(1),Integer(1),Integer(1)])
>>> code = LinearCode(G)
>>> code.generator_matrix()
[1 2 1]
[2 1 1]
class sage.coding.linear_code.LinearCodeGeneratorMatrixEncoder(code)[source]#

Bases: Encoder

Encoder based on generator_matrix for Linear codes.

This is the default encoder of a generic linear code, and should never be used for other codes than LinearCode.

INPUT:

  • code – The associated LinearCode of this encoder.

generator_matrix()[source]#

Return a generator matrix of the associated code of self.

EXAMPLES:

sage: G = Matrix(GF(2), [[1,1,1,0,0,0,0], [1,0,0,1,1,0,0],
....:                    [0,1,0,1,0,1,0], [1,1,0,1,0,0,1]])
sage: C = LinearCode(G)
sage: E = codes.encoders.LinearCodeGeneratorMatrixEncoder(C)
sage: E.generator_matrix()
[1 1 1 0 0 0 0]
[1 0 0 1 1 0 0]
[0 1 0 1 0 1 0]
[1 1 0 1 0 0 1]
>>> from sage.all import *
>>> G = Matrix(GF(Integer(2)), [[Integer(1),Integer(1),Integer(1),Integer(0),Integer(0),Integer(0),Integer(0)], [Integer(1),Integer(0),Integer(0),Integer(1),Integer(1),Integer(0),Integer(0)],
...                    [Integer(0),Integer(1),Integer(0),Integer(1),Integer(0),Integer(1),Integer(0)], [Integer(1),Integer(1),Integer(0),Integer(1),Integer(0),Integer(0),Integer(1)]])
>>> C = LinearCode(G)
>>> E = codes.encoders.LinearCodeGeneratorMatrixEncoder(C)
>>> E.generator_matrix()
[1 1 1 0 0 0 0]
[1 0 0 1 1 0 0]
[0 1 0 1 0 1 0]
[1 1 0 1 0 0 1]
class sage.coding.linear_code.LinearCodeNearestNeighborDecoder(code)[source]#

Bases: Decoder

Construct a decoder for Linear Codes. This decoder will decode to the nearest codeword found.

INPUT:

  • code – A code associated to this decoder

decode_to_code(r)[source]#

Corrects the errors in word and returns a codeword.

INPUT:

  • r – a codeword of self

OUTPUT:

  • a vector of self’s message space

EXAMPLES:

sage: G = Matrix(GF(2), [[1,1,1,0,0,0,0], [1,0,0,1,1,0,0],
....:                    [0,1,0,1,0,1,0], [1,1,0,1,0,0,1]])
sage: C = LinearCode(G)
sage: D = codes.decoders.LinearCodeNearestNeighborDecoder(C)
sage: word = vector(GF(2), (1, 1, 0, 0, 1, 1, 0))
sage: w_err = word + vector(GF(2), (1, 0, 0, 0, 0, 0, 0))
sage: D.decode_to_code(w_err)
(1, 1, 0, 0, 1, 1, 0)
>>> from sage.all import *
>>> G = Matrix(GF(Integer(2)), [[Integer(1),Integer(1),Integer(1),Integer(0),Integer(0),Integer(0),Integer(0)], [Integer(1),Integer(0),Integer(0),Integer(1),Integer(1),Integer(0),Integer(0)],
...                    [Integer(0),Integer(1),Integer(0),Integer(1),Integer(0),Integer(1),Integer(0)], [Integer(1),Integer(1),Integer(0),Integer(1),Integer(0),Integer(0),Integer(1)]])
>>> C = LinearCode(G)
>>> D = codes.decoders.LinearCodeNearestNeighborDecoder(C)
>>> word = vector(GF(Integer(2)), (Integer(1), Integer(1), Integer(0), Integer(0), Integer(1), Integer(1), Integer(0)))
>>> w_err = word + vector(GF(Integer(2)), (Integer(1), Integer(0), Integer(0), Integer(0), Integer(0), Integer(0), Integer(0)))
>>> D.decode_to_code(w_err)
(1, 1, 0, 0, 1, 1, 0)
decoding_radius()[source]#

Return maximal number of errors self can decode.

EXAMPLES:

sage: G = Matrix(GF(2), [[1,1,1,0,0,0,0], [1,0,0,1,1,0,0],
....:                    [0,1,0,1,0,1,0], [1,1,0,1,0,0,1]])
sage: C = LinearCode(G)
sage: D = codes.decoders.LinearCodeNearestNeighborDecoder(C)
sage: D.decoding_radius()                                                   # needs sage.libs.gap
1
>>> from sage.all import *
>>> G = Matrix(GF(Integer(2)), [[Integer(1),Integer(1),Integer(1),Integer(0),Integer(0),Integer(0),Integer(0)], [Integer(1),Integer(0),Integer(0),Integer(1),Integer(1),Integer(0),Integer(0)],
...                    [Integer(0),Integer(1),Integer(0),Integer(1),Integer(0),Integer(1),Integer(0)], [Integer(1),Integer(1),Integer(0),Integer(1),Integer(0),Integer(0),Integer(1)]])
>>> C = LinearCode(G)
>>> D = codes.decoders.LinearCodeNearestNeighborDecoder(C)
>>> D.decoding_radius()                                                   # needs sage.libs.gap
1
class sage.coding.linear_code.LinearCodeSyndromeDecoder(code, maximum_error_weight=None)[source]#

Bases: Decoder

Constructs a decoder for Linear Codes based on syndrome lookup table.

The decoding algorithm works as follows:

  • First, a lookup table is built by computing the syndrome of every error pattern of weight up to maximum_error_weight.

  • Then, whenever one tries to decode a word r, the syndrome of r is computed. The corresponding error pattern is recovered from the pre-computed lookup table.

  • Finally, the recovered error pattern is subtracted from r to recover the original word.

maximum_error_weight need never exceed the covering radius of the code, since there are then always lower-weight errors with the same syndrome. If one sets maximum_error_weight to a value greater than the covering radius, then the covering radius will be determined while building the lookup-table. This lower value is then returned if you query decoding_radius after construction.

If maximum_error_weight is left unspecified or set to a number at least the covering radius of the code, this decoder is complete, i.e. it decodes every vector in the ambient space.

Note

Constructing the lookup table takes time exponential in the length of the code and the size of the code’s base field. Afterwards, the individual decodings are fast.

INPUT:

  • code – A code associated to this decoder

  • maximum_error_weight – (default: None) the maximum number of errors to look for when building the table. An error is raised if it is set greater than \(n-k\), since this is an upper bound on the covering radius on any linear code. If maximum_error_weight is kept unspecified, it will be set to \(n - k\), where \(n\) is the length of code and \(k\) its dimension.

EXAMPLES:

sage: G = Matrix(GF(3), [[1,0,0,1,0,1,0,1,2],
....:                    [0,1,0,2,2,0,1,1,0],
....:                    [0,0,1,0,2,2,2,1,2]])
sage: C = LinearCode(G)
sage: D = codes.decoders.LinearCodeSyndromeDecoder(C)
sage: D
Syndrome decoder for [9, 3] linear code over GF(3) handling errors of weight up to 4
>>> from sage.all import *
>>> G = Matrix(GF(Integer(3)), [[Integer(1),Integer(0),Integer(0),Integer(1),Integer(0),Integer(1),Integer(0),Integer(1),Integer(2)],
...                    [Integer(0),Integer(1),Integer(0),Integer(2),Integer(2),Integer(0),Integer(1),Integer(1),Integer(0)],
...                    [Integer(0),Integer(0),Integer(1),Integer(0),Integer(2),Integer(2),Integer(2),Integer(1),Integer(2)]])
>>> C = LinearCode(G)
>>> D = codes.decoders.LinearCodeSyndromeDecoder(C)
>>> D
Syndrome decoder for [9, 3] linear code over GF(3) handling errors of weight up to 4

If one wants to correct up to a lower number of errors, one can do as follows:

sage: D = codes.decoders.LinearCodeSyndromeDecoder(C, maximum_error_weight=2)
sage: D
Syndrome decoder for [9, 3] linear code over GF(3) handling errors of weight up to 2
>>> from sage.all import *
>>> D = codes.decoders.LinearCodeSyndromeDecoder(C, maximum_error_weight=Integer(2))
>>> D
Syndrome decoder for [9, 3] linear code over GF(3) handling errors of weight up to 2

If one checks the list of types of this decoder before constructing it, one will notice it contains the keyword dynamic. Indeed, the behaviour of the syndrome decoder depends on the maximum error weight one wants to handle, and how it compares to the minimum distance and the covering radius of code. In the following examples, we illustrate this property by computing different instances of syndrome decoder for the same code.

We choose the following linear code, whose covering radius equals to 4 and minimum distance to 5 (half the minimum distance is 2):

sage: G = matrix(GF(5), [[1, 0, 0, 0, 0, 4, 3, 0, 3, 1, 0],
....:                    [0, 1, 0, 0, 0, 3, 2, 2, 3, 2, 1],
....:                    [0, 0, 1, 0, 0, 1, 3, 0, 1, 4, 1],
....:                    [0, 0, 0, 1, 0, 3, 4, 2, 2, 3, 3],
....:                    [0, 0, 0, 0, 1, 4, 2, 3, 2, 2, 1]])
sage: C = LinearCode(G)
>>> from sage.all import *
>>> G = matrix(GF(Integer(5)), [[Integer(1), Integer(0), Integer(0), Integer(0), Integer(0), Integer(4), Integer(3), Integer(0), Integer(3), Integer(1), Integer(0)],
...                    [Integer(0), Integer(1), Integer(0), Integer(0), Integer(0), Integer(3), Integer(2), Integer(2), Integer(3), Integer(2), Integer(1)],
...                    [Integer(0), Integer(0), Integer(1), Integer(0), Integer(0), Integer(1), Integer(3), Integer(0), Integer(1), Integer(4), Integer(1)],
...                    [Integer(0), Integer(0), Integer(0), Integer(1), Integer(0), Integer(3), Integer(4), Integer(2), Integer(2), Integer(3), Integer(3)],
...                    [Integer(0), Integer(0), Integer(0), Integer(0), Integer(1), Integer(4), Integer(2), Integer(3), Integer(2), Integer(2), Integer(1)]])
>>> C = LinearCode(G)

In the following examples, we illustrate how the choice of maximum_error_weight influences the types of the instance of syndrome decoder, alongside with its decoding radius.

We build a first syndrome decoder, and pick a maximum_error_weight smaller than both the covering radius and half the minimum distance:

sage: D = C.decoder("Syndrome", maximum_error_weight=1)
sage: D.decoder_type()
{'always-succeed', 'bounded_distance', 'hard-decision'}
sage: D.decoding_radius()
1
>>> from sage.all import *
>>> D = C.decoder("Syndrome", maximum_error_weight=Integer(1))
>>> D.decoder_type()
{'always-succeed', 'bounded_distance', 'hard-decision'}
>>> D.decoding_radius()
1

In that case, we are sure the decoder will always succeed. It is also a bounded distance decoder.

We now build another syndrome decoder, and this time, maximum_error_weight is chosen to be bigger than half the minimum distance, but lower than the covering radius:

sage: D = C.decoder("Syndrome", maximum_error_weight=3)
sage: D.decoder_type()
{'bounded_distance', 'hard-decision', 'might-error'}
sage: D.decoding_radius()
3
>>> from sage.all import *
>>> D = C.decoder("Syndrome", maximum_error_weight=Integer(3))
>>> D.decoder_type()
{'bounded_distance', 'hard-decision', 'might-error'}
>>> D.decoding_radius()
3

Here, we still get a bounded distance decoder. But because we have a maximum error weight bigger than half the minimum distance, we know it might return a codeword which was not the original codeword.

And now, we build a third syndrome decoder, whose maximum_error_weight is bigger than both the covering radius and half the minimum distance:

sage: D = C.decoder("Syndrome", maximum_error_weight=5)         # long time
sage: D.decoder_type()                                          # long time
{'complete', 'hard-decision', 'might-error'}
sage: D.decoding_radius()                                       # long time
4
>>> from sage.all import *
>>> D = C.decoder("Syndrome", maximum_error_weight=Integer(5))         # long time
>>> D.decoder_type()                                          # long time
{'complete', 'hard-decision', 'might-error'}
>>> D.decoding_radius()                                       # long time
4

In that case, the decoder might still return an unexpected codeword, but it is now complete. Note the decoding radius is equal to 4: it was determined while building the syndrome lookup table that any error with weight more than 4 will be decoded incorrectly. That is because the covering radius for the code is 4.

The minimum distance and the covering radius are both determined while computing the syndrome lookup table. They user did not explicitly ask to compute these on the code C. The dynamic typing of the syndrome decoder might therefore seem slightly surprising, but in the end is quite informative.

decode_to_code(r)[source]#

Correct the errors in word and return a codeword.

INPUT:

  • r – a codeword of self

OUTPUT:

  • a vector of self’s message space

EXAMPLES:

sage: G = Matrix(GF(3), [[1, 0, 0, 0, 2, 2, 1, 1],
....:                    [0, 1, 0, 0, 0, 0, 1, 1],
....:                    [0, 0, 1, 0, 2, 0, 0, 2],
....:                    [0, 0, 0, 1, 0, 2, 0, 1]])
sage: C = LinearCode(G)
sage: D = codes.decoders.LinearCodeSyndromeDecoder(C, maximum_error_weight = 1)
sage: Chan = channels.StaticErrorRateChannel(C.ambient_space(), 1)
sage: c = C.random_element()
sage: r = Chan(c)
sage: c == D.decode_to_code(r)
True
>>> from sage.all import *
>>> G = Matrix(GF(Integer(3)), [[Integer(1), Integer(0), Integer(0), Integer(0), Integer(2), Integer(2), Integer(1), Integer(1)],
...                    [Integer(0), Integer(1), Integer(0), Integer(0), Integer(0), Integer(0), Integer(1), Integer(1)],
...                    [Integer(0), Integer(0), Integer(1), Integer(0), Integer(2), Integer(0), Integer(0), Integer(2)],
...                    [Integer(0), Integer(0), Integer(0), Integer(1), Integer(0), Integer(2), Integer(0), Integer(1)]])
>>> C = LinearCode(G)
>>> D = codes.decoders.LinearCodeSyndromeDecoder(C, maximum_error_weight = Integer(1))
>>> Chan = channels.StaticErrorRateChannel(C.ambient_space(), Integer(1))
>>> c = C.random_element()
>>> r = Chan(c)
>>> c == D.decode_to_code(r)
True
decoding_radius()[source]#

Return the maximal number of errors a received word can have and for which self is guaranteed to return a most likely codeword.

EXAMPLES:

sage: G = Matrix(GF(3), [[1,0,0,1,0,1,0,1,2],
....:                    [0,1,0,2,2,0,1,1,0],
....:                    [0,0,1,0,2,2,2,1,2]])
sage: C = LinearCode(G)
sage: D = codes.decoders.LinearCodeSyndromeDecoder(C)
sage: D.decoding_radius()
4
>>> from sage.all import *
>>> G = Matrix(GF(Integer(3)), [[Integer(1),Integer(0),Integer(0),Integer(1),Integer(0),Integer(1),Integer(0),Integer(1),Integer(2)],
...                    [Integer(0),Integer(1),Integer(0),Integer(2),Integer(2),Integer(0),Integer(1),Integer(1),Integer(0)],
...                    [Integer(0),Integer(0),Integer(1),Integer(0),Integer(2),Integer(2),Integer(2),Integer(1),Integer(2)]])
>>> C = LinearCode(G)
>>> D = codes.decoders.LinearCodeSyndromeDecoder(C)
>>> D.decoding_radius()
4
maximum_error_weight()[source]#

Return the maximal number of errors a received word can have and for which self is guaranteed to return a most likely codeword.

Same as decoding_radius().

EXAMPLES:

sage: G = Matrix(GF(3), [[1,0,0,1,0,1,0,1,2],
....:                    [0,1,0,2,2,0,1,1,0],
....:                    [0,0,1,0,2,2,2,1,2]])
sage: C = LinearCode(G)
sage: D = codes.decoders.LinearCodeSyndromeDecoder(C)
sage: D.maximum_error_weight()
4
>>> from sage.all import *
>>> G = Matrix(GF(Integer(3)), [[Integer(1),Integer(0),Integer(0),Integer(1),Integer(0),Integer(1),Integer(0),Integer(1),Integer(2)],
...                    [Integer(0),Integer(1),Integer(0),Integer(2),Integer(2),Integer(0),Integer(1),Integer(1),Integer(0)],
...                    [Integer(0),Integer(0),Integer(1),Integer(0),Integer(2),Integer(2),Integer(2),Integer(1),Integer(2)]])
>>> C = LinearCode(G)
>>> D = codes.decoders.LinearCodeSyndromeDecoder(C)
>>> D.maximum_error_weight()
4
syndrome_table()[source]#

Return the syndrome lookup table of self.

EXAMPLES:

sage: G = Matrix(GF(2), [[1,1,1,0,0,0,0], [1,0,0,1,1,0,0],
....:                    [0,1,0,1,0,1,0], [1,1,0,1,0,0,1]])
sage: C = LinearCode(G)
sage: D = codes.decoders.LinearCodeSyndromeDecoder(C)
sage: D.syndrome_table()
{(0, 0, 0): (0, 0, 0, 0, 0, 0, 0),
 (1, 0, 0): (1, 0, 0, 0, 0, 0, 0),
 (0, 1, 0): (0, 1, 0, 0, 0, 0, 0),
 (1, 1, 0): (0, 0, 1, 0, 0, 0, 0),
 (0, 0, 1): (0, 0, 0, 1, 0, 0, 0),
 (1, 0, 1): (0, 0, 0, 0, 1, 0, 0),
 (0, 1, 1): (0, 0, 0, 0, 0, 1, 0),
 (1, 1, 1): (0, 0, 0, 0, 0, 0, 1)}
>>> from sage.all import *
>>> G = Matrix(GF(Integer(2)), [[Integer(1),Integer(1),Integer(1),Integer(0),Integer(0),Integer(0),Integer(0)], [Integer(1),Integer(0),Integer(0),Integer(1),Integer(1),Integer(0),Integer(0)],
...                    [Integer(0),Integer(1),Integer(0),Integer(1),Integer(0),Integer(1),Integer(0)], [Integer(1),Integer(1),Integer(0),Integer(1),Integer(0),Integer(0),Integer(1)]])
>>> C = LinearCode(G)
>>> D = codes.decoders.LinearCodeSyndromeDecoder(C)
>>> D.syndrome_table()
{(0, 0, 0): (0, 0, 0, 0, 0, 0, 0),
 (1, 0, 0): (1, 0, 0, 0, 0, 0, 0),
 (0, 1, 0): (0, 1, 0, 0, 0, 0, 0),
 (1, 1, 0): (0, 0, 1, 0, 0, 0, 0),
 (0, 0, 1): (0, 0, 0, 1, 0, 0, 0),
 (1, 0, 1): (0, 0, 0, 0, 1, 0, 0),
 (0, 1, 1): (0, 0, 0, 0, 0, 1, 0),
 (1, 1, 1): (0, 0, 0, 0, 0, 0, 1)}