Skip to content

Simulator

Simulation and sampling manager for color code circuits.

This class handles Monte Carlo simulation, detector/observable sampling, and utility functions for error analysis. It operates on pre-built circuits and integrates with the decoder architecture through dependency injection.

Source code in src/color_code_stim/simulation/simulator.py
class Simulator:
    """
    Simulation and sampling manager for color code circuits.

    This class handles Monte Carlo simulation, detector/observable sampling,
    and utility functions for error analysis. It operates on pre-built circuits
    and integrates with the decoder architecture through dependency injection.
    """

    def __init__(
        self,
        *,
        circuit: stim.Circuit,
        circuit_type: str,
    ):
        """
        Initialize the Simulator with required dependencies.

        Parameters
        ----------
        circuit : stim.Circuit
            The quantum circuit to simulate
        circuit_type : str
            Type of circuit ('tri', 'rec', 'rec_stability', 'growing', 'cult+growing')
        """
        self.circuit = circuit
        self.circuit_type = circuit_type

    def sample(
        self, shots: int, seed: Optional[int] = None
    ) -> Tuple[np.ndarray, np.ndarray]:
        """
        Sample detector outcomes and observables from the quantum circuit.

        Parameters
        ----------
        shots : int
            Number of samples to generate
        seed : int, optional
            Seed value to initialize the random number generator

        Returns
        -------
        det : 2D numpy array of bool
            Detector outcomes. det[i,j] is True if and only if the detector
            with id j in the i-th sample has an outcome of −1.
        obs : 1D or 2D numpy array of bool
            Observable outcomes. If there is only one observable, returns a 1D array;
            otherwise returns a 2D array. obs[i] or obs[i,j] is True if and only if
            the j-th observable (j=0 when 1D) of the i-th sample has an outcome of -1.
        """
        sampler = self.circuit.compile_detector_sampler(seed=seed)
        det, obs = sampler.sample(shots, separate_observables=True)
        if obs.shape[1] == 1:
            obs = obs.ravel()
        return det, obs

    def sample_with_errors(
        self,
        shots: int,
        seed: Optional[int] = None,
    ) -> Tuple[np.ndarray, np.ndarray, np.ndarray]:
        """
        Sample detector outcomes, observables, and errors from the quantum circuit.

        Parameters
        ----------
        shots : int
            Number of samples to generate
        seed : int, optional
            Seed value to initialize the random number generator

        Returns
        -------
        det : 2D numpy array of bool
            Detector outcomes. det[i,j] is True if and only if the detector
            with id j in the i-th sample has an outcome of −1.
        obs : 1D or 2D numpy array of bool
            Observable outcomes. If there is only one observable, returns a 1D array;
            otherwise returns a 2D array. obs[i] or obs[i,j] is True if and only if
            the j-th observable (j=0 when 1D) of the i-th sample has an outcome of -1.
        errors : 2D numpy array of bool
            Errors sampled from the quantum circuit. errors[i,j] is True if and only if
            the j-th error (in the DEM) of the i-th sample has an outcome of -1.
        """
        dem = self.circuit.detector_error_model()
        sampler = dem.compile_sampler(seed=seed)
        det, obs, err = sampler.sample(shots, return_errors=True)
        if obs.shape[1] == 1:
            obs = obs.ravel()

        return det, obs, err

    def simulate(
        self,
        shots: int,
        *,
        decoder_func,
        colors: Union[List[str], str] = "all",
        alpha: float = 0.01,
        confint_method: str = "wilson",
        full_output: bool = False,
        seed: Optional[int] = None,
        verbose: bool = False,
        **kwargs,
    ) -> Tuple[np.ndarray, dict]:
        """
        Monte-Carlo simulation of quantum error correction decoding.

        Parameters
        ----------
        shots : int
            Number of shots to simulate.
        decoder_func : callable
            Decoder function that takes detector outcomes and returns predictions.
            Should have signature: decoder_func(detector_outcomes, **kwargs) -> predictions
        colors : Union[List[str], str], default 'all'
            Colors of the sub-decoding procedures to consider. Can be 'all', one of {'r', 'g', 'b'},
            or a list containing any combination of {'r', 'g', 'b'}.
        alpha : float, default 0.01
            Significance level for the confidence interval calculation.
        confint_method : str, default 'wilson'
            Method to calculate the confidence interval.
            See statsmodels.stats.proportion.proportion_confint for available options.
        full_output: bool = False,
            If True, return additional information.
        seed : Optional[int], default None
            Seed to initialize the random number generator.
        verbose : bool, default False
            If True, print progress information during simulation.
        **kwargs :
            Additional keyword arguments for the decoder function.

        Returns
        -------
        num_fails : numpy.ndarray
            Number of failures for each observable.
        extra_outputs : dict
            Dictionary containing additional information:
            - 'stats': Tuple of (pfail, delta_pfail) where pfail is the estimated failure rate
              and delta_pfail is the half-width of the confidence interval
            - 'fails': Boolean array indicating which samples failed
            - Additional outputs from the decoder function if full_output=True
        """
        if self.circuit_type == "cult+growing":
            raise NotImplementedError(
                "Cult+growing circuit type is not supported for this method."
            )

        if colors == "all":
            colors = ["r", "g", "b"]

        if verbose:
            print("Sampling...")

        shots = round(shots)
        det, obs = self.sample(shots, seed=seed)

        if verbose:
            print("Decoding...")

        preds = decoder_func(
            det,
            colors=colors,
            verbose=verbose,
            full_output=full_output,
            **kwargs,
        )

        if full_output:
            preds, extra_outputs = preds

        if verbose:
            print("Postprocessing...")

        fails = np.logical_xor(obs, preds)
        num_fails = np.sum(fails, axis=0)

        if full_output:
            pfail, delta_pfail = get_pfail(
                shots, num_fails, alpha=alpha, confint_method=confint_method
            )
            extra_outputs["stats"] = (pfail, delta_pfail)
            extra_outputs["fails"] = fails

            return num_fails, extra_outputs

        else:
            return num_fails

__init__(*, circuit, circuit_type)

Initialize the Simulator with required dependencies.

Parameters:

Name Type Description Default
circuit Circuit

The quantum circuit to simulate

required
circuit_type str

Type of circuit ('tri', 'rec', 'rec_stability', 'growing', 'cult+growing')

required
Source code in src/color_code_stim/simulation/simulator.py
def __init__(
    self,
    *,
    circuit: stim.Circuit,
    circuit_type: str,
):
    """
    Initialize the Simulator with required dependencies.

    Parameters
    ----------
    circuit : stim.Circuit
        The quantum circuit to simulate
    circuit_type : str
        Type of circuit ('tri', 'rec', 'rec_stability', 'growing', 'cult+growing')
    """
    self.circuit = circuit
    self.circuit_type = circuit_type

sample(shots, seed=None)

Sample detector outcomes and observables from the quantum circuit.

Parameters:

Name Type Description Default
shots int

Number of samples to generate

required
seed int

Seed value to initialize the random number generator

None

Returns:

Name Type Description
det 2D numpy array of bool

Detector outcomes. det[i,j] is True if and only if the detector with id j in the i-th sample has an outcome of −1.

obs 1D or 2D numpy array of bool

Observable outcomes. If there is only one observable, returns a 1D array; otherwise returns a 2D array. obs[i] or obs[i,j] is True if and only if the j-th observable (j=0 when 1D) of the i-th sample has an outcome of -1.

Source code in src/color_code_stim/simulation/simulator.py
def sample(
    self, shots: int, seed: Optional[int] = None
) -> Tuple[np.ndarray, np.ndarray]:
    """
    Sample detector outcomes and observables from the quantum circuit.

    Parameters
    ----------
    shots : int
        Number of samples to generate
    seed : int, optional
        Seed value to initialize the random number generator

    Returns
    -------
    det : 2D numpy array of bool
        Detector outcomes. det[i,j] is True if and only if the detector
        with id j in the i-th sample has an outcome of −1.
    obs : 1D or 2D numpy array of bool
        Observable outcomes. If there is only one observable, returns a 1D array;
        otherwise returns a 2D array. obs[i] or obs[i,j] is True if and only if
        the j-th observable (j=0 when 1D) of the i-th sample has an outcome of -1.
    """
    sampler = self.circuit.compile_detector_sampler(seed=seed)
    det, obs = sampler.sample(shots, separate_observables=True)
    if obs.shape[1] == 1:
        obs = obs.ravel()
    return det, obs

sample_with_errors(shots, seed=None)

Sample detector outcomes, observables, and errors from the quantum circuit.

Parameters:

Name Type Description Default
shots int

Number of samples to generate

required
seed int

Seed value to initialize the random number generator

None

Returns:

Name Type Description
det 2D numpy array of bool

Detector outcomes. det[i,j] is True if and only if the detector with id j in the i-th sample has an outcome of −1.

obs 1D or 2D numpy array of bool

Observable outcomes. If there is only one observable, returns a 1D array; otherwise returns a 2D array. obs[i] or obs[i,j] is True if and only if the j-th observable (j=0 when 1D) of the i-th sample has an outcome of -1.

errors 2D numpy array of bool

Errors sampled from the quantum circuit. errors[i,j] is True if and only if the j-th error (in the DEM) of the i-th sample has an outcome of -1.

Source code in src/color_code_stim/simulation/simulator.py
def sample_with_errors(
    self,
    shots: int,
    seed: Optional[int] = None,
) -> Tuple[np.ndarray, np.ndarray, np.ndarray]:
    """
    Sample detector outcomes, observables, and errors from the quantum circuit.

    Parameters
    ----------
    shots : int
        Number of samples to generate
    seed : int, optional
        Seed value to initialize the random number generator

    Returns
    -------
    det : 2D numpy array of bool
        Detector outcomes. det[i,j] is True if and only if the detector
        with id j in the i-th sample has an outcome of −1.
    obs : 1D or 2D numpy array of bool
        Observable outcomes. If there is only one observable, returns a 1D array;
        otherwise returns a 2D array. obs[i] or obs[i,j] is True if and only if
        the j-th observable (j=0 when 1D) of the i-th sample has an outcome of -1.
    errors : 2D numpy array of bool
        Errors sampled from the quantum circuit. errors[i,j] is True if and only if
        the j-th error (in the DEM) of the i-th sample has an outcome of -1.
    """
    dem = self.circuit.detector_error_model()
    sampler = dem.compile_sampler(seed=seed)
    det, obs, err = sampler.sample(shots, return_errors=True)
    if obs.shape[1] == 1:
        obs = obs.ravel()

    return det, obs, err

simulate(shots, *, decoder_func, colors='all', alpha=0.01, confint_method='wilson', full_output=False, seed=None, verbose=False, **kwargs)

Monte-Carlo simulation of quantum error correction decoding.

Parameters:

Name Type Description Default
shots int

Number of shots to simulate.

required
decoder_func callable

Decoder function that takes detector outcomes and returns predictions. Should have signature: decoder_func(detector_outcomes, **kwargs) -> predictions

required
colors Union[List[str], str]

Colors of the sub-decoding procedures to consider. Can be 'all', one of {'r', 'g', 'b'}, or a list containing any combination of {'r', 'g', 'b'}.

'all'
alpha float

Significance level for the confidence interval calculation.

0.01
confint_method str

Method to calculate the confidence interval. See statsmodels.stats.proportion.proportion_confint for available options.

'wilson'
full_output bool

If True, return additional information.

False
seed Optional[int]

Seed to initialize the random number generator.

None
verbose bool

If True, print progress information during simulation.

False
**kwargs

Additional keyword arguments for the decoder function.

{}

Returns:

Name Type Description
num_fails ndarray

Number of failures for each observable.

extra_outputs dict

Dictionary containing additional information: - 'stats': Tuple of (pfail, delta_pfail) where pfail is the estimated failure rate and delta_pfail is the half-width of the confidence interval - 'fails': Boolean array indicating which samples failed - Additional outputs from the decoder function if full_output=True

Source code in src/color_code_stim/simulation/simulator.py
def simulate(
    self,
    shots: int,
    *,
    decoder_func,
    colors: Union[List[str], str] = "all",
    alpha: float = 0.01,
    confint_method: str = "wilson",
    full_output: bool = False,
    seed: Optional[int] = None,
    verbose: bool = False,
    **kwargs,
) -> Tuple[np.ndarray, dict]:
    """
    Monte-Carlo simulation of quantum error correction decoding.

    Parameters
    ----------
    shots : int
        Number of shots to simulate.
    decoder_func : callable
        Decoder function that takes detector outcomes and returns predictions.
        Should have signature: decoder_func(detector_outcomes, **kwargs) -> predictions
    colors : Union[List[str], str], default 'all'
        Colors of the sub-decoding procedures to consider. Can be 'all', one of {'r', 'g', 'b'},
        or a list containing any combination of {'r', 'g', 'b'}.
    alpha : float, default 0.01
        Significance level for the confidence interval calculation.
    confint_method : str, default 'wilson'
        Method to calculate the confidence interval.
        See statsmodels.stats.proportion.proportion_confint for available options.
    full_output: bool = False,
        If True, return additional information.
    seed : Optional[int], default None
        Seed to initialize the random number generator.
    verbose : bool, default False
        If True, print progress information during simulation.
    **kwargs :
        Additional keyword arguments for the decoder function.

    Returns
    -------
    num_fails : numpy.ndarray
        Number of failures for each observable.
    extra_outputs : dict
        Dictionary containing additional information:
        - 'stats': Tuple of (pfail, delta_pfail) where pfail is the estimated failure rate
          and delta_pfail is the half-width of the confidence interval
        - 'fails': Boolean array indicating which samples failed
        - Additional outputs from the decoder function if full_output=True
    """
    if self.circuit_type == "cult+growing":
        raise NotImplementedError(
            "Cult+growing circuit type is not supported for this method."
        )

    if colors == "all":
        colors = ["r", "g", "b"]

    if verbose:
        print("Sampling...")

    shots = round(shots)
    det, obs = self.sample(shots, seed=seed)

    if verbose:
        print("Decoding...")

    preds = decoder_func(
        det,
        colors=colors,
        verbose=verbose,
        full_output=full_output,
        **kwargs,
    )

    if full_output:
        preds, extra_outputs = preds

    if verbose:
        print("Postprocessing...")

    fails = np.logical_xor(obs, preds)
    num_fails = np.sum(fails, axis=0)

    if full_output:
        pfail, delta_pfail = get_pfail(
            shots, num_fails, alpha=alpha, confint_method=confint_method
        )
        extra_outputs["stats"] = (pfail, delta_pfail)
        extra_outputs["fails"] = fails

        return num_fails, extra_outputs

    else:
        return num_fails