Template Class StateVectorCudaManaged

Inheritance Relationships

Base Type

Class Documentation

template<class Precision = double>
class StateVectorCudaManaged : public Pennylane::LightningGPU::StateVectorCudaBase<double, StateVectorCudaManaged<double>>

Managed memory CUDA state-vector class using custateVec backed gate-calls.

Template Parameters

Precision – Floating-point precision type.

Public Types

using PrecisionT = Precision
using ComplexT = std::complex<PrecisionT>
using CFP_t = typename StateVectorCudaBase<Precision, StateVectorCudaManaged<Precision>>::CFP_t
using MemoryStorageT = Pennylane::Util::MemoryStorageLocation::Undefined

Public Functions

StateVectorCudaManaged() = delete
inline StateVectorCudaManaged(size_t num_qubits)
inline StateVectorCudaManaged(size_t num_qubits, const DevTag<int> &dev_tag, bool alloc = true, SharedCusvHandle cusvhandle_in = make_shared_cusv_handle(), SharedCublasCaller cublascaller_in = make_shared_cublas_caller(), SharedCusparseHandle cusparsehandle_in = make_shared_cusparse_handle())
inline StateVectorCudaManaged(const CFP_t *gpu_data, size_t length)
inline StateVectorCudaManaged(const CFP_t *gpu_data, size_t length, DevTag<int> dev_tag, SharedCusvHandle handle_in = make_shared_cusv_handle(), SharedCublasCaller cublascaller_in = make_shared_cublas_caller(), SharedCusparseHandle cusparsehandle_in = make_shared_cusparse_handle())
inline StateVectorCudaManaged(const std::complex<Precision> *host_data, size_t length)
inline StateVectorCudaManaged(std::complex<Precision> *host_data, size_t length)
inline StateVectorCudaManaged(const StateVectorCudaManaged &other)
~StateVectorCudaManaged() = default
inline void setBasisState(const std::complex<Precision> &value, const size_t index, const bool async = false)

Set value for a single element of the state-vector on device. This method is implemented by cudaMemcpy.

Parameters
  • value – Value to be set for the target element.

  • index – Index of the target element.

  • async – Use an asynchronous memory copy.

template<class index_type, size_t thread_per_block = 256>
inline void setStateVector(const index_type num_indices, const std::complex<Precision> *values, const index_type *indices, const bool async = false)

Set values for a batch of elements of the state-vector. This method is implemented by the customized CUDA kernel defined in the DataBuffer class.

Parameters
  • num_indices – Number of elements to be passed to the state vector.

  • values – Pointer to values to be set for the target elements.

  • indices – Pointer to indices of the target elements.

  • async – Use an asynchronous memory copy.

template<size_t thread_per_block = 256>
inline void globalPhaseStateVector(const bool adjoint, const Precision param)

Multiplies the state-vector by a global phase.

Parameters
  • adjoint – Indicates whether to use adjoint of gate.

  • param – Complex phase generator.

template<size_t thread_per_block = 256>
inline void cGlobalPhaseStateVector(const bool adjoint, const std::vector<CFP_t> &phase, const bool async = false)

Multiplies the state-vector by a controlled global phase.

Parameters

phase – Controlled complex phase vector.

inline void applyOperation(const std::string &opName, const std::vector<size_t> &wires, bool adjoint, const std::vector<Precision> &params, const std::vector<ComplexT> &matrix)

Apply a single gate to the state-vector. Offloads to custatevec specific API calls if available. If unable, attempts to use prior cached gate values on the device. Lastly, accepts a host-provided matrix if otherwise, and caches on the device for later reuse.

Parameters
  • opName – Name of gate to apply.

  • wires – Wires to apply gate to.

  • adjoint – Indicates whether to use adjoint of gate.

  • params – Optional parameter list for parametric gates.

  • matrix – Gate data (in row-major format).

inline void applyOperation(const std::string &opName, const std::vector<size_t> &wires, bool adjoint = false, const std::vector<Precision> &params = {0.0}, const std::vector<CFP_t> &gate_matrix = {})

Apply a single gate to the state-vector. Offloads to custatevec specific API calls if available. If unable, attempts to use prior cached gate values on the device. Lastly, accepts a host-provided matrix if otherwise, and caches on the device for later reuse.

Parameters
  • opName – Name of gate to apply.

  • wires – Wires to apply gate to.

  • adjoint – Indicates whether to use adjoint of gate.

  • params – Optional parameter list for parametric gates.

  • gate_matrix – Gate data (in row-major format).

template<template<typename...> class complex_t>
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<Precision> &params = {0.0}, const std::vector<complex_t<Precision>> &gate_matrix = {})

Apply a 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.

  • params – Optional parameter list for parametric gates.

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

inline auto applyGenerator(const std::string &opName, const std::vector<size_t> &wires, bool adjoint = 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.

  • adjoint – Indicates whether to use adjoint of gate.

inline void applyMatrix(const std::complex<PrecisionT> *gate_matrix, const std::vector<size_t> &wires, bool adjoint = false)

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

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

  • wires – Wires to apply gate to.

  • adjoint – Indicate whether inverse should be taken.

inline void applyMatrix(const std::vector<std::complex<PrecisionT>> &gate_matrix, const std::vector<size_t> &wires, bool adjoint = false)

Apply a given matrix directly to the statevector using a std vector.

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

  • wires – Wires to apply gate to.

  • adjoint – Indicate whether inverse should be taken.

inline void applyIdentity(const std::vector<std::size_t> &wires, bool adjoint)
inline void applyPauliX(const std::vector<std::size_t> &wires, bool adjoint)
inline void applyPauliY(const std::vector<std::size_t> &wires, bool adjoint)
inline void applyPauliZ(const std::vector<std::size_t> &wires, bool adjoint)
inline void applyHadamard(const std::vector<std::size_t> &wires, bool adjoint)
inline void applyS(const std::vector<std::size_t> &wires, bool adjoint)
inline void applyT(const std::vector<std::size_t> &wires, bool adjoint)
inline void applyRX(const std::vector<std::size_t> &wires, bool adjoint, Precision param)
inline void applyRY(const std::vector<std::size_t> &wires, bool adjoint, Precision param)
inline void applyRZ(const std::vector<std::size_t> &wires, bool adjoint, Precision param)
inline void applyRot(const std::vector<std::size_t> &wires, bool adjoint, Precision param0, Precision param1, Precision param2)
inline void applyRot(const std::vector<std::size_t> &wires, bool adjoint, const std::vector<Precision> &params)
inline void applyPhaseShift(const std::vector<std::size_t> &wires, bool adjoint, Precision param)
inline void applyCNOT(const std::vector<std::size_t> &wires, bool adjoint)
inline void applyCY(const std::vector<std::size_t> &wires, bool adjoint)
inline void applyCZ(const std::vector<std::size_t> &wires, bool adjoint)
inline void applySWAP(const std::vector<std::size_t> &wires, bool adjoint)
inline void applyIsingXX(const std::vector<std::size_t> &wires, bool adjoint, Precision param)
inline void applyIsingYY(const std::vector<std::size_t> &wires, bool adjoint, Precision param)
inline void applyIsingZZ(const std::vector<std::size_t> &wires, bool adjoint, Precision param)
inline void applyIsingXY(const std::vector<std::size_t> &wires, bool adjoint, Precision param)
inline void applyCRot(const std::vector<std::size_t> &wires, bool adjoint, const std::vector<Precision> &params)
inline void applyCRot(const std::vector<std::size_t> &wires, bool adjoint, Precision param0, Precision param1, Precision param2)
inline void applyCRX(const std::vector<std::size_t> &wires, bool adjoint, Precision param)
inline void applyCRY(const std::vector<std::size_t> &wires, bool adjoint, Precision param)
inline void applyCRZ(const std::vector<std::size_t> &wires, bool adjoint, Precision param)
inline void applyControlledPhaseShift(const std::vector<std::size_t> &wires, bool adjoint, Precision param)
inline void applySingleExcitation(const std::vector<std::size_t> &wires, bool adjoint, Precision param)
inline void applySingleExcitationMinus(const std::vector<std::size_t> &wires, bool adjoint, Precision param)
inline void applySingleExcitationPlus(const std::vector<std::size_t> &wires, bool adjoint, Precision param)
inline void applyToffoli(const std::vector<std::size_t> &wires, bool adjoint)
inline void applyCSWAP(const std::vector<std::size_t> &wires, bool adjoint)
inline void applyDoubleExcitation(const std::vector<std::size_t> &wires, bool adjoint, Precision param)
inline void applyDoubleExcitationMinus(const std::vector<std::size_t> &wires, bool adjoint, Precision param)
inline void applyDoubleExcitationPlus(const std::vector<std::size_t> &wires, bool adjoint, Precision param)
inline void applyMultiRZ(const std::vector<std::size_t> &wires, bool adjoint, Precision param)
inline PrecisionT applyGeneratorGlobalPhase([[maybe_unused]] const std::vector<size_t> &wires, [[maybe_unused]] bool adj = false)

Gradient generator function associated with the GlobalPhase gate.

Parameters
  • sv – Statevector

  • wires – Wires to apply operation.

  • adj – Takes adjoint of operation if true. Defaults to false.

inline PrecisionT applyGeneratorRX(const std::vector<size_t> &wires, bool adj = false)

Gradient generator function associated with the RX gate.

Parameters
  • sv – Statevector

  • wires – Wires to apply operation.

  • adj – Takes adjoint of operation if true. Defaults to false.

inline PrecisionT applyGeneratorRY(const std::vector<size_t> &wires, bool adj = false)

Gradient generator function associated with the RY gate.

Parameters
  • sv – Statevector

  • wires – Wires to apply operation.

  • adj – Takes adjoint of operation if true. Defaults to false.

inline PrecisionT applyGeneratorRZ(const std::vector<size_t> &wires, bool adj = false)

Gradient generator function associated with the RZ gate.

Parameters
  • sv – Statevector

  • wires – Wires to apply operation.

  • adj – Takes adjoint of operation if true. Defaults to false.

inline PrecisionT applyGeneratorIsingXX(const std::vector<std::size_t> &wires, bool adjoint)
inline PrecisionT applyGeneratorIsingYY(const std::vector<std::size_t> &wires, bool adjoint)
inline PrecisionT applyGeneratorIsingZZ(const std::vector<std::size_t> &wires, bool adjoint)
inline PrecisionT applyGeneratorIsingXY(const std::vector<std::size_t> &wires, bool adjoint)
inline PrecisionT applyGeneratorPhaseShift(const std::vector<size_t> &wires, bool adj = false)

Gradient generator function associated with the PhaseShift gate.

Parameters
  • wires – Wires to apply operation.

  • adj – Takes adjoint of operation if true. Defaults to false.

inline PrecisionT applyGeneratorCRX(const std::vector<size_t> &wires, bool adj = false)

Gradient generator function associated with the controlled RX gate.

Parameters
  • wires – Wires to apply operation.

  • adj – Takes adjoint of operation if true. Defaults to false.

inline PrecisionT applyGeneratorCRY(const std::vector<size_t> &wires, bool adj = false)

Gradient generator function associated with the controlled RY gate.

Parameters
  • wires – Wires to apply operation.

  • adj – Takes adjoint of operation if true. Defaults to false.

inline PrecisionT applyGeneratorCRZ(const std::vector<size_t> &wires, bool adj = false)

Gradient generator function associated with the controlled RZ gate.

Parameters
  • wires – Wires to apply operation.

  • adj – Takes adjoint of operation if true. Defaults to false.

inline PrecisionT applyGeneratorControlledPhaseShift(const std::vector<size_t> &wires, bool adj = false)

Gradient generator function associated with the controlled PhaseShift gate.

Parameters
  • wires – Wires to apply operation.

  • adj – Takes adjoint of operation if true. Defaults to false.

inline PrecisionT applyGeneratorSingleExcitation(const std::vector<std::size_t> &wires, bool adjoint)
inline PrecisionT applyGeneratorSingleExcitationMinus(const std::vector<std::size_t> &wires, bool adjoint)
inline PrecisionT applyGeneratorSingleExcitationPlus(const std::vector<std::size_t> &wires, bool adjoint)
inline PrecisionT applyGeneratorDoubleExcitation(const std::vector<std::size_t> &wires, bool adjoint)
inline PrecisionT applyGeneratorDoubleExcitationMinus(const std::vector<std::size_t> &wires, bool adjoint)
inline PrecisionT applyGeneratorDoubleExcitationPlus(const std::vector<std::size_t> &wires, bool adjoint)
inline PrecisionT applyGeneratorMultiRZ(const std::vector<std::size_t> &wires, bool adjoint)
inline auto getCublasCaller() const -> const CublasCaller&

Access the CublasCaller the object is using.

Returns

a reference to the object’s CublasCaller object.

inline auto getCusparseHandle() const -> cusparseHandle_t

Get the cuSPARSE handle that the object is using.

Returns

cusparseHandle_t returns the cuSPARSE handle.

inline auto getCusvHandle() const -> custatevecHandle_t

Get the cuStateVec handle that the object is using.

Returns

custatevecHandle_t returns the cuStateVec handle.

inline auto getDataVector() -> std::vector<std::complex<PrecisionT>>

Get a host data copy.

Returns

std::vector<std::complex<PrecisionT>>