Full ECC API (autogenerated)

QuantumClifford.ECC.ConcatType

Concat(c₁, c₂) is a code concatenation of two quantum codes (Knill and Laflamme, 1996).

The inner code c₁ and the outer code c₂. The construction is the following: replace each qubit in code c₂ with logical qubits encoded by code c₁. The resulting code will have n = n₁ × n₂ qubits and k = k₁ × k₂ logical qubits.

source
QuantumClifford.ECC.DistanceMIPAlgorithmType
struct DistanceMIPAlgorithm <: QECCore.AbstractDistanceAlg

A Mixed Integer Programming (MIP) method for computing the code distance of CSS stabilizer codes by finding the minimum-weight non-trivial logical PauliOperator (either X-type or Z-type). Used with distance to select MIP as the method of finding the distance of a code.

Note
  • Requires a JuMP-compatible MIP solver (e.g., HiGHS, SCIP).
  • For some stabilizer CSS codes, the X-distance and Z-distance are equal.
  • logical_qubit: index of the logical qubit to compute code distance for (nothing means compute for all logical qubits)

  • logical_operator_type: type of logical operator to consider (:X or :Z, defaults to :minXZ).

  • solver: JuMP-compatible MIP solver (e.g., HiGHS, SCIP)

  • opt_summary: when true (default=false), prints the MIP solver's solution summary

  • time_limit: time limit (in seconds) for the MIP solver's execution (default=60.0)

source
QuantumClifford.ECC.ShorSyndromeECCSetupType

Configuration for ECC evaluators that simulate the Shor-style syndrome measurement (without a flag qubit).

The simulated circuit includes:

  • perfect noiseless encoding (encoding and its fault tolerance are not being studied here)
  • one round of "memory noise" after the encoding but before the syndrome measurement
  • perfect preparation of entangled ancillary qubits
  • noisy Shor-style syndrome measurement (only two-qubit gate noise)
  • noiseless "logical state measurement" (providing the comparison data when evaluating the decoder)

See also: CommutationCheckECCSetup, NaiveSyndromeECCSetup

source
QuantumClifford.ECC.TableDecoderType

A simple look-up table decoder for error correcting codes.

The lookup table contains only weight=1 errors, thus it is small, but at best it provides only for distance=3 decoding.

The size of the lookup table would grow exponentially quickly for higher distances.

source
QECCore.code_kMethod

The number of logical qubits in a code.

Note that when redundant rows exist in the parity check matrix, the number of logical qubits code_k(c) will be greater than code_n(c) - code_s(c), where the difference equals the redundancy.

source
QuantumClifford.ECC.evaluate_decoderMethod

Evaluate the performance of an error-correcting circuit.

This method requires you give the circuit that performs both syndrome measurements and (probably noiseless) logical state measurements. The faults matrix that translates an error vector into corresponding logical errors is necessary as well.

This is a relatively barebones method that assumes the user prepares necessary circuits, etc. It is a method that is used internally by more user-frienly methods providing automatic conversion of codes and noise models to the necessary noisy circuits.

source
QuantumClifford.ECC.faults_matrixMethod

Error-to-logical-observable map (a.k.a. fault matrix) of a code.

For a code with n physical qubits and k logical qubits this function returns a 2k × 2n binary matrix O such that O[i,j] is true if the logical observable of index i is flipped by the single physical qubit error of index j. Indexing is such that:

  • O[1:k,:] is the error-to-logical-X-observable map (logical X observable, i.e. triggered by logical Z errors)
  • O[k+1:2k,:] is the error-to-logical-Z-observable map
  • O[:,1:n] is the X-physical-error-to-logical-observable map
  • O[n+1:2n,:] is the Z-physical-error-to-logical-observable map

E.g. for k=1, n=10, then if O[2,5] is true, then the logical Z observable is flipped by a X₅ error; and if O[1,12] is true, then the logical X observable is flipped by a Z₂ error.

Of note is that there is a lot of freedom in choosing the logical operations! A logical operator multiplied by a stabilizer operator is still a logical operator. Similarly there is a different fault matrix for each choice of logical operators. But once the logical operators are picked, the fault matrix is fixed.

Below we show an example that uses the Shor code. While it is not the smallest code, it is a convenient choice to showcase the importance of the fault matrix when dealing with degenerate codes where a correction operation and an error do not need to be the same.

First, consider a single-qubit error, potential correction operations, and their effect on the Shor code:

julia> using QuantumClifford.ECC: faults_matrix, Shor9

julia> state = MixedDestabilizer(Shor9())
𝒟ℯ𝓈𝓉𝒶𝒷━━━━━
+ Z________
+ ___Z_____
+ _X_______
+ __X______
+ ____X____
+ _____X___
+ ______X__
+ _______X_
𝒳ₗ━━━━━━━━━
+ ______XXX
𝒮𝓉𝒶𝒷━━━━━━━
+ XXX___XXX
+ ___XXXXXX
+ ZZ_______
+ Z_Z______
+ ___ZZ____
+ ___Z_Z___
+ ______Z_Z
+ _______ZZ
𝒵ₗ━━━━━━━━━
+ Z__Z____Z

julia> err_Z₁ = single_z(9,1) # the error we will simulate
+ Z________

julia> cor_Z₂ = single_z(9,2) # the correction operation we will perform
+ _Z_______

julia> err_Z₁ * state # observe that one of the syndrome bits is now flipped
𝒟ℯ𝓈𝓉𝒶𝒷━━━━━
+ Z________
+ ___Z_____
+ _X_______
+ __X______
+ ____X____
+ _____X___
+ ______X__
+ _______X_
𝒳ₗ━━━━━━━━━
+ ______XXX
𝒮𝓉𝒶𝒷━━━━━━━
- XXX___XXX
+ ___XXXXXX
+ ZZ_______
+ Z_Z______
+ ___ZZ____
+ ___Z_Z___
+ ______Z_Z
+ _______ZZ
𝒵ₗ━━━━━━━━━
+ Z__Z____Z

julia> cor_Z₂ * err_Z₁ * state # we are back to a good code state
𝒟ℯ𝓈𝓉𝒶𝒷━━━━━
+ Z________
+ ___Z_____
- _X_______
+ __X______
+ ____X____
+ _____X___
+ ______X__
+ _______X_
𝒳ₗ━━━━━━━━━
+ ______XXX
𝒮𝓉𝒶𝒷━━━━━━━
+ XXX___XXX
+ ___XXXXXX
+ ZZ_______
+ Z_Z______
+ ___ZZ____
+ ___Z_Z___
+ ______Z_Z
+ _______ZZ
𝒵ₗ━━━━━━━━━
+ Z__Z____Z

julia> bad_Z₆Z₉ = single_z(9,6) * single_z(9,9) # a different "correction" operation
+ _____Z__Z

julia> bad_Z₆Z₉ * err_Z₁ * state # the syndrome is trivial, but now we have a logical error
𝒟ℯ𝓈𝓉𝒶𝒷━━━━━
+ Z________
+ ___Z_____
+ _X_______
+ __X______
+ ____X____
- _____X___
+ ______X__
+ _______X_
𝒳ₗ━━━━━━━━━
- ______XXX
𝒮𝓉𝒶𝒷━━━━━━━
+ XXX___XXX
+ ___XXXXXX
+ ZZ_______
+ Z_Z______
+ ___ZZ____
+ ___Z_Z___
+ ______Z_Z
+ _______ZZ
𝒵ₗ━━━━━━━━━
+ Z__Z____Z

The success of cor_Z₂ and the failure of bad_Z₆Z₉ can be immediately seen through the fault matrix, as the wrong "correction" does not result in the same logical flips ad the error:

julia> O = faults_matrix(Shor9())
2×18 BitMatrix:
 0  0  0  0  0  0  0  0  0  0  0  0  0  0  0  1  1  1
 1  0  0  1  0  0  0  0  1  0  0  0  0  0  0  0  0  0

julia> O * stab_to_gf2(err_Z₁)
2-element Vector{Int64}:
 0
 0

julia> O * stab_to_gf2(cor_Z₂)
2-element Vector{Int64}:
 0
 0

julia> O * stab_to_gf2(bad_Z₆Z₉)
2-element Vector{Int64}:
 1
 0

While its use in this situation is rather contrived, the fault matrix is incredibly useful when running large scale simulations in which we want a separate fast error sampling process, (e.g. with Pauli frames) and a syndrome decoding process, without coupling between them. We just gather all our syndrome measurement and logical observables from the Pauli frame simulations, and then use them with the fault matrix in the syndrome decoding simulation.

source
QuantumClifford.ECC.isdegenerateFunction

Check if the code is degenerate with respect to a given set of error or with respect to all "up to d physical-qubit" errors (defaulting to d=1).

julia> using QuantumClifford.ECC

julia> isdegenerate(Shor9(), [single_z(9,1), single_z(9,2)])
true

julia> isdegenerate(Shor9(), [single_z(9,1), single_x(9,1)])
false

julia> isdegenerate(Steane7(), 1)
false

julia> isdegenerate(Steane7(), 2)
true
source
QuantumClifford.ECC.naive_encoding_circuitMethod

Encoding physical qubits into a larger logical code.

The initial physical qubits to be encoded have to be at indices n-k+1:n.

Encoding circuits are not fault-tolerant

Encoding circuits are not fault-tolerant, and thus should not be used in practice. Instead, you should measure the stabilizers of the code and the logical observables, thus projecting into the code space (which can be fault-tolerant).

The canonicalization operation performed on the code may permute the qubits (see canonicalize_gott!). That permutation is corrected for with SWAP gates by default (controlled by the undoperm keyword argument).

Based on (Cleve and Gottesman, 1997) and (Gottesman, 1997), however it seems the published algorithm has some errors. Consult the erratum, as well as the more recent (Grassl, 2002) and (Grassl, 2011), and be aware that this implementation also uses H instead of Z gates.

source
QuantumClifford.ECC.naive_syndrome_circuitFunction

Generate the non-fault-tolerant stabilizer measurement cicuit for a given code instance or parity check tableau.

Use the ancillary_index and bit_index arguments to offset where the corresponding part the circuit starts.

Returns the circuit, the number of ancillary qubits that were added, and a list of bit indices that will store the measurement results.

See also: shor_syndrome_circuit

source
QuantumClifford.ECC.shor_syndrome_circuitFunction

Generate the Shor fault-tolerant stabilizer measurement cicuit for a given code instance or parity check tableau.

Use the ancillary_index and bit_index arguments to offset where the corresponding part the circuit starts. Ancillary qubits

Returns:

  • The ancillary cat state preparation circuit.
  • The Shor syndrome measurement circuit.
  • The number of ancillary qubits that were added.
  • The list of bit indices that store the final measurement results.

See also: naive_syndrome_circuit

source

Implemented in an extension requiring Hecke.jl

QuantumCliffordHeckeExt.LPCodeType
struct LPCode <: QECCore.AbstractECC

Lifted product codes ((Panteleev and Kalachev, 2021), (Panteleev and Kalachev, Jun 2022))

A lifted product code is defined by the hypergraph product of a base matrices A and the conjugate of another base matrix B'. Here, the hypergraph product is taken over a group algebra, of which the base matrices are consisting.

The binary parity check matrices are obtained by applying A_repr and B_repr representation maps to each element of the base matrices. These linear transformations convert group algebra elements to their matrix representations while preserving the CSS orthogonality condition.

Mathematical Framework

Given classical parity-check matrices:

  • $A \in \mathbb{F}_q^{m_a \times n_a}$

  • $B \in \mathbb{F}_q^{m_b \times n_b}$

The lifted product construction produces quantum CSS codes with parity-check matrices:

\[\begin{aligned} H_X &= [A \otimes I_{m_b}, -I_{m_a} \otimes B] \\ H_Z &= [I_{n_a} \otimes B^*, A^* \otimes I_{n_b}] \end{aligned}\]

Commutative Group Algebra

When R is commutative, a single representation suffices since all elements naturally commute. Here $\rho(a) = \lambda(a)$ for all $a \in R$.

Non-Commutative Group Algebra

When R is non-commutative, distinct representations are essential:

  • A_repr implements the right regular representation: $\rho(a)x = xa$

  • B_repr implements the left regular representation: $\lambda(b)x = bx$

These ensure the critical commutation relation:

\[\begin{aligned} \rho(a)\lambda(b) = \lambda(b)\rho(a) \end{aligned}\]

which follows from the associative property:

\[\begin{aligned} \rho(a)\lambda(b)(x) = b(xa) = (bx)a = \lambda(b)\rho(a)(x) \end{aligned}\]

Constructors

Multiple constructors are available:

  1. Two base matrices of group algebra elements.

  2. Two lifted codes, whose base matrices are for quantum code construction.

  3. Two base matrices of group elements, where each group element will be considered as a group algebra element by assigning a unit coefficient.

  4. Two base matrices of integers, where each integer represent the shift of a cyclic permutation. The order of the cyclic permutation should be specified.

Below is a list of all constructors:

LPCode(A, B; GA, A_repr, B_repr, repr)

defined at /home/runner/work/QuantumClifford.jl/QuantumClifford.jl/ext/QuantumCliffordHeckeExt/lifted_product.jl:158.

LPCode(c₁, c₂; GA, A_repr, B_repr, repr)

defined at /home/runner/work/QuantumClifford.jl/QuantumClifford.jl/ext/QuantumCliffordHeckeExt/lifted_product.jl:175.

LPCode(A, B; GA)

defined at /home/runner/work/QuantumClifford.jl/QuantumClifford.jl/ext/QuantumCliffordHeckeExt/lifted_product.jl:195.

LPCode(group_elem_array1, group_elem_array2; GA)

defined at /home/runner/work/QuantumClifford.jl/QuantumClifford.jl/ext/QuantumCliffordHeckeExt/lifted_product.jl:203.

LPCode(shift_array1, shift_array2, l; GA)

defined at /home/runner/work/QuantumClifford.jl/QuantumClifford.jl/ext/QuantumCliffordHeckeExt/lifted_product.jl:211.

Examples

A [[882, 24, d ≤ 24]] code from Appendix B of (Roffe et al., 2023). We use the 1st constructor to generate the code and check its length and dimension. During the construction, we do arithmetic operations to get the group algebra elements in base matrices A and B. Here x is the generator of the group algebra, i.e., offset-1 cyclic permutation, and GA(1) is the unit element.

julia> import Hecke: group_algebra, GF, abelian_group, gens; import LinearAlgebra: diagind; using QuantumClifford.ECC;

julia> l = 63; GA = group_algebra(GF(2), abelian_group(l)); x = gens(GA)[];

julia> A = zeros(GA, 7, 7);

julia> A[diagind(A)] .= x^27;

julia> A[diagind(A, -1)] .= x^54;

julia> A[diagind(A, 6)] .= x^54;

julia> A[diagind(A, -2)] .= GA(1);

julia> A[diagind(A, 5)] .= GA(1);

julia> B = reshape([1 + x + x^6], (1, 1));

julia> c1 = LPCode(A, B);

julia> code_n(c1), code_k(c1)
(882, 24)

A [[175, 19, d ≤ 10]] code from Eq. (18) in Appendix A of (Raveendran et al., 2022), following the 4th constructor.

julia> import Hecke; using QuantumClifford.ECC;

julia> base_matrix = [0 0 0 0; 0 1 2 5; 0 6 3 1]; l = 7;

julia> c2 = LPCode(base_matrix, l .- base_matrix', l);

julia> code_n(c2), code_k(c2)
(175, 19)

Code subfamilies and convenience constructors for them

  • When the base matrices of the LPCode are 1×1, the code is called a two-block group-algebra code two_block_group_algebra_codes.
  • When the base matrices of the LPCode are 1×1 and their elements are sums of cyclic permutations, the code is called a generalized bicycle code generalized_bicycle_codes.
  • When the two matrices are adjoint to each other, the code is called a bicycle code bicycle_codes.

The representation function

We use the default representation function Hecke.representation_matrix to convert a GF(2)-group algebra element to a binary matrix. The default representation, provided by Hecke, is the permutation representation.

We also accept a custom representation function as detailed in LiftedCode.

See also: LiftedCode, two_block_group_algebra_codes, generalized_bicycle_codes, bicycle_codes, haah_cubic_codes.

All fields:

  • A::Union{LinearAlgebra.Adjoint{<:Hecke.GroupAlgebraElem, <:Matrix{<:Hecke.GroupAlgebraElem}}, Matrix{<:Hecke.GroupAlgebraElem}}: the first base matrix of the code, whose elements are in a group algebra.

  • B::Union{LinearAlgebra.Adjoint{<:Hecke.GroupAlgebraElem, <:Matrix{<:Hecke.GroupAlgebraElem}}, Matrix{<:Hecke.GroupAlgebraElem}}: the second base matrix of the code, whose elements are in the same group algebra as A.

  • GA::Hecke.GroupAlgebra: the group algebra for which elements in A and B are from.

  • A_repr::Function: a function that converts a group algebra element to a binary matrix for A; default to be the right regular representation for GF(2)-algebra.

  • B_repr::Function: a function that converts a group algebra element to a binary matrix for B; default to be the left regular representation for GF(2)-algebra.

source
QuantumCliffordHeckeExt.LaCrossType
struct LaCross <: QECCore.AbstractCSSCode

The La-cross code is a quantum LDPC code constructed using the hypergraph product of two classical seed LDPC codes. It is characterized by its parity check matrix H, which is derived from circulant matrices with specific properties. These codes were introduced in (Pecorari et al., 2025).

The La-cross code has two families, distinguished by their boundary conditions:

  • Periodic boundary (full_rank = false):

\[\begin{aligned} \left[\left[2n^2, 2k^2, d \right]\right] \end{aligned}\]

  • Open boundary (full_rank = true):

\[\begin{aligned} \left[\left[(n - k)^2 + n^2, k^2, d \right]\right] \end{aligned}\]

Note

When H is square and circulant (full_rank=false), classical checks connect opposite endpoints of the length-n classical code and give rise to a quantum code with stabilizers connecting opposite array boundaries, i.e. with periodic boundary conditions. On the contrary, rectangular parity-check matrices in $\mathbb{F}_2^{(n−k)×n}$ give rise to a quantum code with stabilizers stretching up to the array extent, i.e. with open boundary conditions.

Cyclic Codes and Circulant Matrices

A cyclic code is a linear code in which codewords remain valid under cyclic shifts. A circulant matrix is a square matrix where each row is a cyclic shift of the first row. When the parity-check matrix H is circulant, the code is fully determined by its first row:

\[\begin{aligned} H = \text{circ}(c_0, c_1, \dots, c_k, 0, \dots, 0) \in \mathbb{F}_2^{n \times n}. \end{aligned}\]

The elements $c_i$ (for $i = 0, 1, \dots, k$) correspond to the coefficients of a $degree-k$ polynomial:

\[\begin{aligned} h(x) = 1 + \sum_{i=1}^{k} c_i x^i \in \mathbb{F}_2[x]/(x^n - 1). \end{aligned}\]

This establishes a mapping between $\mathbb{F}_2^n$ and the quotient ring $\mathbb{F}_2[x]/(x^n - 1)$, where cyclic shifts in $\mathbb{F}_2^n$ correspond to multiplications by x in the polynomial ring. Since multiplication by x preserves the ideal structure of $\mathbb{F}_2[x]/(x^n - 1)$, cyclic codes correspond to ideals in this ring. These ideals are in one-to-one correspondence with unitary $mod-2$ divisors of $x^n - 1$ with a leading coefficient of 1. Consequently, the fundamental building blocks of cyclic codes correspond to the factorization of $x^n - 1$.

Note

For k = 1, the generator polynomial $h(x) = 1 + x$ defines the repetition code.

Polynomial Representation

The first row of a circulant matrix $H = \text{circ}(c_0, c_1, c_2, \dots, c_{n-1})$ can be mapped to the coefficients of a polynomial $h(x)$. For instance, if the first row is $[1, 1, 0, 1]$, the polynomial is: $h(x) = 1 + x + x^3$. This polynomial-based representation aids in the analysis and design of cyclic codes. For our implementation of La-cross codes, we leverage Hecke.polynomial_ring to work directly with polynomial rings rather than manipulating coefficient arrays explicitly.

Note

The next-to-next-to-nearest neighbor connectivity implies the use of a degree-3 seed polynomial $h(x) = 1 + x + x^2 + x^3$ in the ring $\mathbb{F}_2[x]/(x^n - 1)$ for a specific code length n. Additionally, the condition of low stabilizer weight requires the polynomial $1 + x + x^3$.

[[2n², 2k², d]] La-cross Code

Here is [[98, 18, 4]] La-cross code from with $h(x) = 1 + x + x^3$, n = 7, and k = 3 from Appendix A of (Pecorari et al., 2025).

julia> using QuantumClifford; using QuantumClifford.ECC; # hide

julia> import Hecke: GF, polynomial_ring;

julia> n = 7; k = 3; F = GF(2);

julia> R, x = polynomial_ring(F, "x");

julia> h = 1 + x + x^k;

julia> c = LaCross(n, h, false);

julia> import HiGHS;

julia> code_n(c), code_k(c), distance(c, DistanceMIPAlgorithm(solver=HiGHS))
(98, 18, 4)

[[(n - k)² + n², k², d]] La-cross Code

Here is [[65, 9, 4]] La-cross code from with $h(x) = 1 + x + x^3$, n = 7, k = 3 and full rank seed rectangular circulant matrix from Appendix A of (Pecorari et al., 2025).

julia> c = LaCross(n, h, true);

julia> code_n(c), code_k(c), distance(c, DistanceMIPAlgorithm(solver=HiGHS))
(65, 9, 4)

Here is [[400, 16, 8]] La-cross code from with $h(x) = 1 + x + x^4$, n = 16, k = 4 from (Pecorari et al., 2025).

julia> n = 16; k = 4;

julia> R, x = polynomial_ring(F, "x");

julia> h = 1 + x + x^k;

julia> full_rank = true;

julia> c = LaCross(n, h, full_rank);

julia> code_n(c), code_k(c), distance(c, DistanceMIPAlgorithm(solver=HiGHS))
(400, 16, 8)

The ECC Zoo has an entry for this family.

Fields

  • n::Int64: The block length of the classical seed code.

  • h::Nemo.FqPolyRingElem: The seed vector is represented with a degree-n polynomial of the form $h(x) = \sum_{i=0}^{n-1} c_i x^i$.

  • full_rank::Bool: A flag indicating whether to use the full-rank rectangular matrix (true) or the original circulant matrix (false).

source
QuantumCliffordHeckeExt.LiftedCodeType
struct LiftedCode <: QECCore.AbstractCECC

Classical codes lifted over a group algebra, used for lifted product code construction ((Panteleev and Kalachev, 2021), (Panteleev and Kalachev, Jun 2022))

The parity-check matrix is constructed by applying repr to each element of A, which is mathematically a linear map from a group algebra element to a binary matrix. The size of the parity check matrix will enlarged with each element of A being inflated into a matrix. The procedure is called a lift (Panteleev and Kalachev, Jun 2022).

Constructors

A lifted code can be constructed via the following approaches:

  1. A matrix of group algebra elements.

  2. A matrix of group elements, where a group element will be considered as a group algebra element by assigning a unit coefficient.

  3. A matrix of integers, where each integer represent the shift of a cyclic permutation. The order of the cyclic permutation should be specified.

The default GA is the group algebra of A[1, 1], the default representation repr is the permutation representation.

Below is a list of all constructors:

LiftedCode(A; GA, repr)

defined at /home/runner/work/QuantumClifford.jl/QuantumClifford.jl/ext/QuantumCliffordHeckeExt/lifted.jl:58.

LiftedCode(A; GA)

defined at /home/runner/work/QuantumClifford.jl/QuantumClifford.jl/ext/QuantumCliffordHeckeExt/lifted.jl:68.

LiftedCode(group_elem_array; GA, repr)

defined at /home/runner/work/QuantumClifford.jl/QuantumClifford.jl/ext/QuantumCliffordHeckeExt/lifted.jl:74.

LiftedCode(shift_array, l; GA, repr)

defined at /home/runner/work/QuantumClifford.jl/QuantumClifford.jl/ext/QuantumCliffordHeckeExt/lifted.jl:87.

The representation function repr

We use the default representation function Hecke.representation_matrix to convert a GF(2)-group algebra element to a binary matrix. The default representation, provided by Hecke, is the permutation representation.

We also accept a custom representation function (the repr field of the constructor). Whatever the representation, the matrix elements need to be convertible to Integers (e.g. permit lift(ZZ, ...)). Such a customization would be useful to reduce the number of bits required by the code construction.

For example, if we use a D4 group for lifting, our default representation will be 8×8 permutation matrices, where 8 is the group's order. However, we can find a 4×4 matrix representation for the group, e.g. by using the typical 2×2 representation and converting it into binary representation by replacing "1" with the Pauli I, and "-1" with the Pauli X matrix.

See also: QuantumClifford.ECC.LPCode.

All fields:

  • A::Union{LinearAlgebra.Adjoint{<:Hecke.GroupAlgebraElem, <:Matrix{<:Hecke.GroupAlgebraElem}}, Matrix{<:Hecke.GroupAlgebraElem}}: the base matrix of the code, whose elements are in a group algebra.

  • GA::Hecke.GroupAlgebra: the group algebra for which elements in A are from.

  • repr::Function: a function that converts a group algebra element to a binary matrix; default to be the permutation representation for GF(2)-algebra.

source
QuantumClifford.ECC.bicycle_codesMethod
bicycle_codes(a_shifts::Array{Int64}, l::Int64) -> Any

Bicycle codes are a special case of generalized bicycle codes, where a and b are conjugate to each other. The order of the cyclic group is l, and the shifts a_shifts and b_shifts are reverse to each other. Thus you need to provide only the a_shifts and the rest of the conversions and conjugations are taken care of.

See also: two_block_group_algebra_codes, generalized_bicycle_codes, haah_cubic_codes.

source
QuantumClifford.ECC.generalized_bicycle_codesMethod
generalized_bicycle_codes(
    a_shifts::Array{Int64},
    b_shifts::Array{Int64},
    l::Int64
) -> Any

Generalized bicycle codes, which are a special case of abelian 2GBA codes (and therefore of lifted product codes). Here the group is chosen as the cyclic group of order l, and the base matrices a and b are the sum of the group algebra elements corresponding to the shifts a_shifts and b_shifts.

Behind the scenes, the shifts are converted to the corresponding group algebra elements and piped to two_block_group_algebra_codes.

See also: two_block_group_algebra_codes, bicycle_codes.

Examples

A [[254, 28, 14 ≤ d ≤ 20]] code from (A1) in Appendix B of (Panteleev and Kalachev, 2021).

julia> import Hecke; using QuantumClifford.ECC

julia> c = generalized_bicycle_codes([0, 15, 20, 28, 66], [0, 58, 59, 100, 121], 127);

julia> code_n(c), code_k(c)
(254, 28)

An [[70, 8, 10]] abelian 2BGA code from Table 1 of (Lin and Pryadko, 2024), with cyclic group of order l = 35, illustrates that abelian 2BGA codes can be viewed as GB codes.

julia> import Hecke; using QuantumClifford.ECC

julia> l = 35;

julia> c1 = generalized_bicycle_codes([0, 15, 16, 18], [0, 1, 24, 27], l);

julia> code_n(c1), code_k(c1)
(70, 8)
source
QuantumClifford.ECC.haah_cubic_codesMethod
haah_cubic_codes(
    a_shifts::Array{Int64},
    b_shifts::Array{Int64},
    l::Int64
) -> Any

Haah’s cubic codes (Haah, 2011) can be viewed as generalized bicycle (GB) codes with the group G = Cₗ × Cₗ × Cₗ, where l denotes the lattice size. In particular, a GB code with the group G = ℤ₃ˣ³ corresponds to a cubic code.

Behind the scenes, this function is just a simple shortcut for preparing the group G, before piping the arguments to generalized_bicycle_codes.

The ECC Zoo has an entry for this family.

See also: bicycle_codes, generalized_bicycle_codes, two_block_group_algebra_codes.

Examples

julia> import Hecke; using QuantumClifford.ECC;

julia> c = haah_cubic_codes([0, 15, 20, 28, 66], [0, 58, 59, 100, 121], 6);

julia> code_n(c), code_k(c)
(432, 8)
source
QuantumClifford.ECC.haah_cubic_codesMethod

Haah’s cubic code is defined as $\text{LP}(1 + x + y + z, 1 + xy + xz + yz)$ where $\text{LP}$ is the lifted product code, and x, y, z are elements of the ring $R = \mathbb{F}_2[x, y, z] / (x^L - 1, y^L - 1, z^L - 1)$. Here $\mathbb{F}_2$ is the finite field of order 2 and L is the lattice size. The ring $R$ is the group algebra $\mathbb{F}_qG$ of a finite group G, where $G = (C_L)^3$ and $C_L$ is the cyclic group of order L. This method of Haah's cubic code construction is outlined in Appendix B of (Panteleev and Kalachev, Jun 2022).

Here is an example of a [[1024, 30, 13 ≤ d ≤ 32]] Haah's cubic code from Appendix B, code D of (Panteleev and Kalachev, 2021) on the 8 × 8 × 8 Lattice.

julia> import Hecke; using QuantumClifford.ECC;

julia> l = 8;

julia> c = haah_cubic_codes(l);

julia> code_n(c), code_k(c)
(1024, 30)

See also: bicycle_codes, generalized_bicycle_codes, two_block_group_algebra_codes, honeycomb_color_codes.

source
QuantumClifford.ECC.honeycomb_color_codesMethod

The honeycomb color codes (Eberhardt and Steffan, 2024) are exactly the Bivariate Bicycle (BB) codes defined by the polynomials c = 1 + x + xy and d = 1 + y + xy, provided that both and m are divisible by three. This 6.6.6 code is an example of BB code, as it represents a special case.

The ECC Zoo has an entry for this family.

julia> import Hecke; using QuantumClifford.ECC;

julia> ℓ = 9; m = 6;

julia> c = honeycomb_color_codes(ℓ, m);

julia> code_n(c), code_k(c)
(108, 4)

See also: bicycle_codes, generalized_bicycle_codes, two_block_group_algebra_codes, honeycomb_color_codes.

source
QuantumClifford.ECC.two_block_group_algebra_codesMethod
two_block_group_algebra_codes(
    a::Hecke.GroupAlgebraElem,
    b::Hecke.GroupAlgebraElem
) -> Any

Two-block group algebra (2BGA) codes, which are a special case of lifted product codes from two group algebra elements a and b, used as 1×1 base matrices. To build them, you pick a group and specific generators for that group, then you pick two polynomials made of the group generators, and then, behind the scenes, these two polynomials a and b are piped to the lifted product code constructor as the elements of 1×1 matrices.

See also: LPCode, generalized_bicycle_codes, bicycle_codes, haah_cubic_codes, honeycomb_color_codes.

Examples of 2BGA code subfamilies

C₄ x C₂

Here is an example of a [[56, 28, 2]] 2BGA code from Table 2 of (Lin and Pryadko, 2024) build out of polymonials of generators of the direct product C₄ × C₂.

julia> import Hecke: group_algebra, GF, abelian_group, gens; using QuantumClifford.ECC;

julia> GA = group_algebra(GF(2), abelian_group([14,2]));

julia> x, s = gens(GA);

julia> A = 1 + x^7;

julia> B = 1 + x^7 + s + x^8 + s*x^7 + x;

julia> c = two_block_group_algebra_codes(A,B);

julia> import HiGHS;

julia> code_n(c), code_k(c), distance(c, DistanceMIPAlgorithm(solver=HiGHS))
(56, 28, 2)

Bivariate Bicycle codes

Bivariate Bicycle codes are a class of Abelian 2BGA codes formed by the direct product of two cyclic groups ℤₗ × ℤₘ. The parameters l and m represent the orders of the first and second cyclic groups, respectively.

The ECC Zoo has an entry for this family.

A [[756, 16, ≤ 34]] code from Table 3 of (Bravyi et al., 2024):

julia> import Hecke: group_algebra, GF, abelian_group, gens; using QuantumClifford.ECC;

julia> l=21; m=18;

julia> GA = group_algebra(GF(2), abelian_group([l, m]));

julia> x, y = gens(GA);

julia> A = x^3 + y^10 + y^17;

julia> B = y^5 + x^3  + x^19;

julia> c = two_block_group_algebra_codes(A,B);

julia> code_n(c), code_k(c)
(756, 16)

Multivariate Bicycle code

The group algebra of the qubit multivariate bicycle (MB) code with r variables is 𝔽₂[𝐺ᵣ], where 𝐺ᵣ = ℤ/l₁ × ℤ/l₂ × ... × ℤ/lᵣ.

A [[48, 4, 6]] Weight-6 TB-QLDPC code from Appendix A Table 2 of (Voss et al., 2024).

julia> import Hecke: group_algebra, GF, abelian_group, gens; using QuantumClifford.ECC;

julia> l=4; m=6;

julia> GA = group_algebra(GF(2), abelian_group([l, m]));

julia> x, y = gens(GA);

julia> z = x*y;

julia> A = x^3 + y^5;

julia> B = x + z^5 + y^5 + y^2;

julia> c = two_block_group_algebra_codes(A, B);

julia> import HiGHS

julia> code_n(c), code_k(c), distance(c, DistanceMIPAlgorithm(solver=HiGHS))
(48, 4, 6)

Coprime Bivariate Bicycle code

The coprime bivariate bicycle (BB) codes are defined by two polynomials 𝑎(𝑥,𝑦) and 𝑏(𝑥,𝑦), where 𝑙 and 𝑚 are coprime, and can be expressed as univariate polynomials 𝑎(𝜋) and 𝑏(𝜋), with generator 𝜋 = 𝑥𝑦. They can be viewed as a special case of Lifted Product construction based on abelian group ℤₗ x ℤₘ where ℤⱼ cyclic group of order j.

[[108, 12, 6]] coprime-bivariate bicycle (BB) code from Table 2 of (Wang and Mueller, 2024).

julia> import Hecke: group_algebra, GF, abelian_group, gens; using QuantumClifford.ECC;

julia> l=2; m=27;

julia> GA = group_algebra(GF(2), abelian_group([l*m]));

julia> 𝜋 = gens(GA)[1];

julia> A = 𝜋^2 + 𝜋^5  + 𝜋^44;

julia> B = 𝜋^8 + 𝜋^14 + 𝜋^47;

julia> c = two_block_group_algebra_codes(A, B);

julia> import HiGHS

julia> code_n(c), code_k(c), distance(c, DistanceMIPAlgorithm(solver=HiGHS))
(108, 12, 6)
source

Implemented in an extension requiring Oscar.jl

QuantumCliffordOscarExt.DDimensionalSurfaceCodeType
struct DDimensionalSurfaceCode <: QuantumCliffordOscarExt.DDimensionalCode

Constructs the D-dimensional surface code using chain complexes and $\mathbb{F}_2$-homology.

Homological Algebra Foundations of Quantum Error Correction

The theory of chain complexes over $\mathbb{F}_2$ provides a unified framework for understanding error-correcting codes, where classical $[n, k, d]$ codes correspond to 2-term complexes and quantum CSS codes arise naturally as 3-term complexes satisfying the commutativity condition $H_Z^T H_X = 0$. The homological framework reveals that:

  • Quantum CSS codes arise from chain complexes where boundary operators correspond to parity checks.
  • Logical operators correspond to homology classes.
  • Higher-dimensional codes can be constructed through products of complexes.

Chain Complex Structure

A chain complex $C$ of length $n$ is a sequence of finite-dimensional vector spaces $C_j$ over $\mathbb{F}_2$ connected by boundary operators that are linear transformations $\partial_j \colon C_j \to C_{j-1}$:

\[\begin{aligned} C : \{0\} \longrightarrow C_n \xrightarrow{\partial_n} C_{n-1} \xrightarrow{\partial_{n-1}} \cdots \xrightarrow{\partial_2} C_1 \xrightarrow{\partial_1} C_0 \longrightarrow \{0\} \end{aligned}\]

with boundary operators satisfying $\partial_i \circ \partial_{i+1} = 0$ (equivalently, $\text{im}\, \partial_{j+1} \subseteq \ker \partial_j)$.We define:

  • i-chains: Elements of $C_i$
  • i-cycles: $Z_i(C) := \ker \partial_i$
  • i-boundaries: $B_i(C) := \mathrm{im} \partial_{i+1}$
  • i-th homology: $H_i(C) := Z_i(C)/B_i(C)$

Cohomology of the Dual Complex

Given a chain complex $C$ with boundary operators $\partial_j$, its cochain complex $\widetilde{C}$ is the dual sequence with coboundary maps $\delta^i = \partial_{i+1}^T$:

\[\begin{aligned} \widetilde{C} : \{0\} \longleftarrow C_0^* \xleftarrow{\partial_1^T} C_1^* \xleftarrow{\partial_2^T} \cdots \xleftarrow{\partial_n^T} C_n^* \longleftarrow \{0\} \end{aligned}\]

where

  • i-cocycles: $Z^i(\widetilde{C}) := \ker \partial_{i+1}^T$
  • i-coboundaries: $B^i(\widetilde{C}) := \mathrm{im} \partial_i^T$
  • i-th cohomology: $H^i(\widetilde{C}) := Z^i(\widetilde{C})/B^i(\widetilde{C})$
Note

The cohomology group $H^i(\widetilde{C}) = H(\partial_{i+1}^T, \partial_i^T)$ has the same rank as the homology group $H_i(C)$ and corresponds to $X$-type logical operators in the CSS code $CSS(\partial_i, \partial_{i+1}^T)$, while $H_i(C)$ corresponds to $Z$-type operators.

Classical Codes via Chain Complexes and $\mathbb{F_2}$ Homology

An $[n,k,d]$ classical code corresponds to a 2-term complex:

\[\begin{aligned} \{0\} \longrightarrow C_1 \xrightarrow{\partial_1 = H} C_0 \longrightarrow \{0\} \end{aligned}\]

where

  • $C_1 = \mathbb{F}_2^n$ (codeword space)
  • $C_0 = \mathbb{F}_2^{n-k}$ (syndrome space)
  • $H$ is the parity check matrix

Quantum CSS Codes via Chain Complexes and $\mathbb{F_2}$ Homology

Quantum CSS codes extend this to 3-term complexes:

\[\begin{aligned} \{0\} \longrightarrow C_2 \xrightarrow{\partial_2 = H_Z^T} C_1 \xrightarrow{\partial_1 = H_X} C_0 \longrightarrow \{0\} \end{aligned}\]

where

  • $C_1 = \mathbb{F}_2^n$ (physical qubits)
  • $C_2 = \mathbb{F}_2^{m_Z}$ (Z-stabilizers)
  • $C_0 = \mathbb{F}_2^{m_X}$ (X-stabilizers)

with the condition $\partial_1 \partial_2 = H_Z^TH_X = 0$ ensuring that CSS orthogonality is satisfied.

For any chain complex, selecting two consecutive boundary operators defines a valid CSS code. When qubits are identified with the space $C_i$, the code parameters are:

  • number of physical qubits: $n = \dim C_i$
  • number of logical qubits: $k = \dim H_i(C) = \dim H^i(C)$
  • code distance: $d = \min\{\text{wt}(v) | v \in (H_i(C) \cup H^i(C))\backslash\{0\}\}$
Note

Quantum error-correcting codes, which are represented as 3-term chain complexes, can be constructed by applying the homological or hypergraph product to two 2-term chain complexes.

For a detailed explanation, see the ECC Zoo's writeup on the CSS-to-homology correspondence.

D-dimensional Surface Code ((Berthusen et al., 2024), (Zeng and Pryadko, 2019))

We provide an explicit construction of the D-dimensional surface code within the framework of chain complexes and homology over $\mathbb{F_2}$.

The quantum code is obtained by applying the homological product (or hypergraph product) to two 2-term chain complexes. Our construction relies on taking the hypergraph product of these complexes.

Double Complex

Given chain complexes C and D, we construct a double complex derived from the tensor product of two 2-term chain complexes:

\[\begin{aligned} C \boxtimes D \quad \text{with} \quad \partial_i^v = \partial_i^C \otimes I_{D_i} \quad \text{and} \quad \partial_i^h = I_{C_i} \otimes \partial_i^D \end{aligned}\]

Total Complex

The total complex is derived from a double complex by taking the direct sum of vector spaces and boundary maps that share the same dimension:

\[\begin{aligned} \text{Tot}(C \boxtimes D)_i = \bigoplus_{i=j+k} C_j \otimes D_k = E_i \end{aligned}\]

with boundary maps:

\[\begin{aligned} \partial_i^E = \bigoplus_{i=j+k} \partial_j^v \oplus \partial_k^h \end{aligned}\]

The resulting chain complex, called the tensor product of $C$ and $D$, $C ⊗ D$, enables the construction of a CSS code when selecting any three consecutive terms in its sequence.

Subfamilies

[[L² + (L − 1)², 1, L]] 2D Surface Code

The 2D surface code is constructed using the hypergraph product of two repetition codes.Thus, we obtain a new 3-term chain complex:

\[\begin{aligned} E_2 \xrightarrow{\partial_2^E} E_1 \xrightarrow{\partial_1^E} E_0 \end{aligned}\]

Examples

julia> using Oscar; using QuantumClifford; using QuantumClifford.ECC;

julia> D = 2; L = 2;

julia> c = DDimensionalSurfaceCode(D, L);

julia> code = parity_checks(c)
+ X_X_X
+ _X_XX
+ ZZ__Z
+ __ZZZ

julia> import HiGHS; import JuMP;

julia> code_n(c), code_k(c), distance(c, DistanceMIPAlgorithm(solver=HiGHS))
(5, 1, 2)

When L = 4, we get [[25,1, 4]] 2D surface code from (Berthusen et al., 2024).

julia> using Oscar; using QuantumClifford; using QuantumClifford.ECC;

julia> D = 2; L = 4;

julia> c = DDimensionalSurfaceCode(D, L);

julia> code = parity_checks(c)
+ X___X___________X________
+ _X___X__________XX_______
+ __X___X__________XX______
+ ___X___X__________X______
+ ____X___X__________X_____
+ _____X___X_________XX____
+ ______X___X_________XX___
+ _______X___X_________X___
+ ________X___X_________X__
+ _________X___X________XX_
+ __________X___X________XX
+ ___________X___X________X
+ ZZ______________Z________
+ _ZZ______________Z_______
+ __ZZ______________Z______
+ ____ZZ__________Z__Z_____
+ _____ZZ__________Z__Z____
+ ______ZZ__________Z__Z___
+ ________ZZ_________Z__Z__
+ _________ZZ_________Z__Z_
+ __________ZZ_________Z__Z
+ ____________ZZ________Z__
+ _____________ZZ________Z_
+ ______________ZZ________Z

julia> import HiGHS; import JuMP;

julia> code_n(c), code_k(c), distance(c, DistanceMIPAlgorithm(solver=HiGHS))
(25, 1, 4)

Chain Complex

The construction is as follows:

\[\begin{aligned} C = \left( C_1 \xrightarrow{\partial} C_0 \right) \quad \text{and} \quad D = \left( D_1 \xrightarrow{\partial^T} D_0 \right) \end{aligned}\]

where $\partial$ is the $(L-1) \times L$ parity check matrix:

\[\begin{aligned} H = \begin{pmatrix} 1 & 1 & & \\ & 1 & \ddots & \\ & & \ddots & 1 \\ & & & 1 \end{pmatrix} \end{aligned}\]

[[L³ + 2L(L − 1)², 1, min(L, L²)]] 3D Surface Code

The 3D surface code is obtained by taking the hypergraph product of a 2D surface code with a repetition code. Thus, we obtain a new 4-term chain complex:

\[\begin{aligned} F_3 \xrightarrow{\partial_3^F} F_2 \xrightarrow{\partial_2^F} F_1 \xrightarrow{\partial_1^F} F_0 \end{aligned}\]

Metachecks

  • Z-type metachecks: $M_Z^T = \partial_3^F$

Example

Here is an example of [[12, 1, 2]] 3D Surface code with L = 2 from (Berthusen et al., 2024).

julia> using Oscar; using QuantumClifford; using QuantumClifford.ECC;

julia> D = 3; L = 2;

julia> c = DDimensionalSurfaceCode(D, L);

julia> code = parity_checks(c)
+ XX__X_____X_
+ __XXX______X
+ _____XX__XX_
+ _______XXX_X
+ Z_Z_Z_______
+ _Z_ZZ_______
+ _____Z_Z_Z__
+ ______Z_ZZ__
+ Z____Z____Z_
+ _Z____Z___Z_
+ __Z____Z___Z
+ ___Z____Z__Z
+ ____Z____ZZZ

julia> code_n(c), code_k(c)
(12, 1)
Note

For the 3D surface code, there is an asymmetry between the Z- and X-bases (Berthusen et al., 2024). Specifically, the Z-distance ($d_Z$) is 4, whereas the X-distance ($d_X$) is 2. As a result, the code has the parameters [[12, 1, 2]].

julia> import HiGHS; import JuMP;

julia> dz = distance(c, DistanceMIPAlgorithm(solver=HiGHS, logical_operator_type=:Z))
4

julia> dx = distance(c, DistanceMIPAlgorithm(solver=HiGHS, logical_operator_type=:X))
2

[[6L⁴ − 12L³ + 10L² − 4L + 1, 1, L²]] 4D Surface Code

The 4D surface code is constructed by taking the hypergraph product of a 3D surface code with a repetition code. Thus, we obtain a new 5-term chain complex:

\[\begin{aligned} G_4 \xrightarrow{\partial_4^G} G_3 \xrightarrow{\partial_3^G} G_2 \xrightarrow{\partial_2^G} G_1 \xrightarrow{\partial_1^G} G_0 \end{aligned}\]

[[33, 1, 4]]

Here is an example of [[33, 1, 4]] 4D Surface code with L = 2 from (Berthusen et al., 2024).

julia> using Oscar; using QuantumClifford; using QuantumClifford.ECC;

julia> D = 4; L = 2;

julia> c = DDimensionalSurfaceCode(D, L);

julia> import HiGHS; import JuMP;

julia> code_n(c), code_k(c), distance(c, DistanceMIPAlgorithm(solver=HiGHS))
(33, 1, 4)

Metachecks

Both X and Z-type metachecks available:

  • $M_Z^T = \partial_4^G$
  • $M_X = \partial_1^G$

To obtain surface codes of greater dimensionality, we alternate between C and D and then form a product with the chain complex representing the DDimensionalSurfaceCode (Berthusen et al., 2024).

Note

The procedure described above for the DDimensionalSurfaceCode can alternatively be performed using an L × L repetition code and only the chain complex C. In this case, the result would be the DDimensionalToricCode.

See also: DDimensionalToricCode

Fields

  • D::Int64: Dimension of the Surface code (must be ≥ 2).

  • L::Int64: Size parameter determining the D-dimensional Surface code family, constructed via hypergraph product of (L - 1) × L repetition code chain complexes.

source
QuantumCliffordOscarExt.DDimensionalToricCodeType
struct DDimensionalToricCode <: QuantumCliffordOscarExt.DDimensionalCode

D-dimensional Toric Code ((Berthusen et al., 2024), (Zeng and Pryadko, 2019))

The D-dimensional toric code is obtained by taking the iterated tensor product of a single chain complex:

\[\begin{aligned} C = \left( \mathbb{F}_2^L \xrightarrow{H} \mathbb{F}_2^L \right) \end{aligned}\]

where H is the $L \times L$ parity check matrix of the repetition code. The total complex is built by taking the tensor product $C^{\otimes D}$ and forming the associated total complex via direct sums.

Note

D-dimensional toric code construction differs from the surface code as we use the full $L \times L$ repetition code (rather than $(L-1) \times L$). The tensor products with identical complexes (rather than alternating complexes) are used.

Subfamilies

[[2L², 2, L]] 2D Toric Code

julia> using Oscar; using QuantumClifford; using QuantumClifford.ECC;

julia> D = 2; L = 2;

julia> c = DDimensionalToricCode(D, L);

julia> code = parity_checks(c)
+ X_X_XX__
+ _X_XXX__
+ X_X___XX
+ _X_X__XX
+ ZZ__Z_Z_
+ ZZ___Z_Z
+ __ZZZ_Z_
+ __ZZ_Z_Z

julia> import HiGHS; import JuMP;

julia> code_n(c), code_k(c), distance(c, DistanceMIPAlgorithm(solver=HiGHS))
(8, 2, 2)

[[3L³, 3, min(L, L²)]] 3D Toric Code

julia> using Oscar; using QuantumClifford; using QuantumClifford.ECC;

julia> D = 3; L = 2;

julia> c = DDimensionalToricCode(D, L);

julia> import HiGHS; import JuMP;

julia> code_n(c), code_k(c), distance(c, DistanceMIPAlgorithm(solver=HiGHS))
(24, 3, 2)

[[6L⁴, 6, L²]] 4D Toric Code

julia> using Oscar; using QuantumClifford; using QuantumClifford.ECC;

julia> D = 4; L = 2;

julia> c = DDimensionalToricCode(D, L);

julia> import HiGHS; import JuMP;

julia> code_n(c), code_k(c), distance(c, DistanceMIPAlgorithm(solver=HiGHS))
(96, 6, 4)

See also: DDimensionalSurfaceCode

Fields

  • D::Int64: Dimension of the Toric code (must be ≥ 2).

  • L::Int64: Size parameter determining the D-dimensional Toric code family, constructed via hypergraph product of L × L repetition code chain complexes.

source
QuantumClifford.ECC.boundary_mapsMethod

Returns all boundary maps of the chain complex, including both parity check and metacheck matrices.

Here are the boundary maps of [[12, 1, 2]] 3D Surface code with L = 2 from (Berthusen et al., 2024).

julia> using Oscar; using QuantumClifford; using QECCore;

julia> using QuantumClifford.ECC: DDimensionalSurfaceCode, boundary_maps, metacheck_matrix_z;

julia> D = 3; L = 2;

julia> c = DDimensionalSurfaceCode(D, L);

julia> Mz, Hz, Hx = boundary_maps(c);

The parity check matrices of [[12, 1, 2]] 3D Surface code are

julia> Hx
4×12 Matrix{Int64}:
 1  1  0  0  1  0  0  0  0  0  1  0
 0  0  1  1  1  0  0  0  0  0  0  1
 0  0  0  0  0  1  1  0  0  1  1  0
 0  0  0  0  0  0  0  1  1  1  0  1
Note

For 3D and higher-dimensional codes, Oscar returns Z-type parity check matrix as transpose ($H_Z^T$). We transpose it to convert it back to $H_Z$. See B3, page 11 of (Berthusen et al., 2024).

julia> Hz'
9×12 adjoint(::Matrix{Int64}) with eltype Int64:
 1  0  1  0  1  0  0  0  0  0  0  0
 0  1  0  1  1  0  0  0  0  0  0  0
 0  0  0  0  0  1  0  1  0  1  0  0
 0  0  0  0  0  0  1  0  1  1  0  0
 1  0  0  0  0  1  0  0  0  0  1  0
 0  1  0  0  0  0  1  0  0  0  1  0
 0  0  1  0  0  0  0  1  0  0  0  1
 0  0  0  1  0  0  0  0  1  0  0  1
 0  0  0  0  1  0  0  0  0  1  1  1

Along with the Z- and X-type parity check matrices, we have a metacheck matrix specifically for the Z-type checks. The classical code derived from this metacheck matrix has a distance of d = 2 meaning it can identify (but not correct) a single error in the Z-type syndrome measurements. See page 12 of (Berthusen et al., 2024) for details.

julia> Mz'
2×9 adjoint(::Matrix{Int64}) with eltype Int64:
 1  0  1  0  1  0  1  0  1
 0  1  0  1  0  1  0  1  1

We can use metacheck_matrix_z directly instead of using boundary_maps.

julia> metacheck_matrix_z(c)
2×9 Matrix{Int64}:
 1  0  1  0  1  0  1  0  1
 0  1  0  1  0  1  0  1  1

Metachecks in CSS Codes

The parity-check matrices $M_Z$ and $M_X$ are called metachecks in CSS codes. These matrices emerge from the constraints imposed by boundary maps, which satisfy the condition $\partial_{i+1} \partial_i = 0$. This guarantees that:

\[\begin{aligned} M_Z H_Z = 0 \quad \text{and} \quad M_X H_X = 0, \end{aligned}\]

meaning that:

  • Valid $Z$-type syndromes must be in $\ker M_Z$
  • Valid $X$-type syndromes must be in $\ker M_X$

When measurement errors occur, they distort the syndrome vector $\mathbf{s}$, generating a detectable metasyndrome. By examining $\mathbf{m}$, we can identify and correct errors in $\mathbf{s}$ before proceeding with standard decoding. This technique is called syndrome repair decoding (Higgott and Breuckmann, 2023).

source
QuantumClifford.ECC.twobga_from_direct_productMethod
twobga_from_direct_product(
    a_elts::Vector{Hecke.GroupAlgebraElem{Nemo.FqFieldElem, Hecke.GroupAlgebra{Nemo.FqFieldElem, Oscar.DirectProductGroup, Oscar.BasicGAPGroupElem{Oscar.DirectProductGroup}}}},
    b_elts::Vector{Hecke.GroupAlgebraElem{Nemo.FqFieldElem, Hecke.GroupAlgebra{Nemo.FqFieldElem, Oscar.DirectProductGroup, Oscar.BasicGAPGroupElem{Oscar.DirectProductGroup}}}},
    F2G::Hecke.GroupAlgebra{Nemo.FqFieldElem, Oscar.DirectProductGroup, Oscar.BasicGAPGroupElem{Oscar.DirectProductGroup}}
) -> Any

Constructing two block group algebra codes by specifying the direct product to be used. See also the more general two_block_group_algebra_codes.

Two block group algebra codes are constructed by choosing a group (and specific generators), then choosing two polynomials made out of these generators, then piping these two polynomials as the elements of 1×1 matrices to the lifted product code constructors.

The Hecke library, for which we already have an extension, provides for a fairly easy way to construct such polynomials for many abelian and small groups. See two_block_group_algebra_codes for those capabilities.

However, more esoteric groups can be specified as the direct product of other groups. To support arbitrary direct products we use Oscar, which builds upon Hecke. Oscar supports the direct product operation between two or more arbitrary general groups, including non-abelian groups such as alternating_group, dihedral_group, symmetric_group, and even arbitrary finitely presented groups (e.g., free_group). This capability is not available in Hecke.jl. The 2BGA codes discovered in (Lin and Pryadko, 2024) rely on direct products of two or more general groups, which necessitate the use of Oscar.direct_product.

This particular function is nothing more than a simple wrapper that takes care of argument conversions. Of note, the polynomials here are given as lists of monomials.

Of course, if you are comfortable with Oscar, you can use two_block_group_algebra_codes directly.

See also: two_block_group_algebra_codes, twobga_from_fp_group

Examples

The [[56, 28, 2]] abelian 2BGA code from Appendix C, Table II in (Lin and Pryadko, 2024) can be constructed using the direct product of two cyclic groups. Specifically, the group C₂₈ of order l = 28 can be represented as C₁₄ × C₂, where the first group has order m = 14 and the second group has order n = 2.

julia> import Oscar: cyclic_group, small_group_identification, describe, order

julia> import Hecke: gens, quo, group_algebra, GF, one, direct_product, sub

julia> using QuantumClifford, QuantumClifford.ECC

julia> m = 14; n = 2;

julia> C₁₄ = cyclic_group(m);

julia> C₂ = cyclic_group(n);

julia> G = direct_product(C₁₄, C₂);

julia> GA = group_algebra(GF(2), G);

julia> x, s = gens(GA)[1], gens(GA)[3];

julia> a = [one(GA), x^7];

julia> b = [one(GA), x^7, s, x^8, s * x^7, x];

julia> c = twobga_from_direct_product(a, b, GA);

julia> order(G)
28

julia> code_n(c), code_k(c)
(56, 28)

julia> describe(G), small_group_identification(G)
("C14 x C2", (28, 4))
Danger

When using the direct product, there isn't necessarily a unique set of generators. It is essential to verify that Oscar is providing you with the generators you expect, e.g. for a cycling group that you have the presentation Cₘ = ⟨x, s | xᵐ = s² = xsx⁻¹s⁻¹ = 1⟩. For situations where the generators provided by Oscar are not the ones you want, you can also use twobga_from_fp_group where you specify the group presentation directly.

As a verification that you have the correct generators, Oscar.sub can be used to determine if H is a subgroup of G and to confirm that both C₁₄ and C₂ are subgroups of C₂₈.

julia> order(gens(G)[1])
14

julia> order(gens(G)[3])
2

julia> x^14 == s^2 == x * s * x^-1 * s^-1
true

julia> H, _  = sub(G, [gens(G)[1], gens(G)[3]]);

julia> H == G
true
source
QuantumClifford.ECC.twobga_from_fp_groupMethod
twobga_from_fp_group(
    a_elts::Vector{Oscar.FPGroupElem},
    b_elts::Vector{Oscar.FPGroupElem},
    F2G::Hecke.GroupAlgebra{Nemo.FqFieldElem, Oscar.FPGroup, Oscar.FPGroupElem}
) -> Any

Constructing two block group algebra codes by specifying the group presentation.

Two block group algebra codes are constructed by choosing a group (and specific generators), then choosing two polynomials made out of these generators, then piping these two polynomials as the elements of 1×1 matrices to the lifted product code constructors.

The Hecke library, for which we already have an extension, provides for a fairly easy way to construct such polynomials for many abelian and small groups. See two_block_group_algebra_codes for those capabilities.

However, more esoteric groups are usually specified by a group presentation ⟨S | R⟩, where S is a set of generators and R is the relations those generators obey. To support arbitrary group presentations we use Oscar, which builds upon Hecke. We use Oscar.free_group and quo in order to first prepare the free group generated by S, and then the group obeying also the relations R, i.e. the ⟨S | R⟩ presentation.

After that point we proceed as usual, creating two polynomials of generators and piping them to two_block_group_algebra_codes.

This particular function is nothing more than a simple wrapper that takes care of argument conversions. Of note, the polynomials here are given as lists of monomials.

Of course, if you are comfortable with Oscar, you can use two_block_group_algebra_codes directly.

See also: two_block_group_algebra_codes, twobga_from_direct_product

Examples

The [[96, 12, 10]] 2BGA code from Table I in (Lin and Pryadko, 2024) has the group presentation ⟨r, s | s⁶ = r⁸ = r⁻¹srs = 1⟩ (the group C₂ × (C₃ ⋉ C₈)).

julia> import Oscar: free_group, small_group_identification, describe, order

julia> import Hecke: gens, quo, group_algebra, GF, one

julia> using QuantumClifford, QuantumClifford.ECC

julia> F = free_group(["r", "s"]);

julia> r, s = gens(F); # generators

julia> G, = quo(F, [s^6, r^8, r^(-1) * s * r * s]);  # relations

julia> GA = group_algebra(GF(2), G);

julia> r, s = gens(G);

julia> a = [one(G), r, s^3 * r^2, s^2 * r^3];

julia> b = [one(G), r, s^4 * r^6, s^5 * r^3];

julia> c = twobga_from_fp_group(a, b, GA);

julia> order(G)
48

julia> code_n(c), code_k(c)
(96, 12)

julia> describe(G), small_group_identification(G)
("C2 x (C3 : C8)", (48, 9))

Cyclic Groups

Cyclic groups with specific group presentations, given by Cₘ = ⟨x, s | xᵐ = s² = xsx⁻¹s⁻¹ = 1⟩, where the order is 2m as seen in Table II of (Lin and Pryadko, 2024).

The [[56, 28, 2]] abelian 2BGA code from Appendix C, Table II in (Lin and Pryadko, 2024) is constructed using the group presentation ⟨x, s | xs = sx, xᵐ = s² = 1⟩ (the cyclic group C₂₈ = C₁₄ × C₂).

julia> m = 14;

julia> F = free_group(["x", "s"]);

julia> x, s = gens(F); # generators

julia> G, = quo(F, [x^m, s^2, x * s * x^-1 * s^-1]); # relations

julia> GA = group_algebra(GF(2), G);

julia> x, s = gens(G);

julia> a = [one(G), x^7];

julia> b = [one(G), x^7, s, x^8, s * x^7, x];

julia> c = twobga_from_fp_group(a, b, GA);

julia> order(G)
28

julia> code_n(c), code_k(c)
(56, 28)

julia> describe(G), small_group_identification(G)
("C14 x C2", (28, 4))

Dihedral Groups

Dihedral (non-abelian) groups with group presentations given by Dₘ = ⟨r, s | rᵐ = s² = (rs)² = 1⟩, where the order is 2m.

The [[24, 8, 3]] 2BGA code from Appendix C, Table III in (Lin and Pryadko, 2024) is constructed by specifying a group presentation below (giving the group D₆ = C₆ ⋉ C₂).

julia> m = 6;

julia> F = free_group(["r", "s"]);

julia> r, s = gens(F); # generators

julia> G, = quo(F, [r^m, s^2, (r*s)^2]); # relations

julia> GA = group_algebra(GF(2), G);

julia> r, s = gens(G);

julia> a = [one(G), r^4];

julia> b = [one(G), s*r^4, r^3, r^4, s*r^2, r];

julia> c = twobga_from_fp_group(a, b, GA);

julia> order(G)
12

julia> code_n(c), code_k(c)
(24, 8)

julia> describe(G), small_group_identification(G)
("D12", (12, 4))
Note

Notice how in all of these construction we are specifying a group presentation. We are explicitly not picking a group by name and getting its "canonical" generators, as we do not a priori know whether Oscar would give us the generating set we need (generating sets are not unique).

source

Implemented in an extension requiring JuMP.jl

QECCore.distanceMethod
distance(
    code::QECCore.AbstractECC,
    alg::QuantumClifford.ECC.DistanceMIPAlgorithm
) -> Any

Compute the distance of a code using mixed integer programming. See QuantumClifford.ECC.DistanceMIPAlgorithm for configuration options.

Specifically, it computes the minimum distance of a Calderbank-Shor-Steane code by solving two independent Mixed Integer Programs for X-type ($d_X$) and Z-type ($d_Z$) distances. The code distance is $d = \min(d_X, d_Z)$.

Background on Minimum Distance

For classical codes, the minimum distance, which measures a code's error-correcting capability, is equivalent to its minimum weight. This can be computed by generating all possible codewords from combinations of the generator matrix rows, calculating their weights, and finding the smallest. While accurate, this method takes exponential time. Vardy (Vardy, 1997) demonstrated that computing the minimum distance is NP-hard, and the corresponding decision problem is NP-complete, making polynomial-time algorithms unlikely.

For quantum codes, classical intuition does not always apply. The minimum distance is given by the minimum weight of a non-trivial logical operator. This is generally unrelated to the minimum distance of the corresponding stabilizer code when viewed as a classical, additive code. White and Grassl (White and Grassl, 2006) proposed mapping quantum codes to higher-dimensional classical linear codes. This mapping allows the minimum distance of the quantum additive code to be inferred from that of the classical linear code but increases parameters from n to 3n and d to 2d, adding complexity. Furthermore, once a minimal weight vector is identified, it is essential to verify whether it belongs to the Pauli group 𝒫ₙ over n qubits (Sabo, 2022).

Additionally, to illustrate how classical intuition can be misleading in this context, consider that the [[7, 1, 3]] Steane code has a minimum distance of three, despite all its elements having a weight of four. This discrepancy occurs because stabilizer codes are defined by parity-check matrices, while their minimum distances are determined by the dual (Sabo, 2022).

Note

The minimum distance problem for quantum codes is NP-hard, and this hardness extends to multiplicative and additive approximations, even when restricted to stabilizer or CSS codes, with the result established through a reduction from classical problems in the CWS framework using a 4-cycle free graph (Kapshikar and Kundu, 2023). Despite this, methods that improve on brute-force approaches are actively explored.

For a more in-depth background on minimum distance, see Chapter 3 of (Sabo, 2022).

Quantum Error Correction and Code Distance

The foundation of quantum error correction lies in protecting logical quantum information from physical errors by encoding it across multiple qubits. A quantum code's performance is fundamentally characterized by its distance (d), which quantifies the code's ability to detect and correct errors. Practically, the distance represents the minimum number of physical qubit errors required to cause an undetectable logical error - one that corrupts encoded information while evading the code's error detection mechanisms.

Fundamentals of Quantum Code Distance

The distance d of a quantum error-correcting code represents its robustness against physical errors and is defined as (Ekert et al., undated):

\[\begin{aligned} \begin{equation} d = \min_{P \in N(S)\setminus S} \mathrm{wt}(P) \end{equation} \end{aligned}\]

where:

  • $N(S)$ denotes the normalizer of the stabilizer group S.
  • $\mathrm{wt}(P)$ represents the weight of Pauli operator P.
  • The minimization is taken over all logical operators that are not stabilizers.

The distance reveals the essential property that it equals the smallest number of qubits that must be affected to produce a logical error undetectable by stabilizer measurements. The normalizer condition $P \in N(S)$ ensures the operator commutes with all stabilizers, while $P \notin S$ guarantees it performs a non-trivial logical operation (Ekert et al., undated).

Mixed Integer Programming

We compute the minimum code distance for CSS (Calderbank-Shor-Steane) codes by solving MIPs.

The distance is computed separately for X-type ($d_X$) and Z-type ($d_Z$) logical operators, then combined to give the true code distance: $d = \min(d_X, d_Z)$.

X-type and Z-type Distances

The X-type distance ($d_X$) and Z-type distance ($d_Z$`) are defined as the minimum number of errors required to implement a non-trivial logical operator of the opposite type, subject to the following constraints:

  • For $d_X$ (where $U = X$), the errors are Z-type (phase flips), and the constraints involve the X-stabilizer matrix $\mathbf{H_X}$ and X-logical operators $\mathbf{L_X}$. The error vector is denoted as \mathbf{e}_Z`.
  • For $d_Z$ (where $U = Z$), the errors are X-type (bit flips), with constraints given by the Z-stabilizer matrix $\mathbf{H_Z}$ and Z-logical operators $\mathbf{L_Z}$, and the error vector is $\mathbf{e}_X$.

\[\begin{aligned} \text{Minimize} \quad & \sum_{i=1}^n e_{V,i} \quad \text{(Hamming weight of errors)} \\ \text{Subject to} \quad & \mathbf{H_U} \cdot \mathbf{e}_V \equiv \mathbf{0} \pmod{2} \quad \text{(Commutes with U-stabilizers)} \\ & \mathbf{L_U} \cdot \mathbf{e}_V \equiv 1 \pmod{2} \quad \text{(Anti-commutes with U-logical)} \\ & e_{V,i} \in \{0,1\} \quad \text{(Binary error variables)} \end{aligned}\]

Example

A [[40, 8, 5]] 2BGA code with the minimum distance of 5 from Table 2 of (Lin and Pryadko, 2024).

julia> import Hecke: group_algebra, GF, abelian_group, gens; import HiGHS; import JuMP;

julia> using QuantumClifford, QuantumClifford.ECC

julia> l = 10; m = 2;

julia> GA = group_algebra(GF(2), abelian_group([l,m]));

julia> x, s = gens(GA);

julia> A = 1 + x^6;

julia> B = 1 + x^5 + s + x^6 + x + s*x^2;

julia> c = two_block_group_algebra_codes(A,B);

julia> code_n(c), code_k(c), distance(c, DistanceMIPAlgorithm(solver=HiGHS))
(40, 8, 5)

A [[48, 6, 8]] GB code with the minimum distance of 8 from (A3) in Appendix B of (Panteleev and Kalachev, 2021).

julia> l = 24;

julia> c1 = generalized_bicycle_codes([0, 2, 8, 15], [0, 2, 12, 17], l);

julia> code_n(c1), code_k(c1), distance(c1, DistanceMIPAlgorithm(solver=HiGHS))
(48, 6, 8)

Applications

Mixed-integer programming (MIP) is applied in quantum error correction, notably for decoding and minimum distance computation. Some applications are as follows:

  • The first usecase of the MIP approach was the code capacity Most Likely Error (MLE) decoder for color codes introduced in (Landahl et al., 2011).

  • For all quantum LDPC codes presented in (Panteleev and Kalachev, 2021), the lower and upper bounds on the minimum distance was obtained by reduction to a mixed integer linear program and using the GNU Linear Programming Kit ((Makhorin, 2008)).

  • For all the Bivariate Bicycle (BB) codes presented in (Bravyi et al., 2024), the code distance was calculated using the mixed integer programming approach.

  • (Lacroix et al., 2024) developed a MLE decoder that finds the most likely chain of Pauli errors given the observed error syndrome by solving a mixed-integer program using HiGHS package ((Huangfu and Hall, 2018)).

  • (Cain et al., 2025) formulate maximum-likelihood decoding as a mixed-integer program maximizing $\prod_{j=1}^M p_j^{E_j}(1-p_j)^{1-E_j}$ (where binary variables $E_j \in {0,1}$ indicate error occurrence) subject to syndrome constraints, solved optimally via MIP solvers despite its NP-hard complexity.

source