Template Class StateVectorKokkosMPI

Inheritance Relationships

Base Type

Class Documentation

template<class fp_t = double>
class StateVectorKokkosMPI : public Pennylane::StateVectorBase<double, StateVectorKokkosMPI<double>>

Kokkos state vector class.

Template Parameters

fp_t – Floating-point precision type.

Public Types

using PrecisionT = fp_t
using SVK = StateVectorKokkos<PrecisionT>
using ComplexT = typename SVK::ComplexT
using CFP_t = typename SVK::CFP_t
using KokkosVector = typename SVK::KokkosVector
using UnmanagedComplexHostView = typename SVK::UnmanagedComplexHostView
using UnmanagedConstComplexHostView = typename SVK::UnmanagedConstComplexHostView
using KokkosSizeTVector = typename SVK::KokkosSizeTVector
using UnmanagedSizeTHostView = typename SVK::UnmanagedSizeTHostView
using UnmanagedConstSizeTHostView = typename SVK::UnmanagedConstSizeTHostView
using UnmanagedPrecisionHostView = typename SVK::UnmanagedPrecisionHostView
using KokkosExecSpace = typename SVK::KokkosExecSpace
using HostExecSpace = typename SVK::HostExecSpace
using BaseType = StateVectorBase<fp_t, StateVectorKokkosMPI<fp_t>>

Public Functions

StateVectorKokkosMPI() = delete
inline StateVectorKokkosMPI(MPIManagerKokkos mpi_manager, std::size_t num_global_qubits, std::size_t num_local_qubits, const Kokkos::InitializationSettings &kokkos_args = {})

Create a new state vector with mpi_manager.

Parameters
  • mpi_manager – Kokkos MPIManager

  • num_global_qubits – Number of global qubits

  • num_local_qubits – Number of local qubits

  • kokkos_args – Arguments for Kokkos initialization

inline StateVectorKokkosMPI(MPIManagerKokkos mpi_manager, std::size_t total_num_qubits, const Kokkos::InitializationSettings &kokkos_args = {})

Create a new state vector with mpi_manager.

Parameters
  • mpi_manager – Kokkos MPIManager

  • total_num_qubits – Number of qubits

  • kokkos_args – Arguments for Kokkos initialization

inline StateVectorKokkosMPI(MPI_Comm mpi_communicator, std::size_t total_num_qubits, const Kokkos::InitializationSettings &kokkos_args = {})

Create a new state vector with MPI Comm.

Parameters
  • mpi_communicator – MPI_Comm communicator

  • total_num_qubits – Number of qubits

  • kokkos_args – Arguments for Kokkos initialization

inline StateVectorKokkosMPI(std::size_t total_num_qubits, const Kokkos::InitializationSettings &kokkos_args = {})

Create a new state vector with default MPI_COMM_WORLD.

Parameters
  • total_num_qubits – Number of qubits

  • kokkos_args – Arguments for Kokkos initialization

inline StateVectorKokkosMPI(MPI_Comm mpi_communicator, std::size_t num_global_qubits, std::size_t num_local_qubits, const Kokkos::InitializationSettings &kokkos_args = {})

Create a new state vector with MPI Comm.

Parameters
  • mpi_communicator – MPI_Comm communicator

  • num_global_qubits – Number of global qubits

  • num_local_qubits – Number of local qubits

  • kokkos_args – Arguments for Kokkos initialization

inline StateVectorKokkosMPI(std::size_t num_global_qubits, std::size_t num_local_qubits, const Kokkos::InitializationSettings &kokkos_args = {})

Create a new state vector with default MPI_COMM_WORLD.

Parameters
  • num_global_qubits – Number of global qubits

  • num_local_qubits – Number of local qubits

  • kokkos_args – Arguments for Kokkos initialization

template<class complex>
inline StateVectorKokkosMPI(std::size_t num_global_qubits, std::size_t num_local_qubits, complex *hostdata_, const std::size_t length, const Kokkos::InitializationSettings &kokkos_args = {}, const MPI_Comm &communicator = MPI_COMM_WORLD)

Create a new state vector from data on the host.

Parameters
  • num_global_qubits – Number of global qubits

  • num_local_qubits – Number of local qubits

  • hostdata_ – Host array for state vector

  • length – Length of host array (must be equal to exp2(total qubits))

  • kokkos_args – Arguments for Kokkos initialization

  • communicator – MPI Communicator

inline StateVectorKokkosMPI(std::size_t num_global_qubits, std::size_t num_local_qubits, const ComplexT *hostdata_, const std::size_t length, const Kokkos::InitializationSettings &kokkos_args = {}, const MPI_Comm &communicator = MPI_COMM_WORLD)

Create a new state vector from data on the host.

Parameters
  • num_global_qubits – Number of global qubits

  • num_local_qubits – Number of local qubits

  • hostdata_ – Host array for state vector

  • length – Length of host array (must be equal to exp2(total qubits))

  • kokkos_args – Arguments for Kokkos initialization

  • communicator – MPI Communicator

template<class complex>
inline StateVectorKokkosMPI(std::size_t num_global_qubits, std::size_t num_local_qubits, std::vector<complex> hostdata_, const Kokkos::InitializationSettings &kokkos_args = {}, const MPI_Comm &communicator = MPI_COMM_WORLD)

Create a new state vector from data on the host.

Parameters
  • num_global_qubits – Number of global qubits

  • num_local_qubits – Number of local qubits

  • hostdata_ – Host array for state vector

  • kokkos_args – Arguments for Kokkos initialization

  • communicator – MPI Communicator

inline StateVectorKokkosMPI(const StateVectorKokkosMPI &other, const Kokkos::InitializationSettings &kokkos_args = {})

Copy constructor.

Parameters

other – Another state vector

~StateVectorKokkosMPI() = default
inline SVK &getLocalSV()
inline auto getMPIManager() const
inline auto getSendBuffer() const
inline auto getRecvBuffer() const
inline std::size_t getNumGlobalWires() const
inline std::size_t getNumLocalWires() const
inline void resetIndices()

Reset the indices of the global_wires_, local_wires_ and the mpi_rank_to_global_index_map_.

inline void initZeros()

Set to initial zero state for the state-vector on device.

inline void setBasisState(std::size_t global_index)

Set value for a single element of the state-vector on device.

Parameters

global_index – Index of the target element.

inline void allocateBuffers()

Allocate send and recv buffers for MPI communication.

inline void barrier()
template<typename T>
inline T allReduceSum(const T &data) const

Perform an all-reduce operation with the sum operation.

Parameters

data – Data to be reduced

inline void sendrecvBuffers(const std::size_t send_rank, const std::size_t recv_rank, const std::size_t size, const std::size_t tag)

Perform a send-receive operation between two ranks using the sendbuf_ and recvbuf_.

Parameters
  • send_rank – Rank to send data to

  • recv_rank – Rank to receive data from

  • size – Size of the data to send/receive

  • tag – Tag for the MPI message

inline std::size_t getLocalBlockSize() const

Returns the MPI-distribution block size, or the size of the local state vector data.

inline std::size_t getRevLocalWireIndex(const std::size_t wire) const
inline std::size_t getRevGlobalWireIndex(const std::size_t wire) const
inline const std::vector<std::size_t> &getMPIRankToGlobalIndexMap() const
inline const std::vector<std::size_t> &getGlobalWires() const
inline const std::vector<std::size_t> &getLocalWires() const
inline std::vector<std::size_t> getLocalWireIndices(const std::vector<std::size_t> &wires) const
inline size_t getLocalWireIndex(const std::size_t wire) const
inline std::vector<std::size_t> getGlobalWiresIndices(const std::vector<std::size_t> &wires) const
inline size_t getGlobalWireIndex(const std::size_t wire) const
inline std::vector<std::size_t> findGlobalWires(const std::vector<std::size_t> &wires) const
inline std::vector<std::size_t> findLocalWires(const std::vector<std::size_t> &wires) const
inline std::pair<std::size_t, std::size_t> global2localIndex(const std::size_t index) const

Converts a global state vector index to a local one.

Parameters

index – Pair containing {global index/rank location, local index}

inline std::size_t getGlobalIndexFromMPIRank(const int mpi_rank) const
inline std::size_t getMPIRankFromGlobalIndex(const std::size_t global_index) const
inline bool isWiresLocal(const std::vector<std::size_t> &wires) const
inline bool isWiresGlobal(const std::vector<std::size_t> &wires) const
inline void normalize()

Normalize the state vector.

inline void resetStateVector()

Reset the data back to the \(\ket{0}\) state.

inline void setBasisState(const std::vector<std::size_t> &state, const std::vector<std::size_t> &wires)

Prepares a single computational basis state.

Parameters
  • state – Binary number representing the index

  • wires – Wires.

inline void setStateVector(const std::vector<ComplexT> &state, const std::vector<std::size_t> &wires)

Set values for a batch of elements of the state-vector.

Parameters
  • state – State.

  • wires – Wires.

inline void setStateVector(const ComplexT *state, const std::vector<std::size_t> &wires)

Set values for a batch of elements of the state-vector.

Parameters
  • state – State.

  • wires – Wires.

inline std::vector<std::size_t> localWiresSubsetToSwap(const std::vector<std::size_t> &global_wires, const std::vector<std::size_t> &wires)

Get the local wires that could be used for global/local swaps. It will return a set of local wires of the same size as the global wires to be swapped. The local wires will also not be in wires, i.e. wires used in an operation.

Parameters
  • global_wires – Global wires to swap, obtained from findGlobalWires()

  • wires – Wires used in an operation to be excluded in returned wires

inline void swapGlobalLocalWires(const std::vector<std::size_t> &global_wires_to_swap, const std::vector<std::size_t> &local_wires_to_swap)

Perform a swap between global and local wires.

Example: For Global|Local wires = 0 | 1 2 and swapping global_wires {0} and local_wires {2}, the elements {0|01, 0|11} in global_index 0 will be swap with elements {1|00, 1|10} in global_index 1. This will be done when batch_index = 1:

  • receiving_global_index = local_global_index ^ batch_index

  • the (relevant, i.e. to be swapped) local_wire = (relevant) local_wire ^ batch_index (wire 2 in example)

  • the (irrelevant, not swapped) local_wires are looped over and copied to/from the buffers to send/recv

Note: the batch index loops over all the global wires (including those that are not swapped). W e then check whether the specific batch number requires communication

For each batch, a single pairwise sendrecv is performed. The number of batches = 2^(num_swapping_wires) - 1 . The number of elements sent in each batch is 1/2^(num_swapping_wires) * size_of_subSV.

Parameters
  • global_wires_to_swap – wire indices for global wires to swap

  • local_wires_to_swap – wire indices for local wires to swap

inline void matchWires(const StateVectorKokkosMPI &other_sv)

Match the global/local wires and map from another state vector The following steps are performed (if necessary): Swap Global Local Wires (so that the first n-wires are global, and last wires a local, e.g. G={2,0,1}, L={3,5,4}) Swap Global Global Wires (make sure the global wires are in order, so move G to {0, 1, 2} Swap Local Local Wires (make sure the local wires are in order, so move L to {3, 4, 5}.

Parameters

other_sv – State vector to match

inline void matchLocalWires(const std::vector<std::size_t> &other_local_wires)
inline void matchGlobalWiresAndIndex(const StateVectorKokkosMPI &other_sv)
inline void matchGlobalWiresAndIndex(const std::vector<std::size_t> &global_wires_target, const std::vector<std::size_t> &mpi_rank_to_global_index_map_target)
inline void applyPauliRot(const std::vector<std::size_t> &wires, bool inverse, const std::vector<PrecisionT> &params, const std::string &word)

Apply a PauliRot gate to the state-vector.

Parameters
  • wires – Wires to apply gate to.

  • inverse – Indicates whether to use inverse of gate.

  • params – Rotation angle.

  • word – A Pauli word (e.g. “XYYX”).

inline void applyOperation(const std::string &opName, const std::vector<std::size_t> &wires, bool inverse = false, const std::vector<fp_t> &params = {}, const std::vector<ComplexT> &gate_matrix = {})

Apply a single gate to the state vector.

Parameters
  • opName – Name of gate to apply.

  • wires – Wires to apply gate to.

  • inverse – Indicates whether to use adjoint of gate.

  • params – Optional parameter list for parametric gates.

  • gate_matrix – Optional std gate matrix if opName doesn’t exist.

inline void applyOperation(const std::string &opName, const std::vector<std::size_t> &controlled_wires, const std::vector<bool> &controlled_values, const std::vector<std::size_t> &wires, bool inverse = false, const std::vector<fp_t> &params = {}, const std::vector<ComplexT> &gate_matrix = {})

Apply a controlled-single gate to the state vector.

Parameters
  • opName – Name of gate to apply.

  • controlled_wires – Control wires.

  • controlled_values – Control values (false or true).

  • wires – Wires to apply gate to.

  • inverse – Indicates whether to use adjoint of gate. (Default to false)

  • params – Optional parameter list for parametric gates.

  • gate_matrix – Optional unitary gate matrix if opName doesn’t exist.

inline void applyMatrix(const ComplexT *matrix, const std::vector<std::size_t> &wires, bool inverse = false)

Apply a given matrix directly to the statevector using a raw matrix pointer vector on host memory.

Parameters
  • matrix – Pointer to host matrix to apply to wires (in row-major format).

  • wires – Wires to apply gate to.

  • inverse – Indicate whether inverse should be taken. (Default to false)

inline void applyMatrix(const std::vector<ComplexT> &matrix, const std::vector<std::size_t> &wires, bool inverse = false)

Apply a given matrix as a vector directly to the statevector.

Parameters
  • matrix – Matrix data as a vector (in row-major format).

  • wires – Wires to apply gate to.

  • inverse – Indicate whether inverse should be taken. (Default to false)

inline void applyControlledMatrix(const ComplexT *matrix, const std::vector<std::size_t> &controlled_wires, const std::vector<bool> &controlled_values, const std::vector<std::size_t> &wires, bool inverse = false)

Apply a given matrix directly to the statevector using a raw matrix pointer on host memory.

Parameters
  • matrix – Pointer to the array data (in row-major format).

  • controlled_wires – Controlled wires

  • controlled_values – Controlled values (true or false)

  • wires – Wires to apply gate to.

  • inverse – Indicate whether inverse should be taken. (Default to false)

inline void applyControlledMatrix(const std::vector<ComplexT> &matrix, const std::vector<std::size_t> &controlled_wires, const std::vector<bool> &controlled_values, const std::vector<std::size_t> &wires, bool inverse = false)

Apply a given controlled-matrix as a vector directly to the statevector.

Parameters
  • matrix – Matrix data as a vector to apply to target wires (in row-major format).

  • controlled_wires – Control wires.

  • controlled_values – Control values (false or true).

  • wires – Wires to apply gate to.

  • inverse – Indicate whether inverse should be taken. (Default to false)

inline auto applyGenerator(const std::string &opName, const std::vector<std::size_t> &wires, bool inverse = false) -> PrecisionT

Apply a single generator to the state vector using the given kernel.

Parameters
  • opName – Name of gate to apply.

  • wires – Wires to apply gate to.

  • inverse – Indicates whether to use adjoint of gate. (Default to false)

Returns

PrecisionT Generator scale prefactor

inline auto applyControlledGenerator(const std::string &opName, const std::vector<std::size_t> &controlled_wires, const std::vector<bool> &controlled_values, const std::vector<std::size_t> &wires, bool inverse = false) -> PrecisionT

Apply a single controlled generator to the state vector using the given kernel.

Parameters
  • opName – Name of gate to apply.

  • controlled_wires – Control wires.

  • controlled_values – Control values (true or false).

  • wires – Wires to apply gate to.

  • inverse – Indicates whether to use adjoint of gate. (Default to false)

Returns

PrecisionT Generator scale prefactor

inline auto applyGenerator(const std::string &opName, const std::vector<std::size_t> &controlled_wires, const std::vector<bool> &controlled_values, const std::vector<std::size_t> &wires, bool inverse = false) -> PrecisionT

Apply a single generator to the state vector using the given kernel.

Parameters
  • opName – Name of gate to apply.

  • controlled_wires – Control wires.

  • controlled_values – Control values (true or false).

  • wires – Wires to apply gate to.

  • inverse – Indicates whether to use adjoint of gate. (Default to false)

Returns

PrecisionT Generator scale prefactor

inline void collapse(std::size_t wire, bool branch)

Collapse the state vector after having measured one of the qubits.

The branch parameter imposes the measurement result on the given wire.

Parameters
  • wire – Wire to collapse.

  • branch – Branch 0 or 1.

inline void updateData(const KokkosVector other)

Update data of the class.

Parameters

other – Kokkos View

inline void updateData(const StateVectorKokkos<fp_t> &other)

Update data of the class.

Parameters

other – State vector

inline void updateData(ComplexT *new_data, std::size_t new_size)

Update data of the class.

Parameters
  • new_data – data pointer to new data.

  • new_size – size of underlying data storage.

inline void updateData(std::vector<ComplexT> &other)

Update data of the class.

Parameters

other – STL vector of type ComplexT

inline void updateData(const StateVectorKokkosMPI<PrecisionT> &other)
inline auto getData() -> ComplexT*

Get underlying Kokkos view data on the device.

Returns

ComplexT *

inline auto getData() const -> const ComplexT*
inline auto getView() const -> KokkosVector&

Get the Kokkos data of the state vector.

Returns

The pointer to the data of state vector

inline auto getView() -> KokkosVector&

Get the Kokkos data of the state vector.

Returns

The pointer to the data of state vector

inline auto getDataVector() -> std::vector<ComplexT>

Get the local vector-converted Kokkos view.

Returns

std::vector<ComplexT>

inline auto getDataVector() const -> const std::vector<ComplexT>
inline void HostToDevice(ComplexT *sv, std::size_t length)

Copy data from the host space to the device space.

inline void DeviceToHost(ComplexT *sv, std::size_t length) const

Copy data from the device space to the host space.

inline void DeviceToDevice(KokkosVector vector_to_copy)

Copy data from the device space to the device space.

inline void reorderLocalWires()

Make sure local wires are in ascending order Must be used after reorderGlobalLocalWires()

inline void reorderGlobalLocalWires()

Make sure global wires are lowest numbered wires and in ascending order.

inline void reorderAllWires()