C++ API¶
Pauli¶
-
struct Pauli¶
A class for efficient representation of a 2x2 Pauli matrix \( \sigma_i \in \{ I,X,Y,Z \} \).
Public Functions
-
inline constexpr Pauli()¶
Default constructor, initializes to I.
- template<class T> inline requires constexpr std::convertible_to< T, uint8_t > Pauli (T const code)
Constructor given a numeric code.
- Template Parameters:
T – Any type convertible to uint8_t
- Parameters:
code – 0: I, 1: X, 2: Y, 3: Z
- Returns:
requires
Friends
-
inline friend std::pair<std::complex<double>, Pauli> operator*(Pauli const &lhs, Pauli const &rhs)¶
Returns the product of two pauli matrices and their phase as a pair.
- Parameters:
lhs – left hand side pauli object
rhs – right hand side pauli object
- Returns:
std::pair<std::complex<double>, Pauli> phase and resulting pauli matrix
-
inline constexpr Pauli()¶
PauliString¶
-
struct PauliString¶
A class representation of a Pauli string (i.e. a tensor product of 2x2 pauli matrices) \( $\mathcal{\hat{P}} = \bigotimes_i \sigma_i \) where \( \sigma_i \in \{ I,X,Y,Z \} \).
Public Functions
-
PauliString() noexcept = default¶
Default constructor, initialize weight and empty vector for paulis.
-
inline PauliString(std::vector<Pauli> paulis)¶
Constructs a PauliString from a vector of pauli matrices and calculates the weight.
-
inline PauliString(std::span<fast_pauli::Pauli> paulis)¶
Constructs a PauliString from a span of pauli matrices and calculates the weight.
-
inline PauliString(std::string const &str)¶
Constructs a PauliString from a string and calculates the weight. This is often the most compact way to initialize a PauliString.
-
inline PauliString(char const *str)¶
Allows implicit conversion of string literals to PauliStrings. Ex: std::vector<PauliString> pauli_strings = {“IXYZ”, “IIIII”};.
-
inline size_t n_qubits() const noexcept¶
Return the number of qubits in the PauliString.
- Returns:
size_t
-
inline size_t dim() const noexcept¶
Return the dimension (2^n_qubits) of the PauliString.
Note
this returns 0 if the PauliString is empty.
- Returns:
size_t
-
template<std::floating_point T>
inline void apply(std::mdspan<std::complex<T>, std::dextents<size_t, 1>> new_states, std::mdspan<std::complex<T>, std::dextents<size_t, 1>> states, std::complex<T> const c = 1.0) const¶ Apply a pauli string (using the sparse representation) to a vector. This performs following matrix-vector multiplication \( \mathcal{\hat{P}} \ket{\psi} \).
- Template Parameters:
T – The floating point base to use for all the complex numbers
- Parameters:
new_states – Output state
states – The input vector to apply the PauliString to. Must be the
c – Multiplication factor to apply to the PauliString same size as PauliString.dim().
-
template<std::floating_point T, execution_policy ExecutionPolicy>
inline void apply(ExecutionPolicy&&, std::mdspan<std::complex<T>, std::dextents<size_t, 1>> new_states, std::mdspan<std::complex<T>, std::dextents<size_t, 1>> states, std::complex<T> const c = 1.0) const¶ - Template Parameters:
T –
ExecutionPolicy –
- Parameters:
new_states –
states –
-
template<std::floating_point T>
inline void apply_batch(std::mdspan<std::complex<T>, std::dextents<size_t, 2>> new_states_T, std::mdspan<std::complex<T>, std::dextents<size_t, 2>> const states_T, std::complex<T> const c) const¶ Apply the PauliString to a batch of states. This function takes a different shape of the states than the other apply functions. here all the states (new and old) are transposed so their shape is (n_dims x n_states). All the new_stats are overwritten, no need to initialize.
This performs following matrix-matrix multiplication \( \mathcal{\hat{P}} \hat{\Psi} \) where matrix \( \hat{\Psi} \) has \( \ket{\psi_t} \) as columns
- Template Parameters:
T – The floating point base to use for all the complex numbers
- Parameters:
new_states_T – The output states after applying the PauliString (n_dim x n_states)
states_T – THe original states to apply the PauliString to (n_dim x n_states)
c – Multiplication factor to apply to the PauliString
-
template<std::floating_point T, execution_policy ExecutionPolicy>
inline void apply_batch(ExecutionPolicy&&, std::mdspan<std::complex<T>, std::dextents<size_t, 2>> new_states_T, std::mdspan<std::complex<T>, std::dextents<size_t, 2>> const states_T, std::complex<T> const c) const¶ This performs following matrix-matrix multiplication \( \mathcal{\hat{P}} \hat{\Psi} \) where matrix \( \hat{\Psi} \) has \( \ket{\psi_t} \) as columns
- Template Parameters:
T – The floating point base to use for all the complex numbers
T –
ExecutionPolicy –
- Parameters:
new_states_T – The output states after applying the PauliString (n_dim x n_states)
states_T – THe original states to apply the PauliString to (n_dim x n_states)
c – Multiplication factor to apply to the PauliString
new_states_T –
states_T –
c –
-
template<std::floating_point T>
inline void expectation_value(std::mdspan<std::complex<T>, std::dextents<size_t, 1>> expectation_vals_out, std::mdspan<std::complex<T>, std::dextents<size_t, 2>> states, std::complex<T> const c = 1.0) const¶ Calculate expectation values for a given batch of states. This function takes in transposed states with (n_dims x n_states) shape.
It computes following inner product \( \bra{\psi_t} \mathcal{\hat{P_i}} \ket{\psi_t} \) for each state \( \ket{\psi_t} \) from provided batch.
Note
The expectation values are added to corresponding coordinates in the expectation_vals_out vector.
- Template Parameters:
T – The floating point base to use for all the complex numbers
- Parameters:
expectation_vals_out – accumulator for expectation values for each state in the batch
states – THe original states to apply the PauliString to (n_dim x n_states)
c – Multiplication factor to apply to the PauliString
-
template<std::floating_point T, execution_policy ExecutionPolicy>
inline void expectation_value(ExecutionPolicy&&, std::mdspan<std::complex<T>, std::dextents<size_t, 1>> expectation_vals_out, std::mdspan<std::complex<T>, std::dextents<size_t, 2>> states, std::complex<T> const c = 1.0) const¶ It computes following inner product \( \bra{\psi_t} \mathcal{\hat{P_i}} \ket{\psi_t} \) for each state \( \ket{\psi_t} \) from provided batch.
Note
The expectation values are added to corresponding coordinates in the expectation_vals_out vector.
- Template Parameters:
T – The floating point base to use for all the complex numbers
T –
ExecutionPolicy –
- Parameters:
expectation_vals_out – accumulator for expectation values for each state in the batch
states – THe original states to apply the PauliString to (n_dim x n_states)
c – Multiplication factor to apply to the PauliString
expectation_vals_out –
states –
c –
-
template<std::floating_point T>
inline void to_tensor(std::mdspan<std::complex<T>, std::dextents<size_t, 2>> output) const¶ Get the dense representation of the object as a 2D-array.
- Template Parameters:
T – The floating point base to use for all the complex numbers
- Parameters:
output – The output tensor to fill with the dense representation
Friends
-
inline friend std::pair<std::complex<double>, PauliString> operator*(PauliString const &lhs, PauliString const &rhs)¶
Returns the result of matrix multiplication of two PauliStrings and their phase as a pair.
- Parameters:
lhs – left hand side PauliString
rhs – right hand side PauliString
- Returns:
std::pair<std::complex<double>, PauliString> phase and resulting PauliString
-
PauliString() noexcept = default¶
PauliOp¶
-
template<std::floating_point T, typename H = std::complex<T>>
struct PauliOp¶ A class representation for a Pauli Operator (i.e. a weighted sum of Pauli Strings) \( \big( \sum_i h_i \mathcal{\hat{P}}_i \big) \) where \( \mathcal{\hat{P}}_i \) are composed using \( \sigma_i \in \{ I,X,Y,Z \} \) and \( h_i \) are the coefficients.
Public Functions
-
PauliOp() = default¶
Default constructor, initialize empty vectors for paulis and coefficients.
-
inline PauliOp(std::vector<std::string> const &strings)¶
Construct a PauliOp from a vector of strings and default corresponding coeffs to ones.
- Parameters:
strings – vector of strings
-
inline PauliOp(std::vector<PauliString> strings)¶
Construct a PauliOp from a vector of PauliStrings and default corresponding coeffs to ones.
- Parameters:
strings – vector of PauliString objects
-
inline PauliOp(std::vector<H> coefficients, std::vector<PauliString> strings)¶
Construct a PauliOp from a vector of PauliStrings and corresponding coefficients.
- Parameters:
coefficients – vector of coefficients
strings – vector of PauliString objects
-
inline size_t dim() const¶
Return the dimension (2^n_qubits) of the PauliStrings used to compose PauliOp.
- Returns:
size_t
-
inline size_t n_pauli_strings() const¶
Return the number of PauliStrings in PauliOp.
- Returns:
size_t
-
inline void scale(std::complex<T> factor)¶
Scale each individual term by a factor.
- Parameters:
factors – a factor to scale each term with
-
inline void scale(mdspan<std::complex<T>, std::dextents<size_t, 1>> factors)¶
Scale each individual term by a factor.
- Parameters:
factors – n_pauli_strings length array of factors to scale each term
-
inline PauliOp<T, H> operator-() const¶
Return the copy of PauliOp with the coefficients negated.
- Returns:
PauliOp<T, H>
-
inline void extend(PauliString pauli_str, std::complex<T> coeff, bool dedupe = true)¶
Add a PauliString term with appropriate coefficient to the summation inside PauliOp.
- Parameters:
pauli_str – PauliString to add to the summation
coeff – coefficient to apply to the PauliString
dedupe – whether to deduplicate provided PauliString
-
inline void extend(PauliOp<T, H> const &other_op)¶
Add another PauliOp to the current one by extending the internal summation with new terms.
Note
: for now it’s very sloppy implementation just to have this functionality
- Parameters:
other_op – PauliOp to add to the current one
-
inline void apply(mdspan<std::complex<T>, std::dextents<size_t, 1>> state_out, mdspan<std::complex<T>, std::dextents<size_t, 1>> const state) const¶
Apply the PauliOp to a state.
This performs following matrix-matrix multiplication \( \big( \sum_i h_i \mathcal{\hat{P}}_i \big) \hat{\psi} \)
-
template<execution_policy ExecutionPolicy>
inline void apply(ExecutionPolicy &&policy, mdspan<std::complex<T>, std::dextents<size_t, 1>> state_out, mdspan<std::complex<T>, std::dextents<size_t, 1>> const state) const¶ This performs following matrix-matrix multiplication \( \big( \sum_i h_i \mathcal{\hat{P}}_i \big) \hat{\psi} \)
-
inline void apply(mdspan<std::complex<T>, std::dextents<size_t, 2>> new_states, mdspan<std::complex<T>, std::dextents<size_t, 2>> const states) const¶
Apply the PauliOp to a batch of states. Here all the states (new and old) are transposed so their shape is (n_dims x n_states). All the new_stats are overwritten, no need to initialize.
This performs following matrix-matrix multiplication \( \big( \sum_i h_i \mathcal{\hat{P}}_i \big) \hat{\Psi} \) where matrix \( \hat{\Psi} \) has \( \ket{\psi_t} \) as columns
-
template<execution_policy ExecutionPolicy>
inline void apply(ExecutionPolicy&&, mdspan<std::complex<T>, std::dextents<size_t, 2>> new_states, mdspan<std::complex<T>, std::dextents<size_t, 2>> const states) const¶ This performs following matrix-matrix multiplication \( \big( \sum_i h_i \mathcal{\hat{P}}_i \big) \hat{\Psi} \) where matrix \( \hat{\Psi} \) has \( \ket{\psi_t} \) as columns
-
inline void expectation_value(std::mdspan<std::complex<T>, std::dextents<size_t, 1>> expectation_vals_out, mdspan<std::complex<T>, std::dextents<size_t, 2>> states) const¶
Calculate the expectation value of the PauliOp on a batch of states.
It computes following inner product \( \bra{\psi_t} ( \sum_i h_{ik} \mathcal{\hat{P}}_i ) \ket{\psi_t} \) for each state \( \ket{\psi_t} \) from provided batch.
- Parameters:
expectation_vals_out – expectation values for each state in the batch
states – The states we want to use in our expectation value calculation (n_dim x n_states)
-
template<execution_policy ExecutionPolicy>
inline void expectation_value(ExecutionPolicy&&, std::mdspan<std::complex<T>, std::dextents<size_t, 1>> expectation_vals_out, mdspan<std::complex<T>, std::dextents<size_t, 2>> states) const¶ - Template Parameters:
ExecutionPolicy –
- Parameters:
expectation_vals_out –
states –
Public Static Functions
-
static inline void validate_pauli_strings(std::vector<PauliString> const &pauli_strings)¶
Check that the dims of pauli strings are all the same.
- Parameters:
pauli_strings –
Friends
-
inline friend PauliOp<T, H> operator*(PauliOp<T, H> const &pauli_op_left, PauliString const &pauli_str_right)¶
Matrix multiplication of PauliOp with a PauliString on the right.
- Parameters:
pauli_op_left – left hand side PauliOp
pauli_str_right – right hand side PauliString
- Returns:
PauliOp<T, H> new PauliOp instance containing the result of the multiplication
-
inline friend PauliOp<T, H> operator*(PauliString const &pauli_str_left, PauliOp<T, H> const &pauli_op_right)¶
Matrix multiplication of PauliOp with a PauliString on the left.
- Parameters:
pauli_str_left – left hand side PauliString
pauli_op_right – right hand side PauliOp
- Returns:
PauliOp<T, H> new PauliOp instance containing the result of the multiplication
-
PauliOp() = default¶
SummedPauliOp¶
-
template<std::floating_point T>
struct SummedPauliOp¶ A class representing a sum of Pauli operators \( A = \sum_k A_k = \sum_{ik} h_{ik} \mathcal{\hat{P}}_i \). Where \( \mathcal{\hat{P}}_i \) are Pauli strings and \( h_{ik} \) are complex-valued coefficients.
Public Functions
-
SummedPauliOp() noexcept = default¶
Default constructor.
-
inline SummedPauliOp(std::vector<PauliString> const &pauli_strings, std::vector<std::complex<T>> const &coeffs_raw)¶
Construct a new Summed Pauli Op object from a vector of PauliStrings and a blob of coefficients.
- Parameters:
pauli_strings – The PauliStrings that define the set of PauliStrings used by all operators (n_pauli_strings)
coeffs_raw – A vector of coefficients that define the weights of each PauliString in each operator. The coefficients here are a flattened version of \( h_{ik} \) in \( A_k = \sum_i h_{ik} \mathcal{\hat{P}}_i \) (n_pauli_strings * n_operators,)
-
inline SummedPauliOp(std::vector<PauliString> const &pauli_strings, Tensor<2> const coeffs)¶
Construct a new Summed Pauli Op object from a vector of PauliStrings and an std::mdspan of coefficients.
- Parameters:
pauli_strings – The PauliStrings that define the set of PauliStrings used by all operators (n_pauli_strings,)
coeffs – A 2D std::mdspan of coefficients that define the weights of each PauliString in each operator. The coefficients here are \( h_{ik} \) in \( A_k = \sum_i h_{ik} \mathcal{\hat{P}}_i \). (n_pauli_strings, n_operators)
-
inline SummedPauliOp(std::vector<std::string> const &pauli_strings, Tensor<2> const coeffs)¶
Construct a new Summed Pauli Op object from a vector of strings and a std::mdspan of coefficients.
- Parameters:
pauli_strings – A vector of strings that define the set of PauliStrings used by all operators (n_pauli_strings)
coeffs – A 2D std::mdspan of coefficients that define the weights of each PauliString in each operator. The coefficients here are \( h_{ik} \) in \( A_k = \sum_i h_{ik} \mathcal{\hat{P}}_i \). (n_pauli_strings, n_operators)
-
inline size_t dim() const noexcept¶
Return the number of dimensions of the SummedPauliOp.
- Returns:
size_t
-
inline size_t n_qubits() const noexcept¶
Return the number of qubits in the SummedPauliOp.
- Returns:
size_t
-
inline size_t n_operators() const noexcept¶
Return the number of Pauli operators in the SummedPauliOp.
- Returns:
s
-
inline size_t n_pauli_strings() const noexcept¶
Return the number of PauliStrings in the SummedPauliOp.
- Returns:
size_t
-
inline fast_pauli::SummedPauliOp<T> square() const¶
Square the SummedPauliOp, mathematically \( A_k \rightarrow A_k^2 = \sum_{a,b} T_{ab} A_a A_b \).
We use the sprarse structure of A_k operators and PauliString products to calculate the new coefficients.
- Returns:
SummedPauliOp<T>
-
inline void apply(Tensor<2> new_states, Tensor<2> states) const¶
Apply the SummedPauliOp to a set of states, mathematically \( \big(\sum_k \sum_i h_{ik} \mathcal{\hat{P}}_i \big) \ket{\psi_t} \).
- Parameters:
new_states – The output states after applying the SummedPauliOp (n_dim, n_states)
states – The input states to apply the SummedPauliOp to (n_dim, n_states)
-
template<execution_policy ExecutionPolicy>
inline void apply(ExecutionPolicy&&, Tensor<2> new_states, Tensor<2> states) const¶ Apply the SummedPauliOp to a set of states, mathematically \( \big(\sum_k \sum_i h_{ik} \mathcal{\hat{P}}_i \big) \ket{\psi_t} \).
- Parameters:
new_states – The output states after applying the SummedPauliOp (n_dim, n_states)
states – The input states to apply the SummedPauliOp to (n_dim, n_states)
- Template Parameters:
ExecutionPolicy –
-
template<std::floating_point data_dtype>
inline void apply_weighted(Tensor<2> new_states, Tensor<2> states, std::mdspan<data_dtype, std::dextents<size_t, 2>> data) const¶ Apply the SummedPauliOp to a set of weighted states.
Calculates \( \big(\sum_k x_{tk} \sum_i h_{ik} \mathcal{\hat{P}}_i \big) \ket{\psi_t} \)
- Template Parameters:
data_dtype – The floating point type of the weights \( x_{kt} \) (n_operators, n_states)
- Parameters:
new_states – The output states after applying the SummedPauliOp (n_dim, n_states)
states – The input states to apply the SummedPauliOp to (n_dim, n_states)
data – A 2D std::mdspan of the data \( x_{tk} \) that weights the operators in the expression above (n_operators, n_states)
-
template<execution_policy ExecutionPolicy, std::floating_point data_dtype>
inline void apply_weighted(ExecutionPolicy&&, Tensor<2> new_states, Tensor<2> states, std::mdspan<data_dtype, std::dextents<size_t, 2>> data) const¶ - Template Parameters:
ExecutionPolicy – Execution policy for parallelization
data_dtype – The floating point type of the weights \( x_{kt} \) (n_operators, n_states)
-
inline void expectation_value(Tensor<2> expectation_vals_out, Tensor<2> states) const¶
Calculate the expectation value of the SummedPauliOp on a batch of states. This function returns the expectation values of each operator for each input states, so the output tensor will have shape (n_operators x n_states).
\( \bra{\psi_t} \big(\sum_k \sum_i h_{ik} \mathcal{\hat{P}}_i \big) \ket{\psi_t} \)
- Parameters:
expectation_vals_out – Output tensor for the expectation values (n_operators x n_states)
states – The states used to calculate the expectation values (n_dim, n_states)
-
template<execution_policy ExecutionPolicy>
inline void expectation_value(ExecutionPolicy&&, Tensor<2> expectation_vals_out, Tensor<2> states) const¶ Calculate the expectation value of the SummedPauliOp on a batch of states. This function returns the expectation values of each operator for each input states, so the output tensor will have shape (n_operators x n_states).
\( \bra{\psi_t} \big(\sum_k \sum_i h_{ik} \mathcal{\hat{P}}_i \big) \ket{\psi_t} \)
- Parameters:
expectation_vals_out – Output tensor for the expectation values (n_operators x n_states)
states – The states used to calculate the expectation values (n_dim, n_states)
- Template Parameters:
ExecutionPolicy – Execution policy for parallelization
-
inline void to_tensor(Tensor<3> A_k_out) const¶
Get the dense representation of the SummedPauliOp as a 3D tensor.
- Parameters:
A_k_out – The output tensor to fill with the dense representation
-
SummedPauliOp() noexcept = default¶
Helpers¶
-
std::vector<std::string> fast_pauli::get_nontrivial_paulis(size_t const weight)¶
Get the nontrivial sets of pauli matrices given a weight.
- Parameters:
weight – The Pauli weight to get the nontrivial paulis for.
- Returns:
std::vector<std::string>
-
std::vector<PauliString> fast_pauli::calculate_pauli_strings(size_t const n_qubits, size_t const weight)¶
Calculate all possible PauliStrings for a given number of qubits and weight and return them in lexicographical order.
- Parameters:
n_qubits – The number of qubits.
weight – The Pauli weight.
- Returns:
std::vector<PauliString>
-
std::vector<PauliString> fast_pauli::calculate_pauli_strings_max_weight(size_t n_qubits, size_t weight)¶
Calculate all possible PauliStrings for a given number of qubits and all weights less than or equal to a given weight.
- Parameters:
n_qubits – The number of qubits.
weight – The Pauli weight.
- Returns:
std::vector<PauliString>
-
template<std::floating_point T>
std::tuple<std::vector<size_t>, std::vector<std::complex<T>>> fast_pauli::get_sparse_repr(std::vector<Pauli> const &paulis)¶ Get the sparse representation of the pauli string matrix.
PauliStrings are always sparse and have only single non-zero element per row. It’s N non-zero elements for NxN matrix where N is 2^n_qubits. Therefore k and m will always have N elements.
See Algorithm 1 in https://arxiv.org/pdf/2301.00560.pdf for details about the algorithm.
- Template Parameters:
T – The floating point type (i.e. std::complex<T> for the values of the PauliString matrix)
- Parameters:
k – The column index of the matrix
m – The values of the matrix