Class SimulatorMPI

class SimulatorMPI

Public Functions

SimulatorMPI(uint64_t seed, size_t max_local, size_t max_cluster_size)

Constructor.

Parameters
  • seed: Seed for pseudo-random number generator

  • max_local: Maximum number of local qubits

  • max_cluster_size: Maximum number of qubits in fused multi-qubit gate

SimulatorMPI(mpi::communicator aWorld, uint64_t seed, size_t max_local, size_t max_cluster_size)

Constructor.

With this constructor you can define your own MPI communicator which must be able to work with 2^N processes, where N is the number of qubits.

Parameters
  • aWorld: A MPI communicator

  • seed: Seed for pseudo-random number generator

  • max_local: Maximum number of local qubits

  • max_cluster_size: Maximum number of qubits in fused multi-qubit gate

SimulatorMPI(const SimulatorMPI &other)

Copy constructor.

Note

Actual declaration is SimulatorMPI(const SimulatorMPI &other) = delete.

Copying of the simulator is prohibited.

~SimulatorMPI()

Destructor.

void AllocateQubit(Index id)

Allocate one qubit with a given ID.

Parameters
  • id: ID of the qubit to allocate

void AllocateQureg(const std::vector<Index> &ids, Complex init = 0)

Allocate the quantum register qureg with given qubits IDs.

Parameters
  • ids: Array of qubit IDs

  • init: Value to be assigned to all amplitudes

Exceptions
  • std::runtime_error: if this is not the first allocation and init is not equal 0.0

void DeallocateQubit(Index id)

De-allocate one qubit with a given ID.

Parameters
  • id: ID of the qubit to de-allocate

void Run()

Actually make calculations:

  • fuse applied (used ApplyGate()) gates into a single multi-qubit gate

  • calculate the product of the resulting gate by state vector

    Exceptions
    • std::runtime_error: if gate cannot be applied

void ApplyGate(Matrix m, std::vector<Index> ids, std::vector<Index> ctrl)

Save gate and than apply it to function Run()

Parameters
  • m: Matrix of the gate

  • ids: Array of qubit IDs

  • ctrl: Array of control qubits

Exceptions
  • std::runtime_error: if gate is non-diagonal

std::vector<bool> MeasureQubits(std::vector<Index> const &ids)

Measure qubits state.

Return

Bit array of measurement results (containing either True or False)

Parameters
  • ids: Array of qubit IDs to measure

void SwapQubitsWrapper(const std::vector<Index> &swap_pairs_ids)

See SwapQubits()

Print logs: duration, qubits, bandwidth

void SwapQubits(const std::vector<Index> &swap_pairs_ids)

Swap global and local qubits in pairs.

Parameters
  • swap_pairs_ids: Array of qubit pairs IDs (there is a global qubit ID on the first place, a local - on the second)

Exceptions
  • std::runtime_error: if qubits are not unique

SimulatorMPI::Complex GetAmplitude(const std::vector<bool> &bit_string, const std::vector<Index> &ids) const

Return the amplitude of the supplied amplitude index.

Return

Amplitude of the provided bit_string

Parameters
  • bit_string: Index amplitude

  • ids: Array of qubit IDs

Exceptions
  • std::runtime_error: if the number of IDs does not match the number of qubits or if any of the desired qubits are not already allocated

SimulatorMPI::Float GetProbability(const std::vector<bool> &bit_string, const std::vector<Index> &ids)

Return the probability of the outcome bit_string when measuring the quantum register qureg.

Return

Probability of measuring the provided bit string

Note

Doesn’t actually make the calculations

Parameters
  • bit_string: Measurement outcome

  • ids: Array of qubit IDs

Exceptions
  • std::runtime_error: if the number of IDs does not match the length of bit_string

SimulatorMPI::Float Entropy()

Return the entropy of quantum state.

std::vector<SimulatorMPI::Index> GetQubitsPermutation() const

Get permutation of qubits.

Return

Array of qubit IDs

std::vector<SimulatorMPI::Index> GetLocalQubitsPermutation()

Get permutation of local qubits.

Return

Array of qubit IDs

std::vector<SimulatorMPI::Index> GetGlobalQubitsPermutation()

Get permutation of global qubits.

Return

Array of qubit IDs

void SetQubitsPermutation(std::vector<Index> p)

Set permutation to all qubits without actually reordering quantum state vector.

Parameters
  • p: Array of qubit IDs

std::tuple<std::map<int, int>, SimulatorMPI::StateVector&> cheat_local()

Access the ordering of the qubits and this MPI process’s part of state vector directly.

Return

A tuple where the first entry is a dictionary mapping qubit indices to bit-locations (all qubits) and the second entry is the local part of state vector.

void collapseWaveFunction(const std::vector<Index> &ids, const std::vector<bool> &values)

Collapse a quantum register into a classical basis state.

Parameters
  • ids: Array of qubit IDs

  • values: Measurement outcome for each of the qubits in qureg.

Exceptions
  • std::runtime_error: if the number of IDs does not match the number of values or measured probability less than 1.e-12

template<class F, class QuReg>
void emulate_math(F const&, QuReg, const std::vector<Index>&, unsigned = 1)

This function is not supported.

size_t LocalQubitsCount() const

Return

Number of local qubits

size_t GlobalQubitsCount() const

Return

Number of global qubits

size_t TotalQubitsCount() const

Return

Number of all qubits