- class cudaq_qec.plugin.decoders.tensor_network_decoder.TensorNetworkDecoder
A general class for tensor network decoders for quantum error correction codes.
This decoder constructs a tensor network representation of a quantum code using its parity check matrix, logical observables, and noise model. The tensor network is based on the Tanner graph of the code and can be contracted to compute the probability that a logical observable has flipped, given a syndrome.
The decoder supports both single-syndrome and batch decoding, and can run on CPU or GPU (using cuTensorNet if available).
The Tensor Network Decoder is a Python-only implementation and it requires Python 3.11 or higher. C++ APIs are not available for this decoder.
Note
It is recommended to create decoders using the
cudaq_qec
plugin API:import cudaq_qec as qec import numpy as np # Example: [3,1] repetition code H = np.array([[1, 1, 0], [0, 1, 1]], dtype=np.uint8) logical_obs = np.array([[1, 1, 1]], dtype=np.uint8) noise_model = [0.1, 0.1, 0.1] decoder = qec.get_decoder("tensor_network_decoder", H, logical_obs=logical_obs, noise_model=noise_model) syndrome = [0.0, 1.0] result = decoder.decode(syndrome)
Tensor Network Structure
The tensor network constructed by this decoder is based on the Tanner graph of the code, extended with noise and logical observable tensors. The structure is illustrated below:
open/output index < logical observable -------- | s1 s2 | s3 < syndromes : product of 2D vectors [1 , 1-2pi] (pi is the probability detector i flipped) | | | | ----| c1 c2 l1 c3 < checks / logical | : delta tensors | / | | \ | | H H H H H H < Hadamard matrices | TANNER (bipartite) GRAPH \ | | / | / | e1 e2 e3 < errors | : delta tensors | | / -----| \ / / P(e1, e2, e3) < noise / error model : classical probability density ci, ej, lk are delta tensors represented sparsely as indices.
- Parameters:
H – Parity check matrix (numpy.ndarray), shape (num_checks, num_qubits)
logical_obs – Logical observable matrix (numpy.ndarray), shape (1, num_qubits)
noise_model – Noise model, either a list of probabilities (length = num_qubits) or a quimb.tensor.TensorNetwork
check_inds – (optional) List of check index names
error_inds – (optional) List of error index names
logical_inds – (optional) List of logical index names
logical_tags – (optional) List of logical tags
contract_noise_model – (bool, optional) Whether to contract the noise model at initialization (default: True)
dtype – (str, optional) Data type for tensors (default: “float32”)
device – (str, optional) Device for tensor operations (“cpu”, “cuda”, or “cuda:X”, default: “cuda”)
Methods
- decode(syndrome)
Decode a single syndrome by contracting the tensor network.
- Parameters:
syndrome – List of float values (soft-decision probabilities) for each check.
- Returns:
DecoderResult with the probability that the logical observable flipped.
- decode_batch(syndrome_batch)
Decode a batch of syndromes.
- Parameters:
syndrome_batch – numpy.ndarray of shape (batch_size, num_checks)
- Returns:
List of DecoderResult objects with the probability that the logical observable has flipped for each syndrome.
- optimize_path(optimize=None, batch_size=-1)
Optimize the contraction path for the tensor network.
- Parameters:
optimize – Optimization options or None
batch_size – (int, optional) Batch size for optimization (default: -1, no batching)
- Returns:
Optimizer info object