Full ECC API (autogenerated)
QuantumClifford.ECC.CSS
— TypeAn arbitrary CSS error correcting code defined by its X and Z checks.
julia> CSS([0 1 1 0; 1 1 0 0], [1 1 1 1]) |> parity_checks
+ _XX_
+ XX__
+ ZZZZ
QuantumClifford.ECC.CircuitCode
— TypeCircuitCode
is defined by a given encoding circuit circ
.
n
: qubit numbercirc
: the encoding circuitencode_qubits
: the qubits to be encoded
See also: random_all_to_all_circuit_code
, random_brickwork_circuit_code
QuantumClifford.ECC.Cleve8
— TypeA pedagogical example of a quantum error correcting [8,3] code used in (Cleve and Gottesman, 1997).
QuantumClifford.ECC.CommutationCheckECCSetup
— TypeConfiguration for ECC evaluator that does not simulate any ECC circuits, rather it simply checks the commutation of the parity check and the Pauli error.
This is much faster than any other simulation method, but it is incapable of noisy-circuit simulations and thus useless for fault-tolerance studies.
See also: NaiveSyndromeECCSetup
, ShorSyndromeECCSetup
QuantumClifford.ECC.Concat
— TypeConcat(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.
QuantumClifford.ECC.Gottesman
— TypeThe family of [[2ʲ, 2ʲ - j - 2, 3]]
Gottesman codes, also known as quantum Hamming codes, as described in Gottesman's 1997 PhD thesis and in (Gottesman, 1996).
You might be interested in consulting (Yu et al., 2013) and (Chao and Reichardt, 2017) as well.
The ECC Zoo has an entry for this family
QuantumClifford.ECC.NaiveSyndromeECCSetup
— TypeConfiguration for ECC evaluator that runs the simplest syndrome measurement circuit.
The circuit is being simulated (as opposed to doing only a quick commutation check). This circuit would give poor performance if there is non-zero gate noise.
See also: CommutationCheckECCSetup
, ShorSyndromeECCSetup
QuantumClifford.ECC.QuantumReedMuller
— TypeThe family of [[2ᵐ - 1, 1, 3]]
CSS Quantum-Reed-Muller codes, as discovered by Steane in his 1999 paper (Steane, 1999).
Quantum codes are constructed from shortened Reed-Muller codes RM(1, m)
, by removing the first row and column of the generator matrix Gₘ
. Similarly, we can define truncated dual codes RM(m - 2, m)
using the generator matrix Hₘ
(Anderson et al., 2014). The quantum Reed-Muller codes QRM(m)
derived from RM(1, m)
are CSS codes.
Given that the stabilizers of the quantum code are defined through the generator matrix of the classical code, the minimum distance of the quantum code corresponds to the minimum distance of the dual classical code, which is d = 3
, thus it can correct any single qubit error. Since one stabilizer from the original and one from the dual code are removed in the truncation process, the code parameters are [[2ᵐ - 1, 1, 3]]
.
You might be interested in consulting (Anderson et al., 2014) and (Campbell et al., 2012) as well.
The ECC Zoo has an entry for this family.
QuantumClifford.ECC.ShorSyndromeECCSetup
— TypeConfiguration 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
QuantumClifford.ECC.Surface
— TypeThe planar surface code refers to the code (Kitaev, 2003) in a 2D lattice with open boundaries.
Illustration of a 3×2 surface code, where qubits are located on the edges:
|---1--(Z)--2---|---3---|
| (X) 7 8 o
|---4---|---5---|---6---|
| o o o
| | | |
The surface code has open boundary conditions, unlike the toric code. To this end, we remove qubits (denoted by "o") and parity checks on the right and bottom sides.
Faces like (1,4,7)
have X checks, and crosses like (1,2,7)
have Z checks. Due to the removal of the bottom and right sides, we have some 3-qubit checks on the boundaries.
julia> parity_checks(Surface(3,2))
+ X__X__X_
+ _X__X_XX
+ __X__X_X
+ ZZ____Z_
+ _ZZ____Z
+ ___ZZ_Z_
+ ____ZZ_Z
More information can be seen in (Fowler et al., 2012).
QuantumClifford.ECC.TableDecoder
— TypeA 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.
QuantumClifford.ECC.Toric
— TypeThe Toric code (Kitaev, 2003).
Illustration of a 2x2 toric code, where qubits are located on the edges:
|--1-(Z)-2--|
| (X) 5 6
|--3--|--4--|
| 7 8
| | |
It is important to note that the toric code has periodic boundary conditions, which means that the top and bottom sides are essentially glued together, as are the left and right sides.
Faces like (1,3,5,6)
have X checks, and crosses like (1,2,5,7)
have Z checks.
julia> parity_checks(Toric(2,2))
+ X_X_XX__
+ _X_XXX__
+ X_X___XX
+ ZZ__Z_Z_
+ ZZ___Z_Z
+ __ZZZ_Z_
QuantumClifford.ECC.BeliefPropDecoder
— MethodA simple Belief Propagation decoder built around tools from LDPCDecoders.jl
.
QuantumClifford.ECC.BitFlipDecoder
— MethodAn Iterative Bitflip decoder built around tools from LDPCDecoders.jl
.
QuantumClifford.ECC.LPCode
— MethodLifted product codes ((Panteleev and Kalachev, 2021), (Panteleev and Kalachev, Jun 2022))
Implemented as a package extension with Hecke. Check the QuantumClifford documentation for more details on that extension.
QuantumClifford.ECC.LiftedCode
— MethodClassical codes lifted over a group algebra, used for lifted product code construction ((Panteleev and Kalachev, 2021), (Panteleev and Kalachev, Jun 2022))
Implemented as a package extension with Hecke. Check the QuantumClifford documentation for more details on that extension.
QuantumClifford.ECC.PyBeliefPropDecoder
— MethodA Belief Propagation decoder built around tools from the python package ldpc
available from the julia package PyQDecoders.jl
.
QuantumClifford.ECC.PyBeliefPropOSDecoder
— MethodA Belief Propagation decoder with ordered statistics decoding, built around tools from the python package ldpc
available from the julia package PyQDecoders.jl
.
QuantumClifford.ECC.PyMatchingDecoder
— MethodA perfect matching decoder built around tools from the python package pymatching
available from the julia package PyQDecoders.jl
.
QuantumClifford.ECC.bicycle_codes
— FunctionImplemented in a package extension with Hecke.
QuantumClifford.ECC.code_k
— MethodThe 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.
QuantumClifford.ECC.code_n
— FunctionThe number of physical qubits in a code.
QuantumClifford.ECC.code_s
— FunctionThe number of stabilizer checks in a code. They might not be all linearly independent, thus code_s >= code_n-code_k
. For the number of linearly independent checks you can use LinearAlgebra.rank
.
QuantumClifford.ECC.distance
— FunctionThe distance of a code.
QuantumClifford.ECC.evaluate_decoder
— MethodEvaluate the performance of a given decoder (e.g. TableDecoder
) and a given style of running an ECC code (e.g. ShorSyndromeECCSetup
)
QuantumClifford.ECC.evaluate_decoder
— MethodEvaluate 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.
QuantumClifford.ECC.faults_matrix
— MethodError-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 mapO[:,1:n]
is the X-physical-error-to-logical-observable mapO[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.
QuantumClifford.ECC.generalized_bicycle_codes
— FunctionImplemented in a package extension with Hecke.
QuantumClifford.ECC.haah_cubic_codes
— FunctionImplemented in a package extension with Hecke.
QuantumClifford.ECC.iscss
— MethodCheck if the code is CSS.
Return nothing
if unknown from the type.
QuantumClifford.ECC.isdegenerate
— FunctionCheck 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
QuantumClifford.ECC.naive_encoding_circuit
— MethodEncoding 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, 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.
QuantumClifford.ECC.naive_syndrome_circuit
— FunctionGenerate 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
QuantumClifford.ECC.parity_checks
— FunctionParity check tableau of a code.
See also: parity_checks_x
and parity_checks_z
QuantumClifford.ECC.parity_checks_x
— MethodParity check boolean matrix of a code (only the X entries in the tableau, i.e. the checks for Z errors).
Only CSS codes have this method.
See also: parity_checks
QuantumClifford.ECC.parity_checks_z
— MethodParity check boolean matrix of a code (only the Z entries in the tableau, i.e. the checks for X errors).
Only CSS codes have this method.
See also: parity_checks
QuantumClifford.ECC.random_all_to_all_circuit_code
— FunctionRandom all-to-all Clifford circuit code (Brown and Fawzi, Jul 2013).
The code of n
qubits is generated by an all-to-all random Clifford circuit of ngates
gates that encodes a subset of qubits encode_qubits
into logical qubits.
Because of the random picking, the size of encode_qubits
is the only thing that matters for the code, referred to as k
.
See also: random_all_to_all_clifford_circuit
, CircuitCode
QuantumClifford.ECC.random_brickwork_circuit_code
— FunctionRandom brickwork Clifford circuit code (Brown and Fawzi, Jul 2013).
The code is generated by a brickwork random Clifford circuit of nlayers
layers that encodes a subset of qubits encode_qubits
into logical qubits.
See also: random_brickwork_clifford_circuit
, CircuitCode
QuantumClifford.ECC.rate
— MethodThe rate of a code.
QuantumClifford.ECC.shor_syndrome_circuit
— FunctionGenerate 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
QuantumClifford.ECC.two_block_group_algebra_codes
— FunctionImplemented in a package extension with Hecke.
Implemented in an extension requiring Hecke.jl
QuantumCliffordHeckeExt.LPCode
— Typestruct LPCode <: QuantumClifford.ECC.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 matrix is obtained by applying repr
to each element of the matrix resulted from the hypergraph product, which is mathematically a linear map from each group algebra element to a binary matrix.
Constructors
Multiple constructors are available:
Two base matrices of group algebra elements.
Two lifted codes, whose base matrices are for quantum code construction.
Two base matrices of group elements, where each group element will be considered as a group algebra element by assigning a unit coefficient.
Two base matrices of integers, where each integer represent the shift of a cyclic permutation. The order of the cyclic permutation should be specified.
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 codetwo_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 codegeneralized_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
.
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 asA
.GA::Hecke.GroupAlgebra
: the group algebra for which elements inA
andB
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.
QuantumCliffordHeckeExt.LiftedCode
— Typestruct LiftedCode <: QuantumClifford.ECC.ClassicalCode
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:
A matrix of group algebra elements.
A matrix of group elements, where a group element will be considered as a group algebra element by assigning a unit coefficient.
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.
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: LPCode
.
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 inA
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.
QuantumCliffordHeckeExt.LiftedCode
— MethodLiftedCode
constructor using the default GF(2)
representation (coefficients converted to a permutation matrix by representation_matrix
provided by Hecke).
QuantumClifford.ECC.bicycle_codes
— MethodBicycle 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.
See also: two_block_group_algebra_codes
, generalized_bicycle_codes
, haah_cubic_codes
.
QuantumClifford.ECC.check_repr_commutation_relation
— MethodChecks the commutation relation between the left and right representation matrices for two randomly-sampled elements a
and b
in the group algebra ℱ[G]
with a general group G
. It verifies the commutation relation that states, L(a)·R(b) = R(b)·L(a)
. This property shows that matrices from the left and right representation sets commute with each other, which is an important property related to the CSS orthogonality condition.
QuantumClifford.ECC.generalized_bicycle_codes
— MethodGeneralized 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
.
See also: two_block_group_algebra_codes
, bicycle_codes
.
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)
QuantumClifford.ECC.haah_cubic_codes
— MethodHaah’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.
The ECC Zoo has an entry for this family.
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)
See also: bicycle_codes
, generalized_bicycle_codes
, two_block_group_algebra_codes
.
QuantumClifford.ECC.two_block_group_algebra_codes
— MethodTwo-block group algebra (2BGA) codes, which are a special case of lifted product codes from two group algebra elements a
and b
, used as 1x1
base matrices.
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) with direct product of C₄ x 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> code_n(c), code_k(c)
(56, 28)
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> code_n(c), code_k(c)
(48, 4)
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> code_n(c), code_k(c)
(108, 12)
See also: LPCode
, generalized_bicycle_codes
, bicycle_codes
, haah_cubic_codes
.
QuantumCliffordHeckeExt.group_algebra_conj
— MethodCompute the conjugate of a group algebra element. The conjugate is defined by inversing elements in the associated group.