API Reference

This documentation is divided into two main parts: a High-Level API for general usage and an Advanced & Technical API for developers interested in the library’s internal workings.

High-Level API

This section covers the primary classes needed to build, configure, and run a simulation.

Simulation Control

The main user-facing classes for constructing and simulating a quantum circuit.

template<typename Coefficient_t = coeff_t>
class Circuit

Represents a quantum circuit and provides a high-level simulation interface.

The Circuit class is the main user-facing interface for building and simulating quantum circuits. It allows users to add a sequence of quantum operations and then run a simulation to compute the resulting observable and its expectation value.

This class manages the complexities of the simulation, including the application of gates, noise, truncation, and merging, based on the provided policies. The simulation is performed in the Heisenberg picture, where the observable is evolved backward through the circuit.

Template Parameters:

Coefficient_t – The numeric type for coefficients (e.g., float, double). Must be a floating-point type.

Public Types

using O_t = Observable<Coefficient_t>

Type alias for the Observable used by this circuit instance.

Public Functions

template<typename TruncatorPtr = std::shared_ptr<Truncator<Coefficient_t>>>
inline Circuit(unsigned nb_qubits, TruncatorPtr truncator = std::make_shared<NeverTruncator<Coefficient_t>>(), NoiseModel<Coefficient_t> const &noise_model = {}, std::shared_ptr<SchedulingPolicy> merge_policy = std::make_shared<AlwaysAfterSplittingPolicy>(), std::shared_ptr<SchedulingPolicy> truncate_policy = std::make_shared<AlwaysAfterSplittingPolicy>())

Constructs a new Circuit.

    Circuit qc{ 2 };
    qc.add_operation("H", 0);
    qc.add_operation("Rz", 0, 1.57f);
    qc.add_operation("CX", 0u, 1u);

    auto result = qc.run(Observable{ "ZZ" });
    std::cout << "Expectation value: " << result.expectation_value() << std::endl;
    Circuit qc{ 64,
            combine_truncators(CoefficientTruncator<>{ 0.001f }, // remove terms with coefficient below 0.001
                       WeightTruncator<>{ 6 } // remove terms with pauli weight > 6
                       ) };

    // Apply a layer of Hadamard gates
    for (unsigned i = 0; i < 64; ++i)
        qc.add_operation("H", i);

    // Entangling layer
    for (unsigned i = 0; i < 63; ++i) {
        qc.add_operation("CX", i, i + 1);
    }

    auto result = qc.run(Observable{ std::string(64, 'Z') });
    std::cout << "Expectation value: " << result.expectation_value() << std::endl;
    // Policies are passed to the Circuit constructor to control optimization.
    // This circuit will merge identical Pauli terms only after a splitting gate is applied.
    Circuit qc(4, std::make_shared<NeverTruncator<>>(), {},
           std::make_shared<AlwaysAfterSplittingPolicy>(), // Merge policy
           std::make_shared<NeverPolicy>() // Truncate policy
    );

Template Parameters:

TruncatorPtr – A type convertible to a shared pointer to a Truncator. This is typically auto-deduced.

Parameters:
  • nb_qubits – The number of qubits in the circuit.

  • truncator – A shared pointer to a truncator object that defines how to simplify the observable. The default NeverTruncator performs no truncation, which is only suitable for small circuits.

  • noise_model – A noise model to apply to the gates. By default, an ideal (noiseless) model is used.

  • merge_policy – A policy defining when to merge identical Pauli terms. The default policy merges after any gate that can increase the number of terms.

  • truncate_policy – A policy defining when to apply the truncator. The default policy truncates after any gate that can increase the number of terms.

inline unsigned nb_qubits() const

Gets the number of qubits in the circuit.

Returns:

The number of qubits.

template<typename ...T>
inline void add_operation(std::string op, T&&... args)

Adds a quantum operation to the circuit by name.

See also

add_operation(QGate, T&&…)

Note

This is a convenience function that looks up the QGate enum from a string and forwards to the primary add_operation overload. It is intended for ease of use in user-facing code.

Template Parameters:

T – Parameter pack for the arguments of the operation.

Parameters:
  • op – The name of the operation (e.g., “H”, “CX”, “Rz”). The name is case-insensitive.

  • args – The arguments for the operation. For single-qubit gates, this is (unsigned qubit). For two-qubit gates, (unsigned control, unsigned target). For parameterized gates, (unsigned qubit, Real angle), etc.

Pre:

Qubit indices provided in args must be less than nb_qubits().

template<typename ...T>
inline void add_operation(QGate g, T&&... args)

Adds a quantum operation to the circuit using its enum type.

This is the primary method for adding operations to the circuit. It dispatches to the correct internal implementation based on the gate type and arguments, and applies any noise specified by the NoiseModel after the ideal operation.

Template Parameters:

T – Parameter pack for the arguments of the operation.

Parameters:
  • g – The QGate enum representing the operation.

  • args – The arguments for the operation, such as qubit indices or rotation angles.

Pre:

Qubit indices provided in args must be less than nb_qubits().

template<typename ExecutionPolicy = DefaultExecutionPolicy>
inline Observable<Coefficient_t> run(Observable<Coefficient_t> const &target_observable, ExecutionPolicy &&policy = ExecutionPolicy{})

Runs the simulation on the circuit.

The run method executes the simulation by applying each gate in the circuit to the observable in reverse order (Heisenberg picture). It returns the final observable, from which the expectation value can be calculated. If the observable becomes empty (representing the maximally mixed state, proportional to identity), the simulation may terminate early.

See also

Observable::expectation_value()

    Circuit qc{ 2 };
    qc.add_operation("H", 0);
    qc.add_operation("Rz", 0, 1.57f);
    qc.add_operation("CX", 0u, 1u);

    auto result = qc.run(Observable{ "ZZ" });
    std::cout << "Expectation value: " << result.expectation_value() << std::endl;

Parameters:

target_observable – The initial observable to be propagated backward through the circuit.

Pre:

The number of qubits in target_observable must match nb_qubits().

Returns:

The final, evolved observable after applying all circuit operations in reverse.

inline std::size_t nb_splitting_gates() const

Counts the number of gates in the circuit that can split an observable.

Splitting gates are operations that can increase the number of Pauli terms in the observable, making the simulation more complex. This count is used to initialize the SimulationState.

Returns:

The total number of splitting gates (e.g., Rz, AmplitudeDamping).

inline void reset()

Clears all operations from the circuit.

Resets the circuit to an empty state, ready to be populated with new operations. The number of qubits and configured policies are not affected.

inline void set_truncator(std::shared_ptr<Truncator<Coefficient_t>> truncator)

Sets a new truncator for the circuit.

See also

Truncator

Parameters:

truncator – A shared pointer to the new truncator.

inline void set_merge_policy(std::shared_ptr<SchedulingPolicy> policy)

Sets a new policy for when to merge Pauli terms.

See also

SchedulingPolicy

Parameters:

policy – A shared pointer to the new merge policy.

inline void set_truncate_policy(std::shared_ptr<SchedulingPolicy> policy)

Sets a new policy for when to truncate the observable.

See also

SchedulingPolicy

Parameters:

policy – A shared pointer to the new truncate policy.

template<typename T>
class NoiseModel

A model for applying noise to quantum gates.

A NoiseModel allows you to specify that certain types of noise should be automatically applied after certain gates in a Circuit.

Template Parameters:

T – The numeric type for noise parameters.

Public Functions

inline NoiseModel()

Constructs an empty noise model.

    // Create a noise model
    NoiseModel<float> nm;

    // Add 1% amplitude damping noise to every CX gate
    nm.add_amplitude_damping_on_gate(QGate::Cx, 0.01f);

    // Add 0.5% dephasing noise to every Hadamard gate
    nm.add_unital_noise_on_gate(QGate::H, UnitalNoise::Dephasing, 0.005f);
    // Add 1% depolarizing noise to every Cx gate
    nm.add_unital_noise_on_gate(QGate::Cx, UnitalNoise::Depolarizing, 0.01f);

    // Create a circuit with this noise model
    Circuit qc(4, std::make_shared<NeverTruncator<>>(), nm);

    qc.add_operation("H", 0);
    qc.add_operation("CX", 0, 1); // Noise will be applied automatically

    auto result = qc.run(Observable("ZZZZ"));
inline void add_unital_noise_on_gate(QGate g, UnitalNoise un, T n)

Adds a unital noise channel to be applied after a specific gate type.

Warning

Adding noise to noisy gates (AmplitudeDamping, Depolarizing, Dephasing) is forbidden to prevent infinite loops.

Parameters:
  • g – The gate type after which the noise should be applied.

  • un – The type of unital noise (Depolarizing or Dephasing).

  • n – The strength/probability of the noise.

inline void add_amplitude_damping_on_gate(QGate g, T n)

Adds an amplitude damping channel to be applied after a specific gate type.

Warning

Adding noise to noisy gates (AmplitudeDamping, Depolarizing, Dephasing) is forbidden to prevent infinite loops.

Parameters:
  • g – The gate type after which the noise should be applied.

  • n – The strength/probability of the noise.


Primary Data Structures

The fundamental classes for representing the quantum state and its components.

template<typename T = coeff_t>
class Observable

Represents a quantum observable as a linear combination of Pauli strings.

An Observable is a sum of Pauli terms. This class is the central data structure in the Pauli back-propagation simulation. It is implemented using a memory-efficient PauliTermContainer to store its terms. Quantum gates are applied to this object, evolving it backward in the Heisenberg picture.

Template Parameters:

T – The numeric type for the coefficients (e.g., float, double).

Constructors

Ways to construct an Observable.

inline Observable(std::string_view pauli_string, typename std::enable_if_t<std::is_constructible_v<T, coeff_t>, T> coeff = T{1})

Constructs an observable from a single Pauli string.

Parameters:
  • pauli_string – A string representing the Pauli operators (e.g., “IXYZ”).

  • coeff – The coefficient of this Pauli term.

        // From a single string_view and an explicit coefficient
        auto obs_from_string = Observable<float>("ZIZ", -1.0f);
        std::cout << "Observable from string: " << obs_from_string << std::endl;
    

inline Observable(std::initializer_list<std::string_view> pauli_string_list)

Constructs an observable from a list of Pauli strings.

Parameters:

pauli_string_list – An initializer list of Pauli strings. Each will have a coefficient of 1.

    // From an initializer_list of strings (coefficients default to 1.0)
    auto obs_from_string_list = Observable<double>({"XXI", "IYY", "ZZZ"});
    std::cout << "Observable from string list: " << obs_from_string_list << std::endl;

inline Observable(std::initializer_list<PauliTerm<T>> paulis_list)

Constructs an observable from a list of PauliTerm objects.

Parameters:

paulis_list – An initializer list of PauliTerm objects.

    // From an initializer_list of full PauliTerm objects
    using PauliTermF = PauliTerm<float>;
    auto obs_from_pauli_terms = Observable<float>({
        PauliTermF("X", 0.5f),
        PauliTermF("Y", -0.5f)
    });
    std::cout << "Observable from PauliTerm list: " << obs_from_pauli_terms << std::endl;

template<PauliTermIterator Iter>
inline Observable(Iter &&begin, Iter &&end)

Constructs an observable from a range of PauliTerm-like objects.

Template Parameters:

Iter – An iterator type.

Parameters:
  • begin – An iterator to the beginning of the range.

  • end – An iterator to the end of the range.

        // From a range specified by iterators (e.g., from a std::vector)
        std::vector<PauliTermF> terms = {PauliTermF("IX", 1.0), PauliTermF("YI", 1.0)};
        auto obs_from_iterators = Observable<float>(terms.begin(), terms.end());
        std::cout << "Observable from iterators: " << obs_from_iterators << std::endl;
    

Gate Application

Methods for applying quantum gates to the observable.

template<IsNotVariant ExecutionPolicy = DefaultExecutionPolicy>
inline void apply_pauli(Pauli_gates g, unsigned qubit, ExecutionPolicy &&policy = ExecutionPolicy{})

Applies a single-qubit Pauli gate to the observable.

Parameters:
  • g – The Pauli gate to apply.

  • qubit – The index of the qubit to apply the gate to.

Pre:

qubit must be a valid index less than nb_qubits().

template<IsNotVariant ExecutionPolicy = DefaultExecutionPolicy>
inline void apply_clifford(Clifford_Gates_1Q g, unsigned qubit, ExecutionPolicy &&policy = ExecutionPolicy{})

Applies a single-qubit Clifford gate to the observable.

Parameters:
  • g – The Clifford gate to apply (e.g., Hadamard).

  • qubit – The index of the qubit to apply the gate to.

Pre:

qubit must be a valid index less than nb_qubits().

template<IsNotVariant ExecutionPolicy = DefaultExecutionPolicy>
inline void apply_unital_noise(UnitalNoise n, unsigned qubit, T p, ExecutionPolicy &&policy = ExecutionPolicy{})

Applies a single-qubit unital noise channel to the observable.

Parameters:
  • n – The type of unital noise (e.g., Depolarizing, Dephasing).

  • qubit – The index of the qubit to apply the noise to.

  • p – The noise probability parameter.

Pre:

qubit must be a valid index less than nb_qubits().

template<IsNotVariant ExecutionPolicy = DefaultExecutionPolicy>
inline void apply_cx(unsigned qubit_control, unsigned qubit_target, ExecutionPolicy &&policy = ExecutionPolicy{})

Applies a CNOT (CX) gate to the observable.

Parameters:
  • qubit_control – The index of the control qubit.

  • qubit_target – The index of the target qubit.

Pre:

qubit_control and qubit_target must be valid and distinct qubit indices.

template<IsNotVariant ExecutionPolicy = DefaultExecutionPolicy>
inline void apply_rz(unsigned qubit, T theta, ExecutionPolicy &&policy = ExecutionPolicy{})

Applies a single-qubit Rz rotation gate to the observable.

Note

This is a splitting operation. Applying an Rz gate to a Pauli term containing an X or Y operator on the target qubit will result in two output Pauli terms. This can increase the size of the observable, often necessitating a subsequent merge() or truncate() call.

Parameters:
  • qubit – The index of the qubit to apply the rotation to.

  • theta – The rotation angle in radians.

template<IsNotVariant ExecutionPolicy = DefaultExecutionPolicy>
inline void apply_amplitude_damping(unsigned qubit, T pn, ExecutionPolicy &&policy = ExecutionPolicy{})

Applies an amplitude damping noise channel.

Note

This can be a splitting operation. If a term has a Z operator on the target qubit, it will be split into two. If it has X or Y, its coefficient is simply scaled. If it has I, there is no effect.

Parameters:
  • qubit – The index of the qubit to apply the channel to.

  • pn – The noise probability parameter.

Container Interface

Methods to interact with the Observable like a container.

inline auto operator[](std::size_t idx)

Accesses a specific Pauli term via a non-owning view.

Returns:

A lightweight view of the term. Modifying it modifies the observable directly.

inline auto operator[](std::size_t idx) const

Accesses a specific Pauli term via a read-only non-owning view.

inline decltype(auto) begin()

Returns an iterator to the beginning of the terms.

inline decltype(auto) begin() const

Returns a const iterator to the beginning of the terms.

inline decltype(auto) cbegin() const

Returns a const iterator to the beginning of the terms.

inline decltype(auto) end()

Returns an iterator to the end of the terms.

inline decltype(auto) end() const

Returns a const iterator to the end of the terms.

inline decltype(auto) cend() const

Returns a const iterator to the end of the terms.

inline std::size_t size() const

Gets the number of Pauli terms in the observable.

inline std::size_t nb_qubits() const

Gets the number of qubits in the observable.

inline PauliTerm<T> copy_term(std::size_t idx) const

Creates an owning copy of a term at a specific index.

Note

This is a potentially expensive operation as it involves constructing a new object.

Parameters:

idx – The index of the term to copy.

Returns:

An owning PauliTerm<T> object.

Public Functions

template<IsNotVariant ExecutionPolicy = DefaultExecutionPolicy>
inline T expectation_value(ExecutionPolicy &&policy = ExecutionPolicy{}) const

Calculates the expectation value of the observable.

The expectation value is the sum of coefficients of all Pauli terms composed entirely of I and Z operators. These are the terms that are diagonal in the computational basis.

Returns:

The expectation value.

template<IsNotVariant ExecutionPolicy = DefaultExecutionPolicy>
inline std::size_t merge(ExecutionPolicy &&policy = ExecutionPolicy{})

Merges terms with identical Pauli strings by summing their coefficients.

This is a crucial optimization for reducing the complexity of the simulation. It calls a high-performance, in-place merging algorithm.

Returns:

The number of terms remaining after the merge.

template<typename TruncatorImpl>
inline std::size_t truncate(TruncatorImpl &&truncator)

Truncates the observable based on a given truncation strategy.

See also

Truncator

Template Parameters:

TruncatorImpl – The type of the truncator object, which must satisfy the Truncator interface.

Parameters:

truncator – The truncator object that defines the truncation logic.

Returns:

The number of Pauli terms removed.

template<typename T>
std::ostream &operator<<(std::ostream &os, Observable<T> const &obs)

Prints the observable to an output stream.

Note

This can be an expensive operation for large observables.

template<class T = coeff_t>
class PauliTerm

Represents a single term in an observable, owning its Pauli string and coefficient.

This class stores a Pauli string (e.g., “IXYZ”) and a corresponding numerical coefficient. It provides methods for applying quantum gates, which modify the term in-place or, for certain gates, split it into two terms.

Template Parameters:

T – The numeric type for the coefficient.

Constructors

inline PauliTerm(std::string_view pauli_string, typename std::enable_if_t<std::is_convertible_v<T, coeff_t>, T> coefficient = T{1.f})

Constructs a PauliTerm from a string representation.

Parameters:
  • pauli_string – A string like “IXYZ” representing the Pauli operators.

  • coefficient – The coefficient of the term.

        // From a string_view and a coefficient
        auto pt_from_string = PauliTerm("IXYZ", 0.5f);
        std::cout << "From string: " << pt_from_string << std::endl;
    

inline PauliTerm(std::initializer_list<Pauli> pauli_list, T coefficient = T{1.})

Constructs a PauliTerm from an initializer list of Pauli objects.

Parameters:
  • pauli_list – An initializer list, e.g., {p_i, p_x, p_y, p_z}.

  • coefficient – The coefficient of the term.

        // From an initializer_list of Pauli objects
        using enum Pauli_enum;
        auto pt_from_list = PauliTerm({ I, X, Y, Z }, -1.0f);
        std::cout << "From list: " << pt_from_list << std::endl;
    

inline PauliTerm(Pauli pauli, T coefficient = T{1})

Constructs a single-qubit PauliTerm from a Pauli object.

Parameters:
  • pauli – The Pauli operator.

  • coefficient – The coefficient of the term.

        // From a single Pauli object (for a 1-qubit term)
        auto pt_from_pauli = PauliTerm(p_x, 1.0f);
        std::cout << "From Pauli object: " << pt_from_pauli << std::endl;
    

Conversion Operators

inline explicit operator NonOwningPauliTerm<T>()

Explicitly converts to a mutable non-owning view.

inline explicit operator ReadOnlyNonOwningPauliTerm<T>() const

Explicitly converts to a read-only non-owning view.

Gate Application

inline void apply_pauli(Pauli_gates g, unsigned qubit)

Applies a Pauli gate to a specific qubit of the term.

Parameters:
  • g – The Pauli gate to apply.

  • qubit – The target qubit index.

inline void apply_clifford(Clifford_Gates_1Q g, unsigned qubit)

Applies a Clifford gate to a specific qubit of the term.

Parameters:
  • g – The single-qubit Clifford gate to apply.

  • qubit – The target qubit index.

inline void apply_unital_noise(UnitalNoise n, unsigned qubit, T p)

Applies a unital noise channel to a specific qubit of the term.

Parameters:
  • n – The noise channel type.

  • qubit – The target qubit index.

  • p – The noise probability.

inline void apply_cx(unsigned control, unsigned target)

Applies a CNOT gate to the term.

Parameters:
  • control – The control qubit index.

  • target – The target qubit index.

Pre:

Control and target qubits must be different.

inline PauliTerm<T> apply_rz(unsigned qubit, T theta)

Applies an Rz gate, potentially splitting the term.

Parameters:
  • qubit – The target qubit index.

  • theta – The rotation angle.

Returns:

A new PauliTerm representing the second term created by the split. The original term is modified to become the first term of the split.

Pre:

The operator on the target qubit must be X or Y.

inline void apply_amplitude_damping_xy(unsigned qubit, T p)

Applies the X/Y part of the amplitude damping channel, scaling the coefficient.

Parameters:
  • qubit – The target qubit index.

  • p – The damping probability.

Pre:

The operator on the target qubit must be X or Y.

inline PauliTerm<T> apply_amplitude_damping_z(unsigned qubit, T p)

Applies the Z part of the amplitude damping channel, splitting the term.

Parameters:
  • qubit – The target qubit index.

  • p – The damping probability.

Returns:

The new PauliTerm created by the split. The original term is also modified.

Pre:

The operator on the target qubit must be Z.

Comparison

inline bool equal_bitstring(PauliTerm const &oth) const

Checks if two terms have the same Pauli string, ignoring the coefficient.

Parameters:

oth – The other term to compare with.

Returns:

true if the Pauli strings are identical.

inline friend bool operator==(PauliTerm const &lhs, PauliTerm const &rhs)

Checks for deep equality between two PauliTerm objects.

Returns:

true if both the Pauli strings and coefficients are identical.

Coefficient Manipulation

inline void add_coeff(T c)

Adds a value to the term’s coefficient.

Warning

The resulting coefficient should remain within a valid range.

Parameters:

c – The value to add.

inline T const &coefficient() const noexcept

Gets a const reference to the coefficient.

inline void set_coefficient(T new_c) noexcept

Sets the coefficient of the term.

Public Functions

inline T expectation_value() const

Calculates the expectation value of this single term.

Returns:

The coefficient if the Pauli string contains only I and Z operators; otherwise, 0.

inline std::size_t phash() const noexcept

Computes a hash of the Pauli string part of the term, ignoring the coefficient.

Returns:

A hash value representing the Pauli string.

inline std::size_t pauli_weight() const

Calculates the Pauli weight of the term.

Returns:

The number of non-identity Pauli operators in the string.

template<typename Numeric, typename T = coeff_t>
PauliTerm<T> operator*(Numeric &&x, PauliTerm<T> pt)

Creates a PauliTerm by multiplying a numeric value with a PauliTerm.

template<typename Numeric, typename T = coeff_t>
PauliTerm<T> operator*(PauliTerm<T> pt, Numeric &&x)

Creates a PauliTerm by multiplying a PauliTerm with a numeric value.

template<typename Numeric, typename T = coeff_t>
PauliTerm<T> operator*(Numeric &&x, Pauli p)

Creates a single-qubit PauliTerm from a numeric value and a Pauli object.

template<typename Numeric, typename T = coeff_t>
PauliTerm<T> operator*(Pauli p, Numeric &&x)

Creates a single-qubit PauliTerm from a Pauli object and a numeric value.

template<typename Numeric, typename T = coeff_t>
PauliTerm<T> operator*(Numeric &&x, Pauli_enum p)

Creates a single-qubit PauliTerm from a numeric value and a Pauli_enum.

Advanced & Technical API

This section details the lower-level, performance-oriented components for developers and advanced users.

Performance-Oriented Container

The specialized container that provides the efficient, bit-packed backend for the Observable class.

template<typename T, typename Underlying = std::uint8_t>
class PauliTermContainer

A specialized container for storing Pauli terms in a memory-packed format.

This container stores a collection of Pauli terms, representing an observable. It separates the coefficients and the Pauli strings into two vectors. The Pauli strings are bit-packed, with each Pauli operator using 2 bits. This makes the container highly efficient for storing and processing large numbers of terms.

    // Construct a container from string representations of Pauli terms.
    // The coefficients are implicitly 1.0f.
    PauliTermContainer<float> observable = {"+IX", "-ZY", "+XI"};

    std::cout << "Observable has " << observable.nb_terms() << " terms on "
              << observable.nb_qubits() << " qubits." << std::endl;

    // Access the second term ("-ZY") using a non-owning read-only view.
    auto term_view = observable[1];
    
    // Print the term using the view.
    std::cout << "Second term is: " << term_view << std::endl;
    std::cout << "Its coefficient is: " << term_view.coefficient() << std::endl;
    std::cout << "Its Pauli weight is: " << term_view.pauli_weight() << std::endl;

    // Modify the first term ("+IX") through a mutable view.
    // Let's change its coefficient.
    observable[0].set_coefficient(0.5f);

    std::cout << "Modified observable:" << std::endl;
    for (const auto& term : observable) {
        std::cout << "  " << term << std::endl;
    }
Template Parameters:
  • T – The numeric type for the coefficients (e.g., float, double).

  • Underlying – The unsigned integer type used for packing Pauli operators.

Packing Configuration

These constants define the bit-packing strategy.

Constructors

inline PauliTermContainer(std::size_t nb_qubits)

Constructs an empty container for a given number of qubits.

Parameters:

nb_qubits – The number of qubits.

template<typename It, std::enable_if_t< std::is_convertible_v< decltype((*std::declval< It >()).coefficient()), T >, bool> inline  PauliTermContainer (It &&begin, It &&end)

Constructs the container from a range of PauliTerm-like objects.

Template Parameters:

It – An iterator type that dereferences to an object with size(), coefficient(), and operator[].

inline PauliTermContainer(std::initializer_list<std::string_view> lst)

Constructs the container from an initializer list of string views.

Note

Uses AdapterIt to lazily convert strings to PauliTerm objects during construction.

Raw Data Access

Methods for direct, indexed access to the packed data.

Fast Operations on Packed Data

These methods operate directly on the raw packed bits for maximum performance.

inline void copy_fast(std::size_t index_input, std::size_t index_output)

Copies the packed Pauli string from one index to another.

Parameters:
  • index_input – The source term index.

  • index_output – The destination term index.

inline std::size_t fast_phash(std::size_t index) const noexcept

Computes a hash of a Pauli string directly from the packed data.

Parameters:

index – The index of the term to hash.

Returns:

A hash value for the Pauli string.

inline bool fast_equal_bitstring(std::size_t index_lhs, std::size_t index_rhs) const noexcept

Compares two Pauli strings for equality directly from the packed data.

Parameters:
  • index_lhs – The index of the first term.

  • index_rhs – The index of the second term.

Returns:

true if the Pauli strings are identical.

inline std::size_t fast_pauli_weight(std::size_t index) const noexcept

Calculates the Pauli weight of a term directly from the packed data.

Parameters:

index – The index of the term.

Returns:

The number of non-Identity operators in the term’s Pauli string.

Modifiers

inline NonOwningPauliTermPacked create_pauliterm()

Appends a new, zero-initialized Pauli term to the container.

Returns:

A non-owning view to the newly created term.

inline NonOwningPauliTermPacked duplicate_pauliterm(std::size_t idx)

Appends a copy of an existing Pauli term.

Parameters:

idx – The index of the term to duplicate.

Returns:

A non-owning view to the newly created duplicate term.

inline void remove_pauliterm(std::size_t idx)

Removes a Pauli term from the container.

Note

This operation uses a swap-and-pop technique. It has O(1) complexity but does not preserve the relative order of the remaining elements.

Parameters:

idx – The index of the term to remove.

Comparison

inline friend bool operator==(PauliTermContainer const &lhs, PauliTermContainer const &rhs)

Checks for equality between two containers.

Note

This is an O(N^2) operation as it performs an unordered comparison of all terms.

Public Functions

inline NonOwningPauliTermPacked operator[](std::size_t idx)

Provides access to a Pauli term via a non-owning mutable view.

Note

This is a lightweight operation that returns a view to the internal packed data.

Parameters:

idx – The index of the term.

Returns:

A NonOwningPauliTermPacked view of the term at the given index.

inline ReadOnlyNonOwningPauliTermPacked operator[](std::size_t idx) const

Provides access to a Pauli term via a non-owning read-only view.

Parameters:

idx – The index of the term.

Returns:

A ReadOnlyNonOwningPauliTermPacked view of the term at the given index.

template<typename T, typename Predicate>
size_t erase_if(PauliTermContainer<T> &paulis, Predicate &&predicate)

Erases elements from a PauliTermContainer that satisfy a specific criterion.

Parameters:
  • paulis – The container from which to erase.

  • predicate – A unary predicate that returns true for elements to be erased.

Returns:

The number of elements erased.

Note

Using :members: on PauliTermContainer will automatically document its important nested classes, including the packed non-owning views (ReadOnlyNonOwningPauliTermPacked, NonOwningPauliTermPacked) and the iterators (ReadOnlyNonOwningIterator, NonOwningIterator).


Non-Owning View Classes

Lightweight views for manipulating PauliTerm data without incurring copy overhead. These are the conceptual counterparts to the packed views.

template<typename T>
class ReadOnlyNonOwningPauliTerm

A non-owning, read-only view of a Pauli term.

This class provides a const view of a Pauli term, consisting of a Pauli string (std::span<const Pauli>) and a coefficient (const T&). It does not own the underlying data, so it is cheap to create and copy. The caller must ensure that the data it refers to outlives this object.

    // The data is owned by these variables.
    std::vector<Pauli> pauli_string = {p_x, p_y, p_z};
    float coefficient = 0.707f;

    // Create a read-only view of the data.
    ReadOnlyNonOwningPauliTerm<float> term_view(pauli_string, coefficient);

    // Use the view to inspect the data without copying it.
    std::cout << "Read-only view: " << term_view << std::endl;
    std::cout << "Pauli weight: " << term_view.pauli_weight() << std::endl;
    std::cout << "Pauli at index 1: " << term_view[1] << std::endl;
Template Parameters:

T – The numeric type for the term’s coefficient (e.g., float, double).

Constructors

Construct a view from various sources of Pauli data.

Accessors

Provide read-only access to the underlying Pauli string.

Public Functions

inline T expectation_value() const

Calculates the expectation value of this single term.

Returns:

The coefficient if the Pauli string contains only I and Z operators; otherwise, 0.

inline std::size_t phash() const noexcept

Computes a hash of the Pauli string part of the term, ignoring the coefficient.

Returns:

A hash value representing the Pauli string.

inline bool equal_bitstring(ReadOnlyNonOwningPauliTerm const &oth) const

Checks if two terms have the same Pauli string, ignoring the coefficient.

Parameters:

oth – The other term to compare with.

Returns:

true if the Pauli strings are identical.

inline T const &coefficient() const noexcept

Gets the coefficient of the term.

Returns:

A const reference to the coefficient.

inline std::size_t pauli_weight() const

Calculates the Pauli weight of the term.

Returns:

The number of non-identity Pauli operators in the string.

Friends

inline friend bool operator==(ReadOnlyNonOwningPauliTerm const &lhs, ReadOnlyNonOwningPauliTerm const &rhs)

Checks for deep equality between two Pauli terms.

Returns:

true if both the Pauli strings and coefficients are identical.

inline friend std::ostream &operator<<(std::ostream &os, ReadOnlyNonOwningPauliTerm const &pt)

Enables printing the term to an output stream.

template<typename T>
class NonOwningPauliTerm

A non-owning, mutable view of a Pauli term.

This class provides a mutable view of a Pauli term, allowing for in-place modification of the underlying Pauli string and coefficient data. It is used internally by the simulation engine to apply quantum gates to an observable without incurring copy overhead. The caller must ensure that the data it refers to outlives this object.

    // The data is owned by these variables.
    std::vector<Pauli> mutable_pauli_string = {p_i, p_x, p_i};
    double mutable_coefficient = 1.0;

    // Create a mutable view of the data.
    NonOwningPauliTerm<double> mutable_term_view(mutable_pauli_string, mutable_coefficient);
    
    std::cout << "Original term: " << mutable_term_view << std::endl;

    // Apply a gate using the view. This modifies the original data.
    // Apply Pauli Z to qubit 1 (where the X is).
    mutable_term_view.apply_pauli(Pauli_gates::Z, 1);

    std::cout << "Term after applying Z gate to qubit 1: " << mutable_term_view << std::endl;

    // The original data has been changed.
    std::cout << "Original coefficient is now: " << mutable_coefficient << std::endl; // Will be -1.0
    std::cout << "Original string is now: ";
    for(const auto& p : mutable_pauli_string) {
        std::cout << p;
    }
    std::cout << std::endl; // Will be "IYI"
Template Parameters:

T – The numeric type for the term’s coefficient (e.g., float, double).

Constructors

Construct a mutable view from various sources of Pauli data.

Gate Application

Methods to apply quantum operations in-place.

inline void apply_pauli(Pauli_gates g, unsigned qubit)

Applies a Pauli gate to a specific qubit of the term.

Parameters:
  • g – The Pauli gate to apply.

  • qubit – The target qubit index.

inline void apply_clifford(Clifford_Gates_1Q g, unsigned qubit)

Applies a Clifford gate to a specific qubit of the term.

Parameters:
  • g – The single-qubit Clifford gate to apply.

  • qubit – The target qubit index.

inline void apply_unital_noise(UnitalNoise n, unsigned qubit, T p)

Applies a unital noise channel to a specific qubit of the term.

Parameters:
  • n – The noise channel type.

  • qubit – The target qubit index.

  • p – The noise probability.

inline void apply_cx(unsigned control, unsigned target)

Applies a CNOT gate to the term.

Parameters:
  • control – The control qubit index.

  • target – The target qubit index.

Pre:

Control and target qubits must be different.

inline void apply_rz(unsigned qubit, T theta, NonOwningPauliTerm<T> &output)

Applies an Rz gate, potentially splitting the term.

This operation modifies the current term in-place and, if a split occurs, populates the output term with the new term’s data.

Parameters:
  • qubit – The target qubit index.

  • theta – The rotation angle.

  • output – A mutable view where the new split term will be stored.

Pre:

The operator on the target qubit must be X or Y.

inline void apply_amplitude_damping_xy(unsigned qubit, T p)

Applies the X/Y part of the amplitude damping channel.

Parameters:
  • qubit – The target qubit index.

  • p – The damping probability.

Pre:

The operator on the target qubit must be X or Y.

inline void apply_amplitude_damping_z(unsigned qubit, T p, NonOwningPauliTerm<T> &output)

Applies the Z part of the amplitude damping channel, splitting the term.

Parameters:
  • qubit – The target qubit index.

  • p – The damping probability.

  • output – A mutable view where the new split term will be stored.

Pre:

The operator on the target qubit must be Z.

Accessors and Mutators

Provide read-write access to the underlying Pauli string and coefficient.

inline T const &coefficient() const noexcept

Gets a const reference to the coefficient.

inline void set_coefficient(T new_c) noexcept

Sets the coefficient of the term.

Parameters:

new_c – The new coefficient value.

inline void add_coeff(T c)

Adds a value to the term’s coefficient.

Warning

The resulting coefficient must remain within a valid range.

Parameters:

c – The value to add to the coefficient.

Public Functions

inline void copy_content(NonOwningPauliTerm<T> const &nopt)

Copies the content from another non-owning term into this view’s underlying data.

Parameters:

nopt – The source term to copy from.

template<typename It>
inline void copy_paulis(It &&begin, It &&end)

Copies a Pauli string from an iterator range into this view’s underlying data.

Template Parameters:

It – The type of the iterator.

Parameters:
  • begin – An iterator to the beginning of the source Pauli string.

  • end – An iterator to the end of the source Pauli string.

Pre:

The distance between begin and end must match the size of this view’s span.

inline operator ReadOnlyNonOwningPauliTerm<T>() const

Implicitly converts a mutable view to a read-only view.

inline T expectation_value() const

Calculates the expectation value of this single term.

Returns:

The coefficient if the Pauli string contains only I and Z operators; otherwise, 0.

inline std::size_t phash() const noexcept

Computes a hash of the Pauli string part of the term, ignoring the coefficient.

Returns:

A hash value representing the Pauli string.

inline bool equal_bitstring(NonOwningPauliTerm const &oth) const

Checks if two terms have the same Pauli string, ignoring the coefficient.

Parameters:

oth – The other term to compare with.

Returns:

true if the Pauli strings are identical.

inline std::size_t pauli_weight() const

Calculates the Pauli weight of the term.

Returns:

The number of non-identity Pauli operators in the string.

Friends

inline friend bool operator==(NonOwningPauliTerm const &lhs, NonOwningPauliTerm const &rhs)

Checks for deep equality between two Pauli terms.

Returns:

true if both the Pauli strings and coefficients are identical.

inline friend std::ostream &operator<<(std::ostream &os, NonOwningPauliTerm const &pt)

Enables printing the term to an output stream.


Algorithms & Utilities

Key algorithms and low-level helper functions used throughout the library.

Warning

doxygenfunction: Cannot find function “merge_inplace_noalloc” in doxygen xml output for project “ProPauli” from directory: /home/runner/work/ProPauli/ProPauli/build/doc/doxygen/xml

Provides a collection of constexpr helper functions for low-level bit manipulation.

This header contains various standalone, compile-time-evaluated functions for creating bitmasks and performing other bitwise operations. These utilities are the fundamental building blocks for the memory-efficient bit-packing scheme used by the PauliTermContainer class to store Pauli strings.

Functions

template<typename T>
bool is_power_of_two(T k)
template<typename size_type>
constexpr size_type next_power_of_two(size_type n)
template<typename T>
constexpr T compute_mask(T nb_bits)

Creates a bitmask with a specified number of lower bits set to 1.

Note

For example, compute_mask<uint8_t>(3) returns 7 (0b00000111).

Template Parameters:

T – The integer type for the mask.

Parameters:

nb_bits – The number of low-order bits to set to 1.

Returns:

An integer of type T with the nb_bits least significant bits set.

template<typename Underlying, std::size_t OBJS_PER_UNDERLYING>
constexpr std::array<Underlying, OBJS_PER_UNDERLYING> compute_mask_lut()

Computes a compile-time lookup table (LUT) of bitmasks for packed objects.

Note

This is used by PauliTermContainer to avoid repeated shift calculations at runtime when accessing packed Pauli data.

Template Parameters:
  • Underlying – The integer type in which objects are packed.

  • OBJS_PER_UNDERLYING – The number of objects packed into one Underlying integer.

Returns:

A std::array where each element is a mask for one object slot within the Underlying type.

template<std::unsigned_integral T>
constexpr T create_low_bit_mask()

Creates a mask to select the low bit of every 2-bit pair in an integer.

Template Parameters:

T – An unsigned integral type.

Returns:

A mask of the form ...01010101.

template<std::unsigned_integral T>
constexpr int count_nonzero_pairs(T input)

Efficiently counts the number of non-zero 2-bit pairs in an unsigned integer.

Note

This function is a key optimization for calculating the Pauli weight of a term directly from its packed representation. Since each Pauli operator is stored as a 2-bit value and the Identity operator is 00, counting the non-zero pairs is equivalent to counting the non-Identity operators. The algorithm works by ORing the low and high bits of each pair together, then using std::popcount on the result.

Template Parameters:

T – An unsigned integral type.

Parameters:

input – The integer whose bit pairs are to be counted.

Returns:

The number of 2-bit chunks in input that are not 00.

Provides a generic, lazy-initializing iterator adapter.

This file defines the AdapterIt class template, a forward iterator that wraps an existing iterator and transforms its dereferenced value into a different type on-the-fly. This is useful for creating compatible iterators over containers whose value types do not match a required interface, without needing to create a new, fully populated container.

template<typename TOut, typename ItIn>
class AdapterIt
#include <adapter.hpp>

A forward iterator that lazily adapts an input iterator’s value to a different type.

AdapterIt wraps an input iterator ItIn. When dereferenced for the first time, it constructs an object of type TOut from the value of the input iterator (*it). This converted object is cached in a std::shared_ptr. Subsequent dereferences return the cached object until the iterator is incremented, at which point the cache is cleared.

This lazy conversion mechanism is efficient when the full sequence of converted objects is not needed, or when the conversion itself is an expensive operation.

    /*  A simple struct that can be constructed from an int
    struct NumberString {
        std::string data;
        explicit NumberString(int val) : data("Number: " + std::to_string(val)) {}
    }; */

    // Create a vector of integers.
    std::vector<int> numbers = { 10, 20, 30, 40, 50 };

    // Create AdapterIt begin and end iterators.
    // This will adapt `std::vector<int>::iterator` to an iterator that
    // yields `NumberString` objects on-the-fly.
    AdapterIt<NumberString, decltype(numbers.begin())> begin(numbers.begin());
    AdapterIt<NumberString, decltype(numbers.end())> end(numbers.end());

    // Iterate using the adapter. The conversion from `int` to `NumberString`
    // happens lazily inside the loop upon dereference.
    std::cout << "Using AdapterIt to convert ints to NumberStrings:\n";
    for (auto it = begin; it != end; ++it) {
        // `*it` triggers the conversion from int to NumberString
        std::cout << it->data << std::endl;
    }
Template Parameters:
  • TOut – The target type to which the underlying iterator’s elements are converted.

  • ItIn – The type of the underlying input iterator being adapted.

Standard Iterator Traits

using iterator_category = std::forward_iterator_tag
using difference_type = std::ptrdiff_t
using value_type = TOut
using pointer = value_type*
using reference = value_type&

Public Functions

inline AdapterIt(ItIn it_in)

Constructs an AdapterIt from a given input iterator.

Parameters:

it_in – The input iterator to wrap.

inline reference operator*()

Dereferences the iterator, returning a reference to the converted object.

Note

The conversion from the input iterator’s value to TOut happens lazily on the first call to this operator or operator->().

Returns:

A reference to the converted object of type TOut.

inline pointer operator->()

Dereferences the iterator, returning a pointer to the converted object.

Note

Triggers the lazy conversion if the object has not yet been created.

Returns:

A pointer to the converted object of type TOut.

inline AdapterIt &operator++()

Advances the iterator to the next element (prefix increment).

Note

This invalidates the cached object. The next dereference will trigger a new conversion for the next element.

Returns:

A reference to the incremented iterator.

inline AdapterIt operator++(int)

Advances the iterator to the next element (postfix increment).

Returns:

A copy of the iterator before it was incremented.

Private Functions

inline void init()

Initializes the obj member if it hasn’t been already.

This function is called on the first dereference (* or ->). It constructs the TOut object from the value of the current input iterator and stores it in the obj shared pointer.

Private Members

std::remove_cvref_t<ItIn> it

The underlying input iterator.

std::shared_ptr<TOut> obj

A shared pointer to the lazily-initialized, converted object.

Friends

inline friend bool operator==(AdapterIt const &lhs, AdapterIt const &rhs)

Compares two iterators for equality.

Parameters:
  • lhs – The left-hand side iterator.

  • rhs – The right-hand side iterator.

Returns:

true if the underlying iterators are equal, false otherwise.

inline friend bool operator!=(AdapterIt const &lhs, AdapterIt const &rhs)

Compares two iterators for inequality.

Parameters:
  • lhs – The left-hand side iterator.

  • rhs – The right-hand side iterator.

Returns:

true if the underlying iterators are not equal, false otherwise.

Truncation Framework

A collection of classes for truncating observables to manage simulation complexity.

Base Class

template<typename T>
class Truncator

An abstract base class (interface) for defining truncation strategies.

A truncator is responsible for removing Pauli terms from an observable to manage the simulation’s complexity. All specific truncation strategies should inherit from this class.

Template Parameters:

T – The numeric type of the PauliTerm coefficients.

Subclassed by KeepNTruncator< T >, PredicateTruncator< P, T >, RuntimeMultiTruncators< T >, Truncators< T, Ts >

Public Functions

virtual std::size_t truncate(PauliTermContainer<T> &paulis, T &error_acc) = 0

Applies the truncation strategy to a container of Pauli terms in place.

Parameters:

paulis – The container of Pauli terms to truncate.

Returns:

The number of terms removed.


Provided Truncators

template<typename T = coeff_t>
using CoefficientTruncator = decltype(PredicateTruncator(std::declval<CoefficientPredicate<T>>()))

A truncator that removes Pauli terms with small coefficients.

    // Create a truncator that removes terms with coefficient magnitude < 0.1
    CoefficientTruncator<float> coeff_trunc(0.1f);
    float error = 0;
    auto removed_by_coeff = coeff_trunc.truncate(observable, error);
    // `observable` now contains {"+0.5XX", "-0.8II", "+0.2ZZZ"}
    std::cout << "Removed by coefficient truncator: " << removed_by_coeff << std::endl;
template<typename T = coeff_t>
using WeightTruncator = decltype(PredicateTruncator<WeightPredicate, T>(std::declval<WeightPredicate>()))

A truncator that removes Pauli terms with high Pauli weight.

    // Create a truncator that removes terms with Pauli weight >= 3
    WeightTruncator<> weight_trunc(3);
    float error_w = 0;
    auto removed_by_weight = weight_trunc.truncate(observable, error_w);
    // `observable` now contains {"+0.5XX", "-0.8II"}
    std::cout << "Removed by weight truncator: " << removed_by_weight << std::endl;