Generic structures for linear codes

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
Linear code of length 5, dimension 3 over Finite Field of size 2
sage: C.generator_matrix()
[1 0 0 1 0]
[0 1 0 1 1]
[0 0 1 1 1]

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

Further references

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

If you want to learn more on the design of this library, see http://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): trac ticket ##3893: random_element() should pass on *args and **kwds.
  • Thomas Feulner (2012-11): trac ticket #13723: deprecation of hamming_weight()
  • Thomas Feulner (2013-10): added methods to compute a canonical representative and the automorphism group

TESTS:

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 == loads(dumps(C))
True
class sage.coding.linear_code.AbstractLinearCode(base_field, length, default_encoder_name, default_decoder_name)

Bases: sage.modules.module.Module

Abstract 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(SubclassName, self).__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.
  • 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 AbstractLinearCode is not designed to be implemented, it does not have any representation methods. You should implement _repr_ and _latex_ methods in the subclass.

Note

AbstractLinearCode has generic implementations of the comparison methods __cmp__ and __eq__ which use the generator matrix and are quite slow. In subclasses you are encouraged to override these functions.

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.

add_decoder(name, decoder)

Adds an decoder to the list of registered decoders of self.

Note

This method only adds decoder to self, and not to any member of the class of self. To know how to add an sage.coding.decoder.Decoder, please refer to the documentation of AbstractLinearCode.

INPUT:

  • name – the string name for the decoder
  • decoder – the class name of the decoder

EXAMPLES:

First of all, we create a (very basic) new decoder:

sage: class MyDecoder(sage.coding.decoder.Decoder):
....:   def __init__(self, code):
....:       super(MyDecoder, self).__init__(code)
....:   def _repr_(self):
....:       return "MyDecoder decoder with associated code %s" % self.code()

We now create a new code:

sage: C = codes.HammingCode(GF(2), 3)

We can add our new decoder to the list of available decoders of C:

sage: C.add_decoder("MyDecoder", MyDecoder)
sage: C.decoders_available()
['MyDecoder', 'Syndrome', 'NearestNeighbor']

We can verify that any new code will not know MyDecoder:

sage: C2 = codes.HammingCode(GF(2), 3)
sage: C2.decoders_available()
['Syndrome', 'NearestNeighbor']

TESTS:

It is impossible to use a name which is in the dictionary of available decoders:

sage: C.add_decoder("Syndrome", MyDecoder)
Traceback (most recent call last):
...
ValueError: There is already a registered decoder with this name
add_encoder(name, encoder)

Adds an encoder to the list of registered encoders of self.

Note

This method only adds encoder to self, and not to any member of the class of self. To know how to add an sage.coding.encoder.Encoder, please refer to the documentation of AbstractLinearCode.

INPUT:

  • name – the string name for the encoder
  • encoder – the class name of the encoder

EXAMPLES:

First of all, we create a (very basic) new encoder:

sage: class MyEncoder(sage.coding.encoder.Encoder):
....:   def __init__(self, code):
....:       super(MyEncoder, self).__init__(code)
....:   def _repr_(self):
....:       return "MyEncoder encoder with associated code %s" % self.code()

We now create a new code:

sage: C = codes.HammingCode(GF(2), 3)

We can add our new encoder to the list of available encoders of C:

sage: C.add_encoder("MyEncoder", MyEncoder)
sage: C.encoders_available()
['MyEncoder', 'ParityCheck']

We can verify that any new code will not know MyEncoder:

sage: C2 = codes.HammingCode(GF(2), 3)
sage: C2.encoders_available()
['ParityCheck']

TESTS:

It is impossible to use a name which is in the dictionary of available encoders:

sage: C.add_encoder("ParityCheck", MyEncoder)
Traceback (most recent call last):
...
ValueError: There is already a registered encoder with this name
ambient_space()

Returns the ambient vector space of self.

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: C.ambient_space()
Vector space of dimension 7 over Finite Field of size 2
assmus_mattson_designs(t, mode=None)

Assmus and Mattson Theorem (section 8.4, page 303 of [HP]): 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 [HP])

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.ExtendedBinaryGolayCode()             #  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() for c in C if c.hamming_weight()==8]; len(blocks)  # long time computation
759

REFERENCE:

  • [HP] W. C. Huffman and V. Pless, Fundamentals of ECC, Cambridge Univ. Press, 2003.
automorphism_group_gens(equivalence='semilinear')

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:

sage: C = codes.HammingCode(GF(4, 'z'), 3)
sage: C.automorphism_group_gens()
([((1, 1, 1, z, z + 1, z + 1, z + 1, z, z, 1, 1, 1, z, z, z + 1, z, z, z + 1, z + 1, z + 1, 1); (1,6,12,17)(2,16,4,5,11,8,14,13)(3,21,19,10,20,18,15,9), Ring endomorphism of Finite Field in z of size 2^2
      Defn: z |--> z + 1), ((1, 1, 1, z, z + 1, 1, 1, z, z, z + 1, z, z, z + 1, z + 1, z + 1, 1, z + 1, z, z, 1, 1); (1,6,9,13,15,18)(2,21)(3,16,7)(4,5,11,10,12,14)(17,19), Ring endomorphism of Finite Field in z of size 2^2
      Defn: z |--> z), ((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, z, 1, 1, z + 1, z, z + 1, z, z, z + 1, 1, 1, 1, z + 1, z, z, z + 1, z + 1, 1, 1, z); (1,5,10,9,4,14,11,16,18,20,6,19,12,15,3,8,2,17,7,13,21), Ring endomorphism of Finite Field in z of size 2^2
      Defn: z |--> z), ((z + 1, 1, z, 1, 1, z + 1, z + 1, z, 1, z, z + 1, z, z + 1, z + 1, z, 1, 1, z + 1, z + 1, z + 1, z); (1,17,10)(2,15,13)(4,11,21)(5,18,12)(6,14,19)(7,8,16), 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)
base_field()

Return the base field 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: C.base_field()
Finite Field of size 2
basis()

Returns a basis of self.

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: C.basis()
[(1, 0, 0, 0, 0, 1, 1), (0, 1, 0, 0, 1, 0, 1), (0, 0, 1, 0, 1, 1, 0), (0, 0, 0, 1, 1, 1, 1)]
binomial_moment(i)

Returns 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)
0
sage: C.binomial_moment(4)    # long time
35

Warning

This is slow.

REFERENCE:

  • I. Duursma, “Combinatorics of the two-variable zeta function”, Finite fields and applications, 109-136, Lecture Notes in Comput. Sci., 2948, Springer, Berlin, 2004.
canonical_representative(equivalence='semilinear')

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()

Check that the transporter element is correct:

sage: LinearCode(transp*C.generator_matrix()) == CanRep
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()
sage: CanRep_iso == CanRep
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")
sage: CanRep2, _ = C_iso.canonical_representative("linear")
sage: CanRep2 == CanRep1
True
cardinality()

Return the size of this code.

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: C.cardinality()
16
sage: len(C)
16
characteristic()

Returns the characteristic of the base ring of self.

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: C.characteristic()
2
characteristic_polynomial()

Returns the characteristic polynomial of a linear code, as defined in van Lint’s text [vL].

EXAMPLES:

sage: C = codes.ExtendedBinaryGolayCode()
sage: C.characteristic_polynomial()
-4/3*x^3 + 64*x^2 - 2816/3*x + 4096

REFERENCES:

  • van Lint, Introduction to coding theory, 3rd ed., Springer-Verlag GTM, 86, 1999.
chinen_polynomial()

Returns 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.TernaryGolayCode()
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)

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.
covering_radius()

Wraps Guava’s CoveringRadius command.

The covering radius of a linear code C is the smallest number r with the property that each element v of the ambient vector space of C has at most a distance r to the code C. So for each vector v there must be an element c of C with d(v,c) \leq  r. A binary linear code with reasonable small covering radius is often referred to as a covering code.

For example, if C is a perfect code, the covering radius is equal to t, the number of errors the code can correct, where d = 2t+1, with d the minimum distance of C.

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 5)
sage: C.covering_radius()  # optional - gap_packages (Guava package)
1
decode(right, algorithm='syndrome')

Corrects the errors in right and returns a codeword.

INPUT:

  • right – a vector of the same length as self over the base field of self
  • algorithm – (default: 'syndrome') Name of the decoding algorithm which will be used to decode right. Can be 'syndrome' or 'nearest_neighbor'.

Note

This is a deprecated method which will soon be removed from Sage. Please use decode_to_code() instead.

decode_to_code(word, decoder_name=None, **kwargs)

Corrects the errors in word and returns a codeword.

INPUT:

  • word – a vector of the same length as self over the base field of self
  • decoder_name – (default: None) Name of the decoder which will be used to decode word. The default decoder of self will be used if default value is kept.
  • kwargs – all additional arguments are forwarded to decoder()

OUTPUT:

  • A vector 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: 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: C.decode_to_code(w_err)
(1, 1, 0, 0, 1, 1, 0)

It is possible to manually choose the decoder amongst the list of the available ones:

sage: C.decoders_available()
['Syndrome', 'NearestNeighbor']
sage: C.decode_to_code(w_err, 'NearestNeighbor')
(1, 1, 0, 0, 1, 1, 0)
decode_to_message(word, decoder_name=None, **kwargs)

Correct the errors in word and decodes it to the message space.

INPUT:

  • word – a vector of the same length as self over the base field of self
  • decoder_name – (default: None) Name of the decoder which will be used to decode word. The default decoder of self will be used if default value is kept.
  • kwargs – all additional arguments are forwarded to decoder()

OUTPUT:

  • A vector of the message space 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: word = vector(GF(2), (1, 1, 0, 0, 1, 1, 0))
sage: C.decode_to_message(word)
(0, 1, 1, 0)

It is possible to manually choose the decoder amongst the list of the available ones:

sage: C.decoders_available()
['Syndrome', 'NearestNeighbor']
sage: C.decode_to_message(word, 'NearestNeighbor')
(0, 1, 1, 0)
decoder(decoder_name=None, **kwargs)

Return a decoder of self.

INPUT:

  • decoder_name – (default: None) name of the decoder which will be returned. The default decoder of self will be used if default value is kept.
  • kwargs – all additional arguments will be forwarded to the constructor of the decoder that will be returned by this method

OUTPUT:

  • a decoder object

Besides creating the decoder and returning it, this method also stores the decoder in a cache. With this behaviour, each decoder will be created at most one time for 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: C.decoder()
Syndrome decoder for Linear code of length 7, dimension 4 over Finite Field of size 2 handling errors of weight up to 1

If the name of a decoder which is not known by self is passed, an exception will be raised:

sage: C.decoders_available()
['Syndrome', 'NearestNeighbor']
sage: C.decoder('Try')
Traceback (most recent call last):
...
ValueError: Passed Decoder name not known
decoders_available(classes=False)

Returns a list of the available decoders’ names for self.

INPUT:

  • classes – (default: False) if classes is set to True, it also returns the decoders’ classes associated with the decoders’ names.

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: C.decoders_available()
['Syndrome', 'NearestNeighbor']

sage: C.decoders_available(True)
{'NearestNeighbor': <class 'sage.coding.linear_code.LinearCodeNearestNeighborDecoder'>,
 'Syndrome': <class 'sage.coding.linear_code.LinearCodeSyndromeDecoder'>}
dimension()

Returns the dimension of this code.

EXAMPLES:

sage: G = matrix(GF(2),[[1,0,0],[1,1,0]])
sage: C = LinearCode(G)
sage: C.dimension()
2
direct_sum(other)

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
Linear code of length 14, dimension 8 over Finite Field of size 2
sage: C3 = C1.direct_sum(C2); C3
Linear code of length 21, dimension 12 over Finite Field of size 2
divisor()

Returns the divisor of a code, which is the smallest integer d_0 > 0 such that each A_i > 0 iff i is divisible by d_0.

EXAMPLES:

sage: C = codes.ExtendedBinaryGolayCode()
sage: C.divisor()   # Type II self-dual
4
sage: C = codes.QuadraticResidueCodeEvenPair(17,GF(2))[0]
sage: C.divisor()
2
dual_code()

Returns the dual code C^{\perp} of the code C,

C^{\perp} = \{ v \in V\ |\ v\cdot c = 0,\ \forall c \in C \}.

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: C.dual_code()
Linear code of length 7, dimension 3 over Finite Field of size 2
sage: C = codes.HammingCode(GF(4, 'a'), 3)
sage: C.dual_code()
Linear code of length 21, dimension 3 over Finite Field in a of size 2^2
encode(word, encoder_name=None, **kwargs)

Transforms an element of a message space into a codeword.

INPUT:

  • word – a vector of a message space of the code.
  • encoder_name – (default: None) Name of the encoder which will be used to encode word. The default encoder of self will be used if default value is kept.
  • kwargs – all additional arguments are forwarded to the construction of the encoder that is used.

Note

The default encoder always has F^{k} as message space, with k the dimension of self and F the base ring of self.

OUTPUT:

  • a vector 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: word = vector((0, 1, 1, 0))
sage: C.encode(word)
(1, 1, 0, 0, 1, 1, 0)

It is possible to manually choose the encoder amongst the list of the available ones:

sage: C.encoders_available()
['GeneratorMatrix']
sage: word = vector((0, 1, 1, 0))
sage: C.encode(word, 'GeneratorMatrix')
(1, 1, 0, 0, 1, 1, 0)
encoder(encoder_name=None, **kwargs)

Returns an encoder of self.

The returned encoder provided by this method is cached.

This methods creates a new instance of the encoder subclass designated by encoder_name. While it is also possible to do the same by directly calling the subclass’ constructor, it is strongly advised to use this method to take advantage of the caching mechanism.

INPUT:

  • encoder_name – (default: None) name of the encoder which will be returned. The default encoder of self will be used if default value is kept.
  • kwargs – all additional arguments are forwarded to the constructor of the encoder this method will return.

OUTPUT:

  • an Encoder object.

Note

The default encoder always has F^{k} as message space, with k the dimension of self and F the base ring 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: C.encoder()
Generator matrix-based encoder for Linear code of length 7, dimension 4 over Finite Field of size 2

We check that the returned encoder is cached:

sage: C.encoder.is_in_cache()
True

If the name of an encoder which is not known by self is passed, an exception will be raised:

sage: C.encoders_available()
['GeneratorMatrix']
sage: C.encoder('NonExistingEncoder')
Traceback (most recent call last):
...
ValueError: Passed Encoder name not known
encoders_available(classes=False)

Returns a list of the available encoders’ names for self.

INPUT:

  • classes – (default: False) if classes is set to True, it also returns the encoders’ classes associated with the encoders’ names.

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: C.encoders_available()
['GeneratorMatrix']

sage: C.encoders_available(True)
{'GeneratorMatrix': <class 'sage.coding.linear_code.LinearCodeGeneratorMatrixEncoder'>}
extended_code()

Returns 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 Finite Field in a of size 2^2
sage: Cx = C.extended_code()
sage: Cx
Extended code coming from [21, 18] Hamming Code over Finite Field in a of size 2^2
galois_closure(F0)

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 Finite Field in a of size 2^2
Linear code of length 21, dimension 20 over Finite Field in a of size 2^2
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
generator_matrix(encoder_name=None, **kwargs)

Returns a generator matrix of self.

INPUT:

  • encoder_name – (default: None) name of the encoder which will be used to compute the generator matrix. The default encoder of self will be used 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]
generator_matrix_systematic()

Return a systematic generator matrix of the code.

A generator matrix of a code is called systematic if it contains a set of columns forming an identity matrix.

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]
sage: code.generator_matrix_systematic()
[1 2 0]
[0 0 1]
gens()

Returns the generators of this code as a list of vectors.

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: C.gens()
 [(1, 0, 0, 0, 0, 1, 1), (0, 1, 0, 0, 1, 0, 1), (0, 0, 1, 0, 1, 1, 0), (0, 0, 0, 1, 1, 1, 1)]
genus()

Returns 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 Finite Field of size 2
sage: C1.genus()
1
sage: C2 = codes.HammingCode(GF(4,"a"), 2); C2
[5, 3] Hamming Code over Finite Field in a of size 2^2
sage: C2.genus()
0

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

information_set()

Return an information set of the code.

Return value of this method is cached.

A set of column positions of a generator matrix of a code is called an information set if the corresponding columns form a square matrix of full rank.

OUTPUT:

  • Information set of a systematic generator matrix of the code.

EXAMPLES:

sage: G = matrix(GF(3),2,[1,-1,0,-1,1,1])
sage: code = LinearCode(G)
sage: code.generator_matrix_systematic()
[1 2 0]
[0 0 1]
sage: code.information_set()
(0, 2)
is_galois_closed()

Checks if self is equal to its Galois closure.

EXAMPLES:

sage: C = codes.HammingCode(GF(4,"a"), 3)
sage: C.is_galois_closed()
False
is_permutation_automorphism(g)

Returns 1 if g is an element of S_n (n = length of self) and if g is an automorphism of self.

EXAMPLES:

sage: C = codes.HammingCode(GF(3), 3)
sage: g = SymmetricGroup(13).random_element()
sage: C.is_permutation_automorphism(g)
0
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: S8 = SymmetricGroup(8)
sage: g = S8("(2,3)")
sage: C.is_permutation_automorphism(g)
1
sage: g = S8("(1,2,3,4)")
sage: C.is_permutation_automorphism(g)
0
is_permutation_equivalent(other, algorithm=None)

Returns 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.CyclicCodeFromGeneratingPolynomial(7,g); C1
Linear code of length 7, dimension 4 over Finite Field of size 2
sage: C2 = codes.HammingCode(GF(2), 3); C2
[7, 4] Hamming Code over Finite Field of size 2
sage: C1.is_permutation_equivalent(C2)
True
sage: C1.is_permutation_equivalent(C2,algorithm="verbose")
(True, (3,4)(5,7,6))
sage: C1 = codes.RandomLinearCode(10,5,GF(2))
sage: C2 = codes.RandomLinearCode(10,5,GF(3))
sage: C1.is_permutation_equivalent(C2)
False
is_projective()

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 [BS11] or 9.8.1 from [BH12]).

EXAMPLE:

sage: C = codes.BinaryGolayCode()
sage: C.is_projective()
True
sage: C.dual_code().minimum_distance()
8

A non-projective code:

sage: C = codes.LinearCode(matrix(GF(2),[[1,0,1],[1,1,1]]))
sage: C.is_projective()
False

REFERENCE:

[BS11]E. Byrne and A. Sneyd, On the Parameters of Codes with Two Homogeneous Weights. WCC 2011-Workshop on coding and cryptography, pp. 81-90. 2011. https://hal.inria.fr/inria-00607341/document
is_self_dual()

Returns True if the code is self-dual (in the usual Hamming inner product) and False otherwise.

EXAMPLES:

sage: C = codes.ExtendedBinaryGolayCode()
sage: C.is_self_dual()
True
sage: C = codes.HammingCode(GF(2), 3)
sage: C.is_self_dual()
False
is_self_orthogonal()

Returns True if this code is self-orthogonal and False otherwise.

A code is self-orthogonal if it is a subcode of its dual.

EXAMPLES:

sage: C = codes.ExtendedBinaryGolayCode()
sage: C.is_self_orthogonal()
True
sage: C = codes.HammingCode(GF(2), 3)
sage: C.is_self_orthogonal()
False
sage: C = codes.QuasiQuadraticResidueCode(11)  # optional - gap_packages (Guava package)
sage: C.is_self_orthogonal()             # optional - gap_packages (Guava package)
True
is_subcode(other)

Returns True if self is a subcode of other.

EXAMPLES:

sage: C1 = codes.HammingCode(GF(2), 3)
sage: G1 = C1.generator_matrix()
sage: G2 = G1.matrix_from_rows([0,1,2])
sage: C2 = LinearCode(G2)
sage: C2.is_subcode(C1)
True
sage: C1.is_subcode(C2)
False
sage: C3 = C1.extended_code()
sage: C1.is_subcode(C3)
False
sage: C4 = C1.punctured([1])
sage: C4.is_subcode(C1)
False
sage: C5 = C1.shortened([1])
sage: C5.is_subcode(C1)
False
sage: C1 = codes.HammingCode(GF(9,"z"), 3)
sage: G1 = C1.generator_matrix()
sage: G2 = G1.matrix_from_rows([0,1,2])
sage: C2 = LinearCode(G2)
sage: C2.is_subcode(C1)
True
length()

Returns the length of this code.

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: C.length()
7
list()

Return a list of all elements of this linear code.

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: Clist = C.list()
sage: Clist[5]; Clist[5] in C
(1, 0, 1, 0, 1, 0, 1)
True
minimum_distance(algorithm=None)

Returns the minimum distance of this linear code.

By default, this uses a GAP kernel function (in C and not part of Guava) written by Steve Linton. If algorithm="guava" is set and q is 2 or 3 then this uses a very fast program written in C written by CJ Tjhal. (This is much faster, except in some small examples.)

Raises a ValueError in case there is no non-zero vector in this linear code.

The minimum distance of the code is stored once it has been computed or provided during the initialization of LinearCode. If algorithm is None and the stored value of minimum distance is found, then the stored value will be returned without recomputing the minimum distance again.

INPUT:

  • algorithm - Method to be used, None, "gap", or "guava" (default: None).

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()
3

Once the minimum distance has been computed, it’s value is stored. Hence the following command will return the value instantly, without further computations.:

sage: C.minimum_distance()
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")
3
sage: C.minimum_distance(algorithm="guava")  # optional - gap_packages (Guava package)
3

Another example.:

sage: C = codes.HammingCode(GF(4,"a"), 2); C
[5, 3] Hamming Code over Finite Field in a of size 2^2
sage: C.minimum_distance()
3

TESTS:

sage: C = codes.RandomLinearCode(5, 2, GF(4,"a"))
sage: C.minimum_distance(algorithm='something')
Traceback (most recent call last):
...
ValueError: The algorithm argument must be one of None, 'gap' or 'guava'; got 'something'
module_composition_factors(gp)

Prints the GAP record of the Meataxe composition factors module 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()

Now type “C.module_composition_factors(gp)” to get the record printed.

parity_check_matrix()

Returns the parity check matrix of self.

The parity check matrix of a linear code C corresponds to the generator matrix of the dual code of C.

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: Cperp = C.dual_code()
sage: C; Cperp
[7, 4] Hamming Code over Finite Field of size 2
Linear code of length 7, dimension 3 over Finite Field of size 2
sage: C.generator_matrix()
 [1 0 0 0 0 1 1]
 [0 1 0 0 1 0 1]
 [0 0 1 0 1 1 0]
 [0 0 0 1 1 1 1]
sage: C.parity_check_matrix()
 [1 0 1 0 1 0 1]
 [0 1 1 0 0 1 1]
 [0 0 0 1 1 1 1]
sage: Cperp.parity_check_matrix()
 [1 0 0 0 0 1 1]
 [0 1 0 0 1 0 1]
 [0 0 1 0 1 1 0]
 [0 0 0 1 1 1 1]
sage: Cperp.generator_matrix()
 [1 0 1 0 1 0 1]
 [0 1 1 0 0 1 1]
 [0 0 0 1 1 1 1]
permutation_automorphism_group(algorithm='partition')

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)
sage: C
Linear code of length 8, dimension 4 over Finite Field of size 2
sage: G = C.permutation_automorphism_group()
sage: G.order()
144
sage: GG = C.permutation_automorphism_group("codecan")
sage: 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: C = codes.ExtendedTernaryGolayCode()
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

Other examples:

sage: C = codes.ExtendedBinaryGolayCode()
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 Finite Field of size 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 Finite Field in z of size 2^2
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_packages (Guava package)
Permutation Group with generators [(1,3)(4,5), (1,4)(3,5)]
sage: C = codes.TernaryGolayCode()
sage: C.permutation_automorphism_group(algorithm="gap")  # optional - gap_packages (Guava package)
Permutation Group with generators [(3,4)(5,7)(6,9)(8,11), (3,5,8)(4,11,7)(6,9,10), (2,3)(4,6)(5,8)(7,10), (1,2)(4,11)(5,8)(9,10)]

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").

permuted_code(p)

Returns the permuted code, which is equivalent to self via the column permutation p.

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: G = C.permutation_automorphism_group(); G
Permutation Group with generators [(4,5)(6,7), (4,6)(5,7), (2,3)(6,7), (2,4)(3,5), (1,2)(5,6)]
sage: g = G("(2,3)(6,7)")
sage: Cg = C.permuted_code(g)
sage: Cg
Linear code of length 7, dimension 4 over Finite Field of size 2
sage: C.generator_matrix() == Cg.generator_matrix_systematic()
True
punctured(L)

Returns 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])
Punctured code coming from [7, 4] Hamming Code over Finite Field of size 2 punctured on position(s) [1, 2]
random_element(*args, **kwds)

Returns a random codeword; passes other positional and keyword arguments to random_element() method of vector space.

OUTPUT:

  • Random element of the vector space of this code

EXAMPLES:

sage: C = codes.HammingCode(GF(4,'a'), 3)
sage: C.random_element() # random test
(1, 0, 0, a + 1, 1, a, a, a + 1, a + 1, 1, 1, 0, a + 1, a, 0, a, a, 0, a, a, 1)

Passes extra positional or keyword arguments through:

sage: C.random_element(prob=.5, distribution='1/n') # random test
(1, 0, a, 0, 0, 0, 0, a + 1, 0, 0, 0, 0, 0, 0, 0, 0, a + 1, a + 1, 1, 0, 0)

TESTS:

Test that the codeword returned is immutable (see trac ticket #16469):

sage: c = C.random_element()
sage: c.is_immutable()
True

Test that codeword returned has the same parent as any non-random codeword (see trac ticket #19653):

sage: C = codes.RandomLinearCode(10, 4, GF(16, 'a'))
sage: c1 = C.random_element()
sage: c2 = C[1]
sage: c1.parent() == c2.parent()
True
redundancy_matrix(C)

If C is a linear [n,k,d] code then this function returns a k \times (n-k) matrix A such that G = (I,A) generates a code (in standard form) equivalent to C. If C is already in standard form and G = (I,A) is its generator matrix then this function simply returns that A.

OUTPUT:

  • Matrix, the redundancy matrix

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: C.generator_matrix()
 [1 0 0 0 0 1 1]
 [0 1 0 0 1 0 1]
 [0 0 1 0 1 1 0]
 [0 0 0 1 1 1 1]
sage: C.redundancy_matrix()
 [0 1 1]
 [1 0 1]
 [1 1 0]
 [1 1 1]
sage: C.standard_form()[0].generator_matrix()
 [1 0 0 0 0 1 1]
 [0 1 0 0 1 0 1]
 [0 0 1 0 1 1 0]
 [0 0 0 1 1 1 1]
sage: C = codes.HammingCode(GF(3), 2)
sage: C.generator_matrix()
[1 0 1 1]
[0 1 1 2]
sage: C.redundancy_matrix()
[1 1]
[1 2]
sd_duursma_data(C, i)

Returns the Duursma data v and m of this formally s.d. code C and the type number i in (1,2,3,4). Does not check if this code is actually sd.

INPUT:

  • i - Type number

OUTPUT:

  • Pair (v, m) as in Duursma [D]

REFERENCES:

[D](1, 2, 3) I. Duursma, “Extremal weight enumerators and ultraspherical polynomials”

EXAMPLES:

sage: MS = MatrixSpace(GF(2),2,4)
sage: G = MS([1,1,0,0,0,0,1,1])
sage: C = LinearCode(G)
sage: C == C.dual_code()  # checks that C is self dual
True
sage: for i in [1,2,3,4]: print(C.sd_duursma_data(i))
[2, -1]
[2, -3]
[2, -2]
[2, -1]
sd_duursma_q(C, i, d0)

INPUT:

  • C - sd code; does not check if C is actually an sd code
  • i - Type number, one of 1,2,3,4
  • d0 - Divisor, the smallest integer such that each A_i > 0 iff i is divisible by d0

OUTPUT:

  • Coefficients q_0, q_1, ... of q(T) as in Duursma [D]

REFERENCES:

  • [D] - I. Duursma, “Extremal weight enumerators and ultraspherical polynomials”

EXAMPLES:

sage: C1 = codes.HammingCode(GF(2), 3)
sage: C2 = C1.extended_code(); C2
Extended code coming from [7, 4] Hamming Code over Finite Field of size 2
sage: C2.sd_duursma_q(1,1)
2/5*T^2 + 2/5*T + 1/5
sage: C2.sd_duursma_q(3,1)
3/5*T^4 + 1/5*T^3 + 1/15*T^2 + 1/15*T + 1/15
sd_zeta_polynomial(C, typ=1)

Returns the Duursma zeta function of a self-dual code using the construction in [D].

INPUT:

  • typ - Integer, type of this s.d. code; one of 1,2,3, or 4 (default: 1)

OUTPUT:

  • Polynomial

EXAMPLES:

sage: C1 = codes.HammingCode(GF(2), 3)
sage: C2 = C1.extended_code(); C2
Extended code coming from [7, 4] Hamming Code over Finite Field of size 2
sage: C2.sd_zeta_polynomial()
2/5*T^2 + 2/5*T + 1/5
sage: C2.zeta_polynomial()
2/5*T^2 + 2/5*T + 1/5
sage: P = C2.sd_zeta_polynomial(); P(1)
1
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.sd_zeta_polynomial(4)
1

It is a general fact about Duursma zeta polynomials that P(1) = 1.

REFERENCES:

  • [D] I. Duursma, “Extremal weight enumerators and ultraspherical polynomials”
shortened(L)

Returns 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 [HP], 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])
Linear code of length 5, dimension 2 over Finite Field of size 2
spectrum(algorithm=None)

Returns the spectrum of self as a list.

The default algorithm uses a GAP kernel function (in C) written by Steve Linton.

INPUT:

  • algorithm - None, "gap", "leon", or "binary"; defaults to "gap" except in the binary case. If "gap" then uses the GAP function, if "leon" then uses Jeffrey Leon’s software via Guava, and if "binary" then uses Sage native Cython code
  • List, the spectrum

The optional algorithm ("leon") may create a stack smashing error and a traceback but should return the correct answer. It appears to run much faster than the GAP algorithm in some small examples and much slower than the GAP algorithm in other 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.spectrum()
[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 Finite Field in z of size 2^2
sage: C.spectrum()
[1, 0, 0, 30, 15, 18]
sage: C = codes.HammingCode(GF(2), 3); C
[7, 4] Hamming Code over Finite Field of size 2
sage: C.spectrum(algorithm="leon")   # optional - gap_packages (Guava package)
[1, 0, 0, 7, 7, 0, 0, 1]
sage: C.spectrum(algorithm="gap")
[1, 0, 0, 7, 7, 0, 0, 1]
sage: C.spectrum(algorithm="binary")
[1, 0, 0, 7, 7, 0, 0, 1]
sage: C = codes.HammingCode(GF(3), 3); C
[13, 10] Hamming Code over Finite Field of size 3
sage: C.spectrum() == C.spectrum(algorithm="leon")   # optional - gap_packages (Guava package)
True
sage: C = codes.HammingCode(GF(5), 2); C
[6, 4] Hamming Code over Finite Field of size 5
sage: C.spectrum() == C.spectrum(algorithm="leon")   # optional - gap_packages (Guava package)
True
sage: C = codes.HammingCode(GF(7), 2); C
[8, 6] Hamming Code over Finite Field of size 7
sage: C.spectrum() == C.spectrum(algorithm="leon")   # optional - gap_packages (Guava package)
True
standard_form()

Returns the standard form of this linear code.

An [n,k] linear code with generator matrix G in standard form is the row-reduced echelon form of G is (I,A), where I denotes the k \times k identity matrix and A is a k \times (n-k) block. This method returns a pair (C,p) where C is a code permutation equivalent to self and p in S_n, with n the length of C, is the permutation sending self to C. This does not call GAP.

Thanks to Frank Luebeck for (the GAP version of) this code.

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: C.generator_matrix()
[1 0 0 0 0 1 1]
[0 1 0 0 1 0 1]
[0 0 1 0 1 1 0]
[0 0 0 1 1 1 1]
sage: Cs,p = C.standard_form()
sage: p
()
sage: MS = MatrixSpace(GF(3),3,7)
sage: G = MS([[1,0,0,0,1,1,0],[0,1,0,1,0,1,0],[0,0,0,0,0,0,1]])
sage: C = LinearCode(G)
sage: Cs, p = C.standard_form()
sage: p
(3,7)
sage: Cs.generator_matrix()
 [1 0 0 0 1 1 0]
 [0 1 0 1 0 1 0]
 [0 0 1 0 0 0 0]
support()

Returns the set of indices j where A_j is nonzero, where spectrum(self) = [A_0,A_1,...,A_n].

OUTPUT:

  • List of integers

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: C.spectrum()
[1, 0, 0, 7, 7, 0, 0, 1]
sage: C.support()
[0, 3, 4, 7]
syndrome(r)

Returns the syndrome of r.

The syndrome of r is the result of H \times r where H is the parity check matrix of self. If r belongs to self, its syndrome equals to the zero vector.

INPUT:

  • r – a vector of the same length as self

OUTPUT:

  • a column vector

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: r = vector(GF(2), (1,0,1,0,1,0,1))
sage: r in C
True
sage: C.syndrome(r)
(0, 0, 0)

If r is not a codeword, its syndrome is not equal to zero:

sage: r = vector(GF(2), (1,0,1,0,1,1,1))
sage: r in C
False
sage: C.syndrome(r)
(0, 1, 1)

Syndrome computation works fine on bigger fields:

sage: C = codes.RandomLinearCode(12, 4, GF(59))
sage: r = C.random_element()
sage: C.syndrome(r)
(0, 0, 0, 0, 0, 0, 0, 0)
unencode(c, encoder_name=None, nocheck=False, **kwargs)

Returns the message corresponding to c.

This is the inverse of encode().

INPUT:

  • c – a codeword of self.
  • encoder_name – (default: None) name of the decoder which will be used to decode word. The default decoder of self will be used if default value is kept.
  • nocheck – (default: False) checks if c is in self. You might set this to True to disable the check for saving computation. Note that if c is not in self and nocheck = True, then the output of unencode() is not defined (except that it will be in the message space of self).
  • kwargs – all additional arguments are forwarded to the construction of the encoder that is used.

OUTPUT:

  • an element of the message space of encoder_name 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: c = vector(GF(2), (1, 1, 0, 0, 1, 1, 0))
sage: C.unencode(c)
(0, 1, 1, 0)
weight_distribution(algorithm=None)

Returns the spectrum of self as a list.

The default algorithm uses a GAP kernel function (in C) written by Steve Linton.

INPUT:

  • algorithm - None, "gap", "leon", or "binary"; defaults to "gap" except in the binary case. If "gap" then uses the GAP function, if "leon" then uses Jeffrey Leon’s software via Guava, and if "binary" then uses Sage native Cython code
  • List, the spectrum

The optional algorithm ("leon") may create a stack smashing error and a traceback but should return the correct answer. It appears to run much faster than the GAP algorithm in some small examples and much slower than the GAP algorithm in other 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.spectrum()
[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 Finite Field in z of size 2^2
sage: C.spectrum()
[1, 0, 0, 30, 15, 18]
sage: C = codes.HammingCode(GF(2), 3); C
[7, 4] Hamming Code over Finite Field of size 2
sage: C.spectrum(algorithm="leon")   # optional - gap_packages (Guava package)
[1, 0, 0, 7, 7, 0, 0, 1]
sage: C.spectrum(algorithm="gap")
[1, 0, 0, 7, 7, 0, 0, 1]
sage: C.spectrum(algorithm="binary")
[1, 0, 0, 7, 7, 0, 0, 1]
sage: C = codes.HammingCode(GF(3), 3); C
[13, 10] Hamming Code over Finite Field of size 3
sage: C.spectrum() == C.spectrum(algorithm="leon")   # optional - gap_packages (Guava package)
True
sage: C = codes.HammingCode(GF(5), 2); C
[6, 4] Hamming Code over Finite Field of size 5
sage: C.spectrum() == C.spectrum(algorithm="leon")   # optional - gap_packages (Guava package)
True
sage: C = codes.HammingCode(GF(7), 2); C
[8, 6] Hamming Code over Finite Field of size 7
sage: C.spectrum() == C.spectrum(algorithm="leon")   # optional - gap_packages (Guava package)
True
weight_enumerator(names='xy', name2=None)

Returns the weight enumerator of the code.

INPUT:

  • names - String of length 2, containing two variable names (default: "xy"). Alternatively, it can be a variable name or a string, or a tuple of variable names or strings.
  • name2 - string or symbolic variable (default: None). If name2 is provided then it is assumed that names contains only one variable.

OUTPUT:

  • Polynomial over \QQ

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: (var1, var2) = var('var1, var2')
sage: C.weight_enumerator((var1, var2))
var1^7 + 7*var1^4*var2^3 + 7*var1^3*var2^4 + var2^7
sage: C.weight_enumerator(var1, var2)
var1^7 + 7*var1^4*var2^3 + 7*var1^3*var2^4 + var2^7
zero()

Returns the zero vector of self.

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3)
sage: C.zero()
(0, 0, 0, 0, 0, 0, 0)
sage: C.sum(()) # indirect doctest
(0, 0, 0, 0, 0, 0, 0)
sage: C.sum((C.gens())) # indirect doctest
(1, 1, 1, 1, 1, 1, 1)
zeta_function(name='T')

Returns 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()
(2/5*T^2 + 2/5*T + 1/5)/(2*T^2 - 3*T + 1)
zeta_polynomial(name='T')

Returns 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()
2/5*T^2 + 2/5*T + 1/5
sage: C = best_known_linear_code(6,3,GF(2))  # optional - gap_packages (Guava package)
sage: C.minimum_distance()                   # optional - gap_packages (Guava package)
3
sage: C.zeta_polynomial()                    # optional - gap_packages (Guava package)
2/5*T^2 + 2/5*T + 1/5
sage: C = codes.HammingCode(GF(2), 4)
sage: C.zeta_polynomial()
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()
1

REFERENCES:

  • I. Duursma, “From weight enumerators to zeta functions”, in Discrete Applied Mathematics, vol. 111, no. 1-2, pp. 55-73, 2001.
class sage.coding.linear_code.LinearCode(generator, d=None)

Bases: sage.coding.linear_code.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 – (optional, 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
Linear code of length 7, dimension 4 over Finite Field of size 2
sage: C.base_ring()
Finite Field of size 2
sage: C.dimension()
4
sage: C.length()
7
sage: C.minimum_distance()
3
sage: C.spectrum()
[1, 0, 0, 7, 7, 0, 0, 1]
sage: 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()
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)
sage: C
Linear code of length 7, dimension 4 over Finite Field of size 5

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

sage: C = codes.GeneralizedReedSolomonCode(GF(23).list(), 12)
sage: LinearCode(C)
Linear code of length 23, dimension 12 over Finite Field of size 23

Another example:

sage: C = codes.HammingCode(GF(7), 3)
sage: C
[57, 54] Hamming Code over Finite Field of size 7
sage: LinearCode(C)
Linear code of length 57, dimension 54 over Finite Field of size 7

AUTHORS:

  • David Joyner (11-2005)
  • Charles Prior (03-2016): trac ticket #20198, LinearCode from a code
generator_matrix(encoder_name=None, **kwargs)

Returns 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]
sage.coding.linear_code.LinearCodeFromVectorSpace(V, d=None)

Simply converts a vector subspace V of GF(q)^n into a LinearCode.

INPUT:

  • V – The vector space
  • d – (Optional, default: None) the minimum distance of the code, if known. This is an optional parameter.

Note

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

EXAMPLES:

sage: V = VectorSpace(GF(2), 8)
sage: L = V.subspace([[1,1,1,1,0,0,0,0],[0,0,0,0,1,1,1,1]])
sage: C = LinearCodeFromVectorSpace(L)
sage: C.generator_matrix()
[1 1 1 1 0 0 0 0]
[0 0 0 0 1 1 1 1]
sage: C.minimum_distance()
4

Here, we provide the minimum distance of the code.:

sage: C = LinearCodeFromVectorSpace(L, d=4)
sage: C.minimum_distance()
4
class sage.coding.linear_code.LinearCodeGeneratorMatrixEncoder(code)

Bases: sage.coding.encoder.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()

Returns 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]
class sage.coding.linear_code.LinearCodeNearestNeighborDecoder(code)

Bases: sage.coding.decoder.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)

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)
decoding_radius()

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()
1
class sage.coding.linear_code.LinearCodeParityCheckEncoder(code)

Bases: sage.coding.encoder.Encoder

Encoder based on parity_check_matrix() for Linear codes.

It constructs the generator matrix through the parity check matrix.

INPUT:

  • code – The associated code of this encoder.
generator_matrix()

Returns 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.LinearCodeParityCheckEncoder(C)
sage: E.generator_matrix()
[1 0 0 0 0 1 1]
[0 1 0 0 1 0 1]
[0 0 1 0 1 1 0]
[0 0 0 1 1 1 1]
class sage.coding.linear_code.LinearCodeSyndromeDecoder(code, maximum_error_weight=None)

Bases: sage.coding.decoder.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 Linear code of length 9, dimension 3 over Finite Field of size 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 Linear code of length 9, dimension 3 over Finite Field of size 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)

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', 'unique'}
sage: 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', 'unique'}
sage: 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)
sage: D.decoder_type()
{'complete', 'hard-decision', 'might-error', 'unique'}
sage: D.decoding_radius()
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)

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(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 = 2)
sage: Chan = channels.StaticErrorRateChannel(C.ambient_space(), 2)
sage: c = C.random_element()
sage: r = Chan(c)
sage: c == D.decode_to_code(r)
True
decoding_radius()

Returns 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
maximum_error_weight()

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

Same as self.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
syndrome_table()

Returns 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)}
sage.coding.linear_code.best_known_linear_code(n, k, F)

Returns the best known (as of 11 May 2006) linear code of length n, dimension k over field F. The function uses the tables described in bounds_minimum_distance to construct this code.

This does not require an internet connection.

EXAMPLES:

sage: best_known_linear_code(10,5,GF(2))    # long time; optional - gap_packages (Guava package)
Linear code of length 10, dimension 5 over Finite Field of size 2
sage: gap.eval("C:=BestKnownLinearCode(10,5,GF(2))")     # long time; optional - gap_packages (Guava package)
'a linear [10,5,4]2..4 shortened code'

This means that best possible binary linear code of length 10 and dimension 5 is a code with minimum distance 4 and covering radius somewhere between 2 and 4. Use bounds_minimum_distance(10,5,GF(2)) for further details.

sage.coding.linear_code.best_known_linear_code_www(n, k, F, verbose=False)

Explains the construction of the best known linear code over GF(q) with length n and dimension k, courtesy of the www page http://www.codetables.de/.

INPUT:

  • n - Integer, the length of the code
  • k - Integer, the dimension of the code
  • F - Finite field, of order 2, 3, 4, 5, 7, 8, or 9
  • verbose - Bool (default: False)

OUTPUT:

  • Text about why the bounds are as given

EXAMPLES:

sage: L = best_known_linear_code_www(72, 36, GF(2)) # optional - internet
sage: print(L)                                      # optional - internet
Construction of a linear code
[72,36,15] over GF(2):
[1]:  [73, 36, 16] Cyclic Linear Code over GF(2)
     CyclicCode of length 73 with generating polynomial x^37 + x^36 + x^34 +
x^33 + x^32 + x^27 + x^25 + x^24 + x^22 + x^21 + x^19 + x^18 + x^15 + x^11 +
x^10 + x^8 + x^7 + x^5 + x^3 + 1
[2]:  [72, 36, 15] Linear Code over GF(2)
     Puncturing of [1] at 1

last modified: 2002-03-20

This function raises an IOError if an error occurs downloading data or parsing it. It raises a ValueError if the q input is invalid.

AUTHORS:

  • Steven Sivek (2005-11-14)
  • David Joyner (2008-03)
sage.coding.linear_code.bounds_minimum_distance(n, k, F)

Calculates a lower and upper bound for the minimum distance of an optimal linear code with word length n and dimension k over the field F.

The function returns a record with the two bounds and an explanation for each bound. The function Display can be used to show the explanations.

The values for the lower and upper bound are obtained from a table constructed by Cen Tjhai for GUAVA, derived from the table of Brouwer. See http://www.codetables.de/ for the most recent data. These tables contain lower and upper bounds for q=2 (when n <= 257), q=3 (when n <= 243), q=4 (n <= 256). (Current as of 11 May 2006.) For codes over other fields and for larger word lengths, trivial bounds are used.

This does not require an internet connection. The format of the output is a little non-intuitive. Try bounds_minimum_distance(10,5,GF(2)) for an example.

This function requires optional GAP package (Guava).

EXAMPLES:

sage: print(bounds_minimum_distance(10,5,GF(2))) # optional - gap_packages (Guava package)
rec(
  construction :=
   [ <Operation "ShortenedCode">,
      [
          [ <Operation "UUVCode">,
              [
                  [ <Operation "DualCode">,
                      [ [ <Operation "RepetitionCode">, [ 8, 2 ] ] ] ],
                  [ <Operation "UUVCode">,
                      [
                          [ <Operation "DualCode">,
                              [ [ <Operation "RepetitionCode">, [ 4, 2 ] ] ] ]
                            , [ <Operation "RepetitionCode">, [ 4, 2 ] ] ] ]
                 ] ], [ 1, 2, 3, 4, 5, 6 ] ] ],
  k := 5,
  lowerBound := 4,
  lowerBoundExplanation := ...
  n := 10,
  q := 2,
  references := rec(
       ),
  upperBound := 4,
  upperBoundExplanation := ... )
sage.coding.linear_code.code2leon(C)

Writes a file in Sage’s temp directory representing the code C, returning the absolute path to the file. This is the Sage translation of the GuavaToLeon command in Guava’s codefun.gi file.

INPUT:

  • C - a linear code (over GF(p), p < 11)

OUTPUT:

  • Absolute path to the file written

EXAMPLES:

sage: C = codes.HammingCode(GF(2), 3); C
[7, 4] Hamming Code over Finite Field of size 2
sage: file_loc = sage.coding.linear_code.code2leon(C)
sage: f = open(file_loc); print(f.read())
LIBRARY code;
code=seq(2,4,7,seq(
1,0,0,0,0,1,1,
0,1,0,0,1,0,1,
0,0,1,0,1,1,0,
0,0,0,1,1,1,1
));
FINISH;
sage: f.close()
sage.coding.linear_code.min_wt_vec_gap(Gmat, n, k, F, algorithm=None)

Returns a minimum weight vector of the code generated by Gmat.

Uses C programs written by Steve Linton in the kernel of GAP, so is fairly fast. The option algorithm="guava" requires Guava. The default algorithm requires GAP but not Guava.

INPUT:

  • Gmat - String representing a GAP generator matrix G of a linear code
  • n - Length of the code generated by G
  • k - Dimension of the code generated by G
  • F - Base field

OUTPUT:

  • Minimum weight vector of the code generated by Gmat

REMARKS:

  • The code in the default case allows one (for free) to also compute the message vector m such that m\*G = v, and the (minimum) distance, as a triple. however, this output is not implemented.
  • The binary case can presumably be done much faster using Robert Miller’s code (see the docstring for the spectrum method). This is also not (yet) implemented.

EXAMPLES:

sage: Gstr = "Z(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: sage.coding.linear_code.min_wt_vec_gap(Gstr,7,4,GF(2))
(0, 1, 0, 1, 0, 1, 0)

This output is different but still a minimum weight vector:

sage: sage.coding.linear_code.min_wt_vec_gap(Gstr,7,4,GF(2),algorithm="guava")    # optional - gap_packages (Guava package)
(0, 0, 1, 0, 1, 1, 0)

Here Gstr is a generator matrix of the Hamming [7,4,3] binary code.

TESTS:

We check that trac ticket #18480 is fixed:

sage: codes.HammingCode(GF(2), 2).minimum_distance()
3

AUTHORS:

  • David Joyner (11-2005)
sage.coding.linear_code.self_orthogonal_binary_codes(n, k, b=2, parent=None, BC=None, equal=False, in_test=None)

Returns a Python iterator which generates a complete set of representatives of all permutation equivalence classes of self-orthogonal binary linear codes of length in [1..n] and dimension in [1..k].

INPUT:

  • n - Integer, maximal length
  • k - Integer, maximal dimension
  • b - Integer, requires that the generators all have weight divisible by b (if b=2, all self-orthogonal codes are generated, and if b=4, all doubly even codes are generated). Must be an even positive integer.
  • parent - Used in recursion (default: None)
  • BC - Used in recursion (default: None)
  • equal - If True generates only [n, k] codes (default: False)
  • in_test - Used in recursion (default: None)

EXAMPLES:

Generate all self-orthogonal codes of length up to 7 and dimension up to 3:

sage: for B in self_orthogonal_binary_codes(7,3):
....:    print(B)
Linear code of length 2, dimension 1 over Finite Field of size 2
Linear code of length 4, dimension 2 over Finite Field of size 2
Linear code of length 6, dimension 3 over Finite Field of size 2
Linear code of length 4, dimension 1 over Finite Field of size 2
Linear code of length 6, dimension 2 over Finite Field of size 2
Linear code of length 6, dimension 2 over Finite Field of size 2
Linear code of length 7, dimension 3 over Finite Field of size 2
Linear code of length 6, dimension 1 over Finite Field of size 2

Generate all doubly-even codes of length up to 7 and dimension up to 3:

sage: for B in self_orthogonal_binary_codes(7,3,4):
....:    print(B); print(B.generator_matrix())
Linear code of length 4, dimension 1 over Finite Field of size 2
[1 1 1 1]
Linear code of length 6, dimension 2 over Finite Field of size 2
[1 1 1 1 0 0]
[0 1 0 1 1 1]
Linear code of length 7, dimension 3 over Finite Field of size 2
[1 0 1 1 0 1 0]
[0 1 0 1 1 1 0]
[0 0 1 0 1 1 1]

Generate all doubly-even codes of length up to 7 and dimension up to 2:

sage: for B in self_orthogonal_binary_codes(7,2,4):
....:    print(B); print(B.generator_matrix())
Linear code of length 4, dimension 1 over Finite Field of size 2
[1 1 1 1]
Linear code of length 6, dimension 2 over Finite Field of size 2
[1 1 1 1 0 0]
[0 1 0 1 1 1]

Generate all self-orthogonal codes of length equal to 8 and dimension equal to 4:

sage: for B in self_orthogonal_binary_codes(8, 4, equal=True):
....:     print(B); print(B.generator_matrix())
Linear code of length 8, dimension 4 over Finite Field of size 2
[1 0 0 1 0 0 0 0]
[0 1 0 0 1 0 0 0]
[0 0 1 0 0 1 0 0]
[0 0 0 0 0 0 1 1]
Linear code of length 8, dimension 4 over Finite Field of size 2
[1 0 0 1 1 0 1 0]
[0 1 0 1 1 1 0 0]
[0 0 1 0 1 1 1 0]
[0 0 0 1 0 1 1 1]

Since all the codes will be self-orthogonal, b must be divisible by 2:

sage: list(self_orthogonal_binary_codes(8, 4, 1, equal=True))
Traceback (most recent call last):
...
ValueError: b (1) must be a positive even integer.
sage.coding.linear_code.wtdist_gap(Gmat, n, F)