Maestro 0.2.5
Unified interface for quantum circuit simulation
Loading...
Searching...
No Matches
QCSimState.h
Go to the documentation of this file.
1
13#pragma once
14
15#ifndef _QCSIMSTATE_H_
16#define _QCSIMSTATE_H_
17
18#ifdef INCLUDED_BY_FACTORY
19
20#include <algorithm>
21
22#include "Simulator.h"
23
24#include "Clifford.h"
25#include "MPSSimulator.h"
26#include "QubitRegister.h"
28
29#include "../TensorNetworks/ForestContractor.h"
30#include "../TensorNetworks/TensorNetwork.h"
31
32#include "../Utils/Alias.h"
33
34#include "MPSDummySimulator.h"
35
36namespace Simulators {
37// TODO: Maybe use the pimpl idiom
38// https://en.cppreference.com/w/cpp/language/pimpl to hide the implementation
39// for good but during development this should be good enough
40namespace Private {
41
52class QCSimState : public ISimulator {
53 public:
54 QCSimState() : rng(std::random_device{}()), uniformZeroOne(0, 1) {}
55
63 void Initialize() override {
64 if (nrQubits != 0) {
65 if (simulationType == SimulationType::kMatrixProductState) {
66 mpsSimulator =
67 std::make_unique<QC::TensorNetworks::MPSSimulator>(nrQubits);
68 if (limitEntanglement && singularValueThreshold > 0.)
69 mpsSimulator->setLimitEntanglement(singularValueThreshold);
70 if (limitSize && chi > 0) mpsSimulator->setLimitBondDimension(chi);
71 } else if (simulationType == SimulationType::kStabilizer)
72 cliffordSimulator =
73 std::make_unique<QC::Clifford::StabilizerSimulator>(nrQubits);
74 else if (simulationType == SimulationType::kTensorNetwork) {
75 tensorNetwork =
76 std::make_unique<TensorNetworks::TensorNetwork>(nrQubits);
77 // for now the only used contractor is the forest one, but we'll use
78 // more in the future
79 const auto tensorContractor =
80 std::make_shared<TensorNetworks::ForestContractor>();
81 tensorNetwork->SetContractor(tensorContractor);
82 } else if (simulationType == SimulationType::kPauliPropagator) {
83 pp = std::make_unique<Simulators::QcsimPauliPropagator>();
84 pp->SetNrQubits(static_cast<int>(nrQubits));
85 } else
86 state = std::make_unique<QC::QubitRegister<>>(nrQubits);
87
88 SetMultithreading(enableMultithreading);
89 }
90 }
91
103 void InitializeState(size_t num_qubits,
104 std::vector<std::complex<double>> &amplitudes) override {
105 if (num_qubits == 0) return;
106 Clear();
107 nrQubits = num_qubits;
108 Initialize();
109 if (simulationType != SimulationType::kStatevector)
110 throw std::runtime_error(
111 "QCSimState::InitializeState: Invalid "
112 "simulation type for initializing the state.");
113
114 Eigen::VectorXcd amplitudesEigen(
115 Eigen::Map<Eigen::VectorXcd, Eigen::Unaligned>(amplitudes.data(),
116 amplitudes.size()));
117 state->setRegisterStorageFastNoNormalize(amplitudesEigen);
118 }
119
131 /*
132 void InitializeState(size_t num_qubits, std::vector<std::complex<double>,
133 avoid_init_allocator<std::complex<double>>>& amplitudes) override
134 {
135 Clear();
136 nrQubits = num_qubits;
137 Initialize();
138 Eigen::VectorXcd amplitudesEigen(Eigen::Map<Eigen::VectorXcd,
139 Eigen::Unaligned>(amplitudes.data(), amplitudes.size()));
140 state->setRegisterStorageFastNoNormalize(amplitudesEigen);
141 }
142 */
143
155#ifndef NO_QISKIT_AER
156 void InitializeState(size_t num_qubits,
157 AER::Vector<std::complex<double>> &amplitudes) override {
158 if (num_qubits == 0) return;
159 Clear();
160 nrQubits = num_qubits;
161 Initialize();
162 if (simulationType != SimulationType::kStatevector)
163 throw std::runtime_error(
164 "QCSimState::InitializeState: Invalid "
165 "simulation type for initializing the state.");
166
167 Eigen::VectorXcd amplitudesEigen(
168 Eigen::Map<Eigen::VectorXcd, Eigen::Unaligned>(amplitudes.data(),
169 amplitudes.size()));
170 state->setRegisterStorageFastNoNormalize(amplitudesEigen);
171 }
172#endif
173
185 void InitializeState(size_t num_qubits,
186 Eigen::VectorXcd &amplitudes) override {
187 if (num_qubits == 0) return;
188 Clear();
189 nrQubits = num_qubits;
190 Initialize();
191
192 if (simulationType != SimulationType::kStatevector)
193 throw std::runtime_error(
194 "QCSimState::InitializeState: Invalid "
195 "simulation type for initializing the state.");
196
197 state = std::make_unique<QC::QubitRegister<>>(nrQubits, amplitudes);
198 state->SetMultithreading(enableMultithreading);
199 }
200
207 void Reset() override {
208 if (mpsSimulator)
209 mpsSimulator->Clear();
210 else if (cliffordSimulator)
211 cliffordSimulator->Reset();
212 else if (tensorNetwork)
213 tensorNetwork->Clear();
214 else if (state)
215 state->Reset();
216 else if (pp)
217 pp->ClearOperations();
218 }
219
227 bool SupportsMPSSwapOptimization() const override { return true; }
228
237 void SetInitialQubitsMap(
238 const std::vector<long long int> &initialMap) override {
239 if (mpsSimulator) {
240 mpsSimulator->SetInitialQubitsMap(initialMap);
241 if (!dummySim || dummySim->getNrQubits() != initialMap.size()) {
242 dummySim = std::make_unique<Simulators::MPSDummySimulator>(initialMap.size());
243 dummySim->SetMaxBondDimension(
244 limitSize ? static_cast<long long int>(chi) : 0);
245 }
246 dummySim->SetInitialQubitsMap(initialMap);
247 }
248 }
249
250 void SetUseOptimalMeetingPosition(bool enable) override {
251 useOptimalMeetingPositionOnly = enable;
252 if (mpsSimulator)
253 mpsSimulator->SetUseOptimalMeetingPosition(enable);
254 }
255
256 void SetLookaheadDepth(int depth) override {
257 lookaheadDepth = depth;
258 if (mpsSimulator && depth > 0 && !useOptimalMeetingPositionOnly)
259 mpsSimulator->SetUseOptimalMeetingPosition(true);
260 }
261
262 void SetLookaheadDepthWithHeuristic(int depth) override
263 {
264 lookaheadDepthWithHeuristic = depth;
265 if (lookaheadDepth < depth) SetLookaheadDepth(depth);
266 }
267
268 void SetUpcomingGates(
269 const std::vector<std::shared_ptr<Circuits::IOperation<double>>> &gates) override {
270 upcomingGates = gates;
271 upcomingGateIndex = 0;
272
273 if (!mpsSimulator || lookaheadDepth <= 0) return;
274
275 // Register an observer that advances the gate index
276 ClearObservers(); // for now we only have this observer, so this should be fine
277 gateCounterObserver = std::make_shared<GateCounterObserver>(upcomingGateIndex);
278 RegisterObserver(gateCounterObserver);
279
280 // Set up a meeting position callback that uses MPSDummySimulator
281 // for lookahead evaluation with actual bond dimensions
282 // the callback is called only for two qubits gates and only if executing them would require a swap
283 mpsSimulator->SetMeetingPositionCallback(
284 [this](/*const auto &qMap,*/ const auto &bondDims)
285 -> QC::TensorNetworks::MPSSimulatorInterface::IndexType {
286 const size_t nQ = bondDims.size() + 1;
287
288 if (!dummySim || dummySim->getNrQubits() != nQ) {
289 dummySim = std::make_unique<Simulators::MPSDummySimulator>(nQ);
290 dummySim->SetMaxBondDimension(limitSize ? static_cast<long long int>(chi) : 0);
291 }
292
293 // Seed dummy with current real simulator state
294 //std::vector<long long int> map64(qMap.begin(), qMap.end());
295 //dummySim->SetInitialQubitsMap(map64);
296 dummySim->setTotalSwappingCost(0);
297
298
299 // check qubits map:
300 /*
301 auto qbitmMap = dummySim->getQubitsMap();
302 for (size_t i = 0; i < nQ; ++i) {
303 if (qbitmMap[i] != qMap[i]) {
304 std::cerr << "Error: qubits map mismatch at index " << i
305 << ": dummySim has " << qbitmMap[i]
306 << " but real sim has " << qMap[i] << std::endl;
307 exit(0);
308 }
309 }
310 */
311
312 // check them, they should be the same, otherwise something is wrong
313
314
315 // Convert actual bond dims to doubles
316 std::vector<double> bondDimsD(bondDims.begin(), bondDims.end());
317 dummySim->SetCurrentBondDimensions(bondDimsD);
318
319 // display bond dimensions for debugging
320 /*
321 std::cout << "Bond dimensions before swapping and applying the gate: ";
322 for (size_t i = 0; i < bondDims.size(); ++i) {
323 std::cout << bondDims[i] << " ";
324 }
325 std::cout << std::endl;
326 */
327
328 const auto &op = upcomingGates[upcomingGateIndex];
329 const auto qbits = op->AffectedQubits();
330
331 /*
332 const auto &qmap = dummySim->getQubitsMap();
333
334 std::cout << "Applying 2-qubit gate on physical qubits " << qmap[qbits[0]] << " and "
335 << qmap[qbits[1]]
336 << std::endl;
337 */
338 /*
339 std::cout << "Finding best meeting position for upcoming gates starting at index "
340 << upcomingGateIndex << " with lookahead depth " << lookaheadDepth << " and heuristic depth "
341 << lookaheadDepthWithHeuristic << std::endl;
342
343
344
345 std::cout << "Affected qubits: ";
346 for (const auto &q : qbits) std::cout << q << " ";
347 std::cout << std::endl;
348
349 std::cout << "Current qubits map: ";
350 for (size_t i = 0; i < qMap.size(); ++i) std::cout << qMap[i] << " ";
351 std::cout << std::endl;
352
353 std::cout << "Current inverse qubits map: ";
354 for (size_t i = 0; i < qMapInv.size(); ++i) std::cout << qMapInv[i] << " ";
355 std::cout << std::endl;
356 */
357
358 double bestCost = std::numeric_limits<double>::infinity();
359 auto res = dummySim->FindBestMeetingPosition(upcomingGates, upcomingGateIndex, lookaheadDepth, lookaheadDepthWithHeuristic, 0, bestCost);
360
361 //std::cout << "Swapping the two qubits on position: " << res << " and " << (res + 1) << std::endl;
362
363 dummySim->SwapQubitsToPosition(qbits[0], qbits[1], res);
364 dummySim->ApplyGate(op);
365
366 // display the expected bond dimensions after applying the gate for
367 // debugging
368
369 /*
370 const auto &expectedBondDims = dummySim->getCurrentBondDimensions();
371 std::cout << "Expected bond dimensions after swapping and applying "
372 "the gate: ";
373 for (size_t i = 0; i < expectedBondDims.size(); ++i) {
374 std::cout << expectedBondDims[i] << " ";
375 }
376 std::cout << std::endl;
377 */
378
379 //std::cout << "Best meeting position: " << res
380 // << " with estimated cost: " << bestCost << std::endl;
381
382 return res;
383 });
384 }
385
394 long long int GetGatesCounter() const override { return upcomingGateIndex; }
395
405 void SetGatesCounter(long long int counter) override {
406 upcomingGateIndex = counter;
407 }
408
417 void IncrementGatesCounter() override {
418 ++upcomingGateIndex;
419 }
420
429 void Configure(const char *key, const char *value) override {
430 if (std::string("method") == key) {
431 if (std::string("statevector") == value)
432 simulationType = SimulationType::kStatevector;
433 else if (std::string("matrix_product_state") == value)
434 simulationType = SimulationType::kMatrixProductState;
435 else if (std::string("stabilizer") == value)
436 simulationType = SimulationType::kStabilizer;
437 else if (std::string("tensor_network") == value)
438 simulationType = SimulationType::kTensorNetwork;
439 else if (std::string("pauli_propagator") == value)
440 simulationType = SimulationType::kPauliPropagator;
441 } else if (std::string("matrix_product_state_truncation_threshold") ==
442 key) {
443 singularValueThreshold = std::stod(value);
444 if (singularValueThreshold > 0.) {
445 limitEntanglement = true;
446 if (mpsSimulator)
447 mpsSimulator->setLimitEntanglement(singularValueThreshold);
448 } else
449 limitEntanglement = false;
450 } else if (std::string("matrix_product_state_max_bond_dimension") == key) {
451 chi = std::stoi(value);
452 if (chi > 0) {
453 limitSize = true;
454 if (mpsSimulator) mpsSimulator->setLimitBondDimension(chi);
455 if (dummySim) dummySim->SetMaxBondDimension(static_cast<long long int>(chi));
456 } else {
457 limitSize = false;
458 if (mpsSimulator) mpsSimulator->setLimitBondDimension(0);
459 if (dummySim) dummySim->SetMaxBondDimension(0);
460 }
461 } else if (std::string("mps_sample_measure_algorithm") == key)
462 useMPSMeasureNoCollapse = std::string("mps_probabilities") == value;
463 }
464
472 std::string GetConfiguration(const char *key) const override {
473 if (std::string("method") == key) {
474 switch (simulationType) {
475 case SimulationType::kStatevector:
476 return "statevector";
477 case SimulationType::kMatrixProductState:
478 return "matrix_product_state";
479 case SimulationType::kStabilizer:
480 return "stabilizer";
481 case SimulationType::kTensorNetwork:
482 return "tensor_network";
483 case SimulationType::kPauliPropagator:
484 return "pauli_propagator";
485 default:
486 return "other";
487 }
488 } else if (std::string("matrix_product_state_truncation_threshold") ==
489 key) {
490 if (limitEntanglement && singularValueThreshold > 0.)
491 return std::to_string(singularValueThreshold);
492 } else if (std::string("matrix_product_state_max_bond_dimension") == key) {
493 if (limitSize && chi > 0) return std::to_string(chi);
494 } else if (std::string("mps_sample_measure_algorithm") == key) {
495 return useMPSMeasureNoCollapse ? "mps_probabilities"
496 : "mps_apply_measure";
497 }
498
499 return "";
500 }
501
509 size_t AllocateQubits(size_t num_qubits) override {
510 if ((simulationType == SimulationType::kStatevector && state) ||
511 (simulationType == SimulationType::kMatrixProductState &&
512 mpsSimulator) ||
513 (simulationType == SimulationType::kStabilizer && cliffordSimulator) ||
514 (simulationType == SimulationType::kTensorNetwork && tensorNetwork))
515 return 0;
516
517 const size_t oldNrQubits = nrQubits;
518 nrQubits += num_qubits;
519 if (simulationType == SimulationType::kPauliPropagator)
520 if (pp) pp->SetNrQubits(static_cast<int>(nrQubits));
521
522 return oldNrQubits;
523 }
524
531 size_t GetNumberOfQubits() const override { return nrQubits; }
532
540 void Clear() override {
541 state = nullptr;
542 mpsSimulator = nullptr;
543 cliffordSimulator = nullptr;
544 tensorNetwork = nullptr;
545 pp = nullptr;
546 nrQubits = 0;
547 }
548
559 size_t Measure(const Types::qubits_vector &qubits) override {
560 // TODO: this is inefficient, maybe implement it better in qcsim
561 // for now it has the possibility of measuring a qubits interval, but not a
562 // list of qubits
563 if (qubits.size() > sizeof(size_t) * 8)
564 std::cerr
565 << "Warning: The number of qubits to measure is larger than the "
566 "number of bits in the size_t type, the outcome will be undefined"
567 << std::endl;
568
569 size_t res = 0;
570 size_t mask = 1ULL;
571
572 DontNotify();
573 if (simulationType == SimulationType::kStatevector) {
574 for (size_t qubit : qubits) {
575 if (state->MeasureQubit(static_cast<unsigned int>(qubit))) res |= mask;
576 mask <<= 1;
577 }
578 } else if (simulationType == SimulationType::kStabilizer) {
579 for (size_t qubit : qubits) {
580 if (cliffordSimulator->MeasureQubit(static_cast<unsigned int>(qubit)))
581 res |= mask;
582 mask <<= 1;
583 }
584 } else if (simulationType == SimulationType::kTensorNetwork) {
585 for (size_t qubit : qubits) {
586 if (tensorNetwork->Measure(static_cast<unsigned int>(qubit)))
587 res |= mask;
588 mask <<= 1;
589 }
590 } else if (simulationType == SimulationType::kPauliPropagator) {
591 std::vector<int> qubitsInt(qubits.begin(), qubits.end());
592 const auto res = pp->Measure(qubitsInt);
593 Types::qubit_t result = 0;
594 for (size_t i = 0; i < res.size(); ++i) {
595 if (res[i]) result |= mask;
596 mask <<= 1;
597 }
598 return result;
599 } else {
600 /*
601 for (size_t qubit : qubits)
602 {
603 if (mpsSimulator->MeasureQubit(static_cast<unsigned int>(qubit)))
604 res |= mask;
605 mask <<= 1;
606 }
607 */
608 const std::set<Eigen::Index> qubitsSet(qubits.begin(), qubits.end());
609 auto measured = mpsSimulator->MeasureQubits(qubitsSet);
610 for (Types::qubit_t qubit : qubits) {
611 if (measured[qubit]) res |= mask;
612 mask <<= 1;
613 }
614 }
615 Notify();
616
617 NotifyObservers(qubits);
618
619 return res;
620 }
621
628 std::vector<bool> MeasureMany(const Types::qubits_vector &qubits) override {
629 std::vector<bool> res(qubits.size(), false);
630 DontNotify();
631
632 if (simulationType == SimulationType::kStatevector) {
633 for (size_t q = 0; q < qubits.size(); ++q)
634 if (state->MeasureQubit(static_cast<unsigned int>(qubits[q])))
635 res[q] = true;
636 } else if (simulationType == SimulationType::kStabilizer) {
637 for (size_t q = 0; q < qubits.size(); ++q)
638 if (cliffordSimulator->MeasureQubit(
639 static_cast<unsigned int>(qubits[q])))
640 res[q] = true;
641 } else if (simulationType == SimulationType::kTensorNetwork) {
642 for (size_t q = 0; q < qubits.size(); ++q)
643 if (tensorNetwork->Measure(static_cast<unsigned int>(qubits[q])))
644 res[q] = true;
645 } else if (simulationType == SimulationType::kPauliPropagator) {
646 std::vector<int> qubitsInt(qubits.begin(), qubits.end());
647 res = pp->Measure(qubitsInt);
648 } else {
649 const std::set<Eigen::Index> qubitsSet(qubits.begin(), qubits.end());
650 auto measured = mpsSimulator->MeasureQubits(qubitsSet);
651 for (size_t q = 0; q < qubits.size(); ++q)
652 if (measured[qubits[q]]) res[q] = true;
653 }
654 Notify();
655 NotifyObservers(qubits);
656
657 return res;
658 }
659
666 void ApplyReset(const Types::qubits_vector &qubits) override {
667 QC::Gates::PauliXGate xGate;
668
669 DontNotify();
670 if (simulationType == SimulationType::kStatevector) {
671 for (size_t qubit : qubits)
672 if (state->MeasureQubit(static_cast<unsigned int>(qubit)))
673 state->ApplyGate(xGate, static_cast<unsigned int>(qubit));
674 } else if (simulationType == SimulationType::kStabilizer) {
675 for (size_t qubit : qubits)
676 if (cliffordSimulator->MeasureQubit(static_cast<unsigned int>(qubit)))
677 cliffordSimulator->ApplyX(static_cast<unsigned int>(qubit));
678 } else if (simulationType == SimulationType::kTensorNetwork) {
679 for (size_t qubit : qubits)
680 if (tensorNetwork->Measure(static_cast<unsigned int>(qubit)))
681 tensorNetwork->AddGate(xGate, static_cast<unsigned int>(qubit));
682 } else if (simulationType == SimulationType::kPauliPropagator) {
683 std::vector<int> qubitsInt(qubits.begin(), qubits.end());
684 const auto res = pp->Measure(qubitsInt);
685 for (size_t i = 0; i < res.size(); ++i) {
686 if (res[i]) pp->ApplyX(qubitsInt[i]);
687 }
688 } else {
689 for (size_t qubit : qubits)
690 if (mpsSimulator->MeasureQubit(static_cast<unsigned int>(qubit)))
691 mpsSimulator->ApplyGate(xGate, static_cast<unsigned int>(qubit));
692 }
693 Notify();
694
695 NotifyObservers(qubits);
696 }
697
709 double Probability(Types::qubit_t outcome) override {
710 if (simulationType == SimulationType::kMatrixProductState)
711 return mpsSimulator->getBasisStateProbability(
712 static_cast<unsigned int>(outcome));
713 else if (simulationType == SimulationType::kStabilizer)
714 return cliffordSimulator->getBasisStateProbability(
715 static_cast<unsigned int>(outcome));
716 else if (simulationType == SimulationType::kTensorNetwork)
717 return tensorNetwork->getBasisStateProbability(outcome);
718 else if (simulationType == SimulationType::kPauliPropagator)
719 return pp->Probability(outcome);
720
721 return state->getBasisStateProbability(static_cast<unsigned int>(outcome));
722 }
723
734 std::complex<double> Amplitude(Types::qubit_t outcome) override {
735 if (simulationType == SimulationType::kMatrixProductState)
736 return mpsSimulator->getBasisStateAmplitude(
737 static_cast<unsigned int>(outcome));
738 else if (simulationType == SimulationType::kStabilizer)
739 throw std::runtime_error(
740 "QCSimState::Amplitude: Invalid simulation type for obtaining the "
741 "amplitude of the specified outcome.");
742 else if (simulationType == SimulationType::kTensorNetwork)
743 throw std::runtime_error(
744 "QCSimState::Amplitude: Not supported for the "
745 "tensor network simulator.");
746 else if (simulationType == SimulationType::kPauliPropagator)
747 throw std::runtime_error(
748 "QCSimState::Amplitude: Invalid simulation type for obtaining the "
749 "amplitude of the specified outcome.");
750
751 return state->getBasisStateAmplitude(static_cast<unsigned int>(outcome));
752 }
753
767 std::complex<double> ProjectOnZero() override {
768 if (simulationType == SimulationType::kMatrixProductState)
769 return mpsSimulator->ProjectOnZero();
770 return Amplitude(0);
771 }
772
783 std::vector<double> AllProbabilities() override {
784 // TODO: In principle this could be done, but why? It should be costly.
785 if (simulationType == SimulationType::kTensorNetwork)
786 throw std::runtime_error(
787 "QCSimState::AllProbabilities: Invalid "
788 "simulation type for obtaining probabilities.");
789 else if (simulationType == SimulationType::kStabilizer)
790 return cliffordSimulator->AllProbabilities();
791 else if (simulationType == SimulationType::kPauliPropagator) {
792 const size_t nrBasisStates = 1ULL << GetNumberOfQubits();
793 std::vector<double> result(nrBasisStates);
794 for (size_t i = 0; i < nrBasisStates; ++i) result[i] = pp->Probability(i);
795 return result;
796 }
797
798 const Eigen::VectorXcd probs =
799 simulationType == SimulationType::kMatrixProductState
800 ? mpsSimulator->getRegisterStorage().cwiseAbs2()
801 : state->getRegisterStorage().cwiseAbs2();
802
803 std::vector<double> result(probs.size());
804
805 for (int i = 0; i < probs.size(); ++i) result[i] = probs[i].real();
806
807 return result;
808 }
809
821 std::vector<double> Probabilities(
822 const Types::qubits_vector &qubits) override {
823 if (simulationType == SimulationType::kStabilizer)
824 throw std::runtime_error(
825 "QCSimState::Probabilities: Invalid simulation "
826 "type for obtaining probabilities.");
827 else if (simulationType == SimulationType::kTensorNetwork) {
828 // TODO: Implement this!!!
829 throw std::runtime_error(
830 "QCSimState::Probabilities: Not implemented yet "
831 "for the tensor network simulator.");
832 }
833
834 std::vector<double> result(qubits.size());
835
836 if (simulationType == SimulationType::kMatrixProductState) {
837 for (int i = 0; i < static_cast<int>(qubits.size()); ++i)
838 result[i] = mpsSimulator->getBasisStateProbability(qubits[i]);
839 } else if (simulationType == SimulationType::kPauliPropagator) {
840 for (int i = 0; i < static_cast<int>(qubits.size()); ++i)
841 result[i] = pp->Probability(qubits[i]);
842 } else {
843 const Eigen::VectorXcd &reg = state->getRegisterStorage();
844
845 for (int i = 0; i < static_cast<int>(qubits.size()); ++i)
846 result[i] = std::norm(reg[qubits[i]]);
847 }
848
849 return result;
850 }
851
868 std::unordered_map<Types::qubit_t, Types::qubit_t> SampleCounts(
869 const Types::qubits_vector &qubits, size_t shots = 1000) override {
870 if (qubits.empty() || shots == 0) return {};
871
872 if (qubits.size() > sizeof(size_t) * 8)
873 std::cerr
874 << "Warning: The number of qubits to measure is larger than the "
875 "number of bits in the size_t type, the outcome will be undefined"
876 << std::endl;
877
878 // TODO: this is inefficient, maybe implement it better in qcsim
879 // for now it has the possibility of measuring a qubits interval, but not a
880 // list of qubits
881 std::unordered_map<Types::qubit_t, Types::qubit_t> result;
882
883 DontNotify();
884
885 if (simulationType == SimulationType::kMatrixProductState) {
886 bool normal = true;
887 if (useMPSMeasureNoCollapse) {
888 // check to see if it can be used
889 const std::set<Eigen::Index> qset(qubits.begin(), qubits.end());
890 if (qset.size() == GetNumberOfQubits()) {
891 // it can!
892 normal = false;
893 for (size_t shot = 0; shot < shots; ++shot) {
894 const size_t measRaw = MeasureNoCollapse();
895 size_t meas = 0;
896 size_t mask = 1ULL;
897
898 // translate the measurement
899 for (auto q : qubits) {
900 const size_t qubitMask = 1ULL << q;
901 if (measRaw & qubitMask) meas |= mask;
902 mask <<= 1ULL;
903 }
904
905 ++result[meas];
906 }
907 } else if (qset.size() > 1) {
908 mpsSimulator->MoveAtBeginningOfChain(qset);
909 // now sample
910 normal = false;
911 for (size_t shot = 0; shot < shots; ++shot) {
912 const auto measRaw = mpsSimulator->MeasureNoCollapse(qset);
913 size_t meas = 0;
914 size_t mask = 1ULL;
915
916 // might not be in the requested order
917 // translate the measurement
918 for (auto q : qubits) {
919 const size_t qubitMask = 1ULL << q;
920 if (measRaw.at(q)) meas |= mask;
921 mask <<= 1ULL;
922 }
923
924 ++result[meas];
925 }
926
927 } else if (qset.size() == 1) {
928 // if only one qubit is measured, we can use the probability
929 normal = false;
930 const auto prob0 = mpsSimulator->GetProbability(qubits[0]);
931 for (size_t shot = 0; shot < shots; ++shot) {
932 const size_t meas = uniformZeroOne(rng) < prob0 ? 0ULL : 1ULL;
933 size_t m = meas;
934 // why would somebody set more than one time?
935 for (size_t i = 1; i < qubits.size(); ++i) {
936 m <<= 1ULL;
937 m |= meas;
938 }
939 ++result[m];
940 }
941 }
942 }
943
944 if (normal) {
945 auto savedState = mpsSimulator->getState();
946 for (size_t shot = 0; shot < shots; ++shot) {
947 const size_t meas = Measure(qubits);
948 ++result[meas];
949 mpsSimulator->setState(savedState);
950 }
951 }
952 } else if (simulationType == SimulationType::kStabilizer) {
953 cliffordSimulator->SaveState();
954 for (size_t shot = 0; shot < shots; ++shot) {
955 const size_t meas = Measure(qubits);
956 ++result[meas];
957 cliffordSimulator->RestoreState();
958 }
959 cliffordSimulator->ClearSavedState();
960 } else if (simulationType == SimulationType::kTensorNetwork) {
961 tensorNetwork->SaveState();
962 for (size_t shot = 0; shot < shots; ++shot) {
963 const size_t meas = Measure(qubits);
964 ++result[meas];
965 tensorNetwork->RestoreState();
966 }
967 tensorNetwork->ClearSavedState();
968 } else if (simulationType == SimulationType::kPauliPropagator) {
969 std::vector<int> qubitsInt(qubits.begin(), qubits.end());
970 for (size_t shot = 0; shot < shots; ++shot) {
971 const auto res = pp->Sample(qubitsInt);
972
973 size_t meas = 0;
974 for (size_t i = 0; i < qubits.size(); ++i) {
975 if (res[i]) meas |= (1ULL << i);
976 }
977
978 ++result[meas];
979 }
980 } else {
981 if (shots > 1) {
982 const auto &statev = state->getRegisterStorage();
983
984 const Utils::Alias alias(statev);
985
986 for (size_t shot = 0; shot < shots; ++shot) {
987 const double prob = 1. - uniformZeroOne(rng);
988 const size_t measRaw = alias.Sample(prob);
989
990 size_t meas = 0;
991 size_t mask = 1ULL;
992 for (auto q : qubits) {
993 const size_t qubitMask = 1ULL << q;
994 if ((measRaw & qubitMask) != 0) meas |= mask;
995 mask <<= 1ULL;
996 }
997
998 ++result[meas];
999 }
1000 } else {
1001 for (size_t shot = 0; shot < shots; ++shot) {
1002 const size_t measRaw = MeasureNoCollapse();
1003 size_t meas = 0;
1004 size_t mask = 1ULL;
1005
1006 for (auto q : qubits) {
1007 const size_t qubitMask = 1ULL << q;
1008 if ((measRaw & qubitMask) != 0) meas |= mask;
1009 mask <<= 1ULL;
1010 }
1011
1012 ++result[meas];
1013 }
1014 }
1015 }
1016
1017 Notify();
1018 NotifyObservers(qubits);
1019
1020 return result;
1021 }
1022
1036 std::unordered_map<std::vector<bool>, Types::qubit_t> SampleCountsMany(
1037 const Types::qubits_vector &qubits, size_t shots = 1000) override {
1038 if (qubits.empty() || shots == 0) return {};
1039
1040 std::unordered_map<std::vector<bool>, Types::qubit_t> result;
1041
1042 DontNotify();
1043
1044 if (simulationType == SimulationType::kMatrixProductState) {
1045 bool normal = true;
1046 if (useMPSMeasureNoCollapse) {
1047 // check to see if it can be used
1048 const std::set<Eigen::Index> qset(qubits.begin(), qubits.end());
1049 if (qset.size() == GetNumberOfQubits()) {
1050 // it can!
1051 normal = false;
1052 for (size_t shot = 0; shot < shots; ++shot) {
1053 const auto meas = MeasureNoCollapseMany();
1054
1055 // might not be in the requested order
1056 // translate the measurement
1057 std::vector<bool> measVec(qubits.size());
1058 for (size_t i = 0; i < qubits.size(); ++i)
1059 measVec[i] = meas[qubits[i]];
1060
1061 ++result[measVec];
1062 }
1063 } else if (qset.size() > 1) {
1064 mpsSimulator->MoveAtBeginningOfChain(qset);
1065 // now sample
1066 normal = false;
1067 for (size_t shot = 0; shot < shots; ++shot) {
1068 const auto meas = mpsSimulator->MeasureNoCollapse(qset);
1069
1070 // might not be in the requested order
1071 // translate the measurement
1072 std::vector<bool> measVec(qubits.size());
1073 for (size_t i = 0; i < qubits.size(); ++i)
1074 measVec[i] = meas.at(qubits[i]);
1075
1076 ++result[measVec];
1077 }
1078 } else if (qset.size() == 1) {
1079 // if only one qubit is measured, we can use the probability
1080 normal = false;
1081 const auto prob0 = mpsSimulator->GetProbability(qubits[0]);
1082 for (size_t shot = 0; shot < shots; ++shot) {
1083 const size_t meas = uniformZeroOne(rng) < prob0 ? 0ULL : 1ULL;
1084 const std::vector<bool> m(qubits.size(), meas);
1085 ++result[m];
1086 }
1087 }
1088 }
1089
1090 if (normal) {
1091 auto savedState = mpsSimulator->getState();
1092 for (size_t shot = 0; shot < shots; ++shot) {
1093 const auto meas = MeasureMany(qubits);
1094
1095 ++result[meas];
1096 mpsSimulator->setState(savedState);
1097 }
1098 }
1099 } else if (simulationType == SimulationType::kStabilizer) {
1100 cliffordSimulator->SaveState();
1101 for (size_t shot = 0; shot < shots; ++shot) {
1102 const auto meas = MeasureMany(qubits);
1103 ++result[meas];
1104 cliffordSimulator->RestoreState();
1105 }
1106 cliffordSimulator->ClearSavedState();
1107 } else if (simulationType == SimulationType::kTensorNetwork) {
1108 tensorNetwork->SaveState();
1109 for (size_t shot = 0; shot < shots; ++shot) {
1110 const auto meas = MeasureMany(qubits);
1111 ++result[meas];
1112 tensorNetwork->RestoreState();
1113 }
1114 tensorNetwork->ClearSavedState();
1115 } else if (simulationType == SimulationType::kPauliPropagator) {
1116 std::vector<int> qubitsInt(qubits.begin(), qubits.end());
1117 for (size_t shot = 0; shot < shots; ++shot) {
1118 const auto meas = pp->Sample(qubitsInt);
1119 ++result[meas];
1120 }
1121 } else {
1122 if (shots > 1) {
1123 const auto &statev = state->getRegisterStorage();
1124
1125 const Utils::Alias alias(statev);
1126
1127 for (size_t shot = 0; shot < shots; ++shot) {
1128 const double prob = 1. - uniformZeroOne(rng);
1129 const size_t measRaw = alias.Sample(prob);
1130
1131 std::vector<bool> meas(qubits.size(), false);
1132
1133 for (size_t i = 0; i < qubits.size(); ++i)
1134 if (((measRaw >> qubits[i]) & 1) == 1) meas[i] = true;
1135
1136 ++result[meas];
1137 }
1138 } else {
1139 for (size_t shot = 0; shot < shots; ++shot) {
1140 const auto measRaw = MeasureNoCollapseMany();
1141 std::vector<bool> meas(qubits.size(), false);
1142
1143 for (size_t i = 0; i < qubits.size(); ++i)
1144 if (measRaw[qubits[i]]) meas[i] = true;
1145
1146 ++result[meas];
1147 }
1148 }
1149 }
1150
1151 Notify();
1152 NotifyObservers(qubits);
1153
1154 return result;
1155 }
1156
1168 double ExpectationValue(const std::string &pauliStringOrig) override {
1169 if (pauliStringOrig.empty()) return 1.0;
1170
1171 std::string pauliString = pauliStringOrig;
1172 if (pauliString.size() > GetNumberOfQubits()) {
1173 for (size_t i = GetNumberOfQubits(); i < pauliString.size(); ++i) {
1174 const auto pauliOp = toupper(pauliString[i]);
1175 if (pauliOp != 'I' && pauliOp != 'Z') return 0.0;
1176 }
1177
1178 pauliString.resize(GetNumberOfQubits());
1179 }
1180
1181 if (simulationType == SimulationType::kStabilizer)
1182 return cliffordSimulator->ExpectationValue(pauliString);
1183 else if (simulationType == SimulationType::kTensorNetwork)
1184 return tensorNetwork->ExpectationValue(pauliString);
1185 else if (simulationType == SimulationType::kPauliPropagator)
1186 return pp->ExpectationValue(pauliString);
1187
1188 // statevector or mps
1189 static const QC::Gates::PauliXGate<> xgate;
1190 static const QC::Gates::PauliYGate<> ygate;
1191 static const QC::Gates::PauliZGate<> zgate;
1192
1193 std::vector<QC::Gates::AppliedGate<Eigen::MatrixXcd>> pauliStringVec;
1194 pauliStringVec.reserve(pauliString.size());
1195
1196 for (size_t q = 0; q < pauliString.size(); ++q) {
1197 switch (toupper(pauliString[q])) {
1198 case 'X': {
1199 QC::Gates::AppliedGate<Eigen::MatrixXcd> ag(
1200 xgate.getRawOperatorMatrix(), static_cast<Types::qubit_t>(q));
1201 pauliStringVec.emplace_back(std::move(ag));
1202 } break;
1203 case 'Y': {
1204 QC::Gates::AppliedGate<Eigen::MatrixXcd> ag(
1205 ygate.getRawOperatorMatrix(), static_cast<Types::qubit_t>(q));
1206 pauliStringVec.emplace_back(std::move(ag));
1207 } break;
1208 case 'Z': {
1209 QC::Gates::AppliedGate<Eigen::MatrixXcd> ag(
1210 zgate.getRawOperatorMatrix(), static_cast<Types::qubit_t>(q));
1211 pauliStringVec.emplace_back(std::move(ag));
1212 } break;
1213 case 'I':
1214 [[fallthrough]];
1215 default:
1216 break;
1217 }
1218 }
1219
1220 if (pauliStringVec.empty()) return 1.0;
1221
1222 if (simulationType == SimulationType::kMatrixProductState)
1223 return mpsSimulator->ExpectationValue(pauliStringVec).real();
1224
1225 return state->ExpectationValue(pauliStringVec).real();
1226 }
1227
1235 SimulatorType GetType() const override { return SimulatorType::kQCSim; }
1236
1245 SimulationType GetSimulationType() const override { return simulationType; }
1246
1255 void Flush() override {}
1256
1266 void SaveStateToInternalDestructive() override {}
1267
1274 void RestoreInternalDestructiveSavedState() override {}
1275
1285 void SaveState() override {
1286 if (simulationType == SimulationType::kMatrixProductState)
1287 mpsSimulator->SaveState();
1288 else if (simulationType == SimulationType::kStabilizer)
1289 cliffordSimulator->SaveState();
1290 else if (simulationType == SimulationType::kTensorNetwork)
1291 tensorNetwork->SaveState();
1292 else if (simulationType == SimulationType::kPauliPropagator)
1293 pp->SaveState();
1294 else
1295 state->SaveState();
1296 }
1297
1306 void RestoreState() override {
1307 if (simulationType == SimulationType::kMatrixProductState)
1308 mpsSimulator->RestoreState();
1309 else if (simulationType == SimulationType::kStabilizer)
1310 cliffordSimulator->RestoreState();
1311 else if (simulationType == SimulationType::kTensorNetwork)
1312 tensorNetwork->RestoreState();
1313 else if (simulationType == SimulationType::kPauliPropagator)
1314 pp->RestoreState();
1315 else
1316 state->RestoreState();
1317 }
1318
1326 std::complex<double> AmplitudeRaw(Types::qubit_t outcome) override {
1327 return Amplitude(outcome);
1328 }
1329
1338 void SetMultithreading(bool multithreading = true) override {
1339 enableMultithreading = multithreading;
1340 if (state) state->SetMultithreading(multithreading);
1341 if (cliffordSimulator) cliffordSimulator->SetMultithreading(multithreading);
1342 if (tensorNetwork) tensorNetwork->SetMultithreading(multithreading);
1343 if (pp) {
1344 if (multithreading)
1345 pp->EnableParallel();
1346 else
1347 pp->DisableParallel();
1348 }
1349 }
1350
1358 bool GetMultithreading() const override { return enableMultithreading; }
1359
1370 bool IsQcsim() const override { return true; }
1371
1389 if (GetNumberOfQubits() > sizeof(Types::qubit_t) * 8)
1390 std::cerr
1391 << "Warning: The number of qubits to measure is larger than the "
1392 "number of bits in the Types::qubit_t type, the outcome will be "
1393 "undefined"
1394 << std::endl;
1395
1396 if (simulationType == SimulationType::kStatevector)
1397 return state->MeasureNoCollapse();
1398 else if (simulationType == SimulationType::kMatrixProductState) {
1399 const auto measured = mpsSimulator->MeasureNoCollapse();
1400 Types::qubit_t result = 0;
1401 Types::qubit_t mask = 1;
1402 for (Types::qubit_t q = 0; q < measured.size(); ++q) {
1403 if (measured.at(q)) result |= mask;
1404 mask <<= 1;
1405 }
1406 return result;
1407 } else if (simulationType == SimulationType::kPauliPropagator) {
1408 std::vector<int> qubitsInt(GetNumberOfQubits());
1409 std::iota(qubitsInt.begin(), qubitsInt.end(), 0);
1410 const auto res = pp->Sample(qubitsInt);
1411 Types::qubit_t result = 0;
1412 for (size_t i = 0; i < res.size(); ++i) {
1413 if (res[i]) result |= (1ULL << i);
1414 }
1415 return result;
1416 }
1417
1418 throw std::runtime_error(
1419 "QCSimState::MeasureNoCollapse: Invalid simulation type for measuring "
1420 "all the qubits without collapsing the state.");
1421
1422 return 0;
1423 }
1424
1439 std::vector<bool> MeasureNoCollapseMany() override {
1440 if (simulationType == SimulationType::kStatevector) {
1441 auto state = MeasureNoCollapse();
1442 std::vector<bool> res(nrQubits);
1443 for (size_t i = 0; i < nrQubits; ++i) res[i] = ((state >> i) & 1) == 1;
1444 return res;
1445 } else if (simulationType == SimulationType::kMatrixProductState) {
1446 const auto measured = mpsSimulator->MeasureNoCollapse();
1447 std::vector<bool> res(nrQubits);
1448 for (size_t i = 0; i < nrQubits; ++i) res[i] = measured.at(i);
1449 return res;
1450 } else if (simulationType == SimulationType::kPauliPropagator) {
1451 std::vector<int> qubitsInt(GetNumberOfQubits());
1452 std::iota(qubitsInt.begin(), qubitsInt.end(), 0);
1453 return pp->Sample(qubitsInt);
1454 }
1455
1456 throw std::runtime_error(
1457 "QCSimState::MeasureNoCollapseMany: Invalid simulation type for "
1458 "measuring all the qubits without collapsing the state.");
1459
1460 return {};
1461 }
1462
1463 protected:
1464 SimulationType simulationType =
1465 SimulationType::kStatevector;
1467 std::unique_ptr<QC::QubitRegister<>> state;
1468 std::unique_ptr<QC::TensorNetworks::MPSSimulator>
1469 mpsSimulator;
1470 std::unique_ptr<QC::Clifford::StabilizerSimulator>
1471 cliffordSimulator;
1472 std::unique_ptr<TensorNetworks::TensorNetwork>
1473 tensorNetwork;
1474 std::unique_ptr<QcsimPauliPropagator> pp;
1476 size_t nrQubits = 0;
1477 bool limitSize = false;
1478 bool limitEntanglement = false;
1479 Eigen::Index chi = 10; // if limitSize is true
1480 double singularValueThreshold = 0.; // if limitEntanglement is true
1481 bool enableMultithreading = true;
1482 bool useMPSMeasureNoCollapse =
1483 true;
1485 int lookaheadDepth = 0;
1486 int lookaheadDepthWithHeuristic = 0;
1487 bool useOptimalMeetingPositionOnly = false;
1488 std::vector<std::shared_ptr<Circuits::IOperation<>>> upcomingGates;
1489 long long int upcomingGateIndex = 0;
1490 std::unique_ptr<Simulators::MPSDummySimulator> dummySim;
1491
1492 // Observer that counts applied gates to track position in upcomingGates
1493 class GateCounterObserver : public ISimulatorObserver {
1494 public:
1495 GateCounterObserver(long long int &indexRef) : index(indexRef) {}
1496 void Update(const Types::qubits_vector &) override {
1497 ++index;
1498 }
1499 private:
1500 long long int &index;
1501 };
1502 std::shared_ptr<GateCounterObserver> gateCounterObserver;
1503
1504 std::mt19937_64 rng;
1505 std::uniform_real_distribution<double> uniformZeroOne;
1506};
1507
1508} // namespace Private
1509} // namespace Simulators
1510
1511#endif
1512
1513#endif // !_QCSIMSTATE_H_
double Probability(void *sim, unsigned long long int outcome)
char * GetConfiguration(void *sim, const char *key)
int RestoreState(void *sim)
int ApplyReset(void *sim, const unsigned long int *qubits, unsigned long int nrQubits)
unsigned long int AllocateQubits(void *sim, unsigned long int nrQubits)
unsigned long int GetNumberOfQubits(void *sim)
double * AllProbabilities(void *sim)
unsigned long long int MeasureNoCollapse(void *sim)
int GetMultithreading(void *sim)
unsigned long long int Measure(void *sim, const unsigned long int *qubits, unsigned long int nrQubits)
double * Amplitude(void *sim, unsigned long long int outcome)
double * Probabilities(void *sim, const unsigned long long int *qubits, unsigned long int nrQubits)
int SetMultithreading(void *sim, int multithreading)
int SaveStateToInternalDestructive(void *sim)
int GetSimulationType(void *sim)
unsigned long long int * SampleCounts(void *sim, const unsigned long long int *qubits, unsigned long int nrQubits, unsigned long int shots)
int RestoreInternalDestructiveSavedState(void *sim)
int IsQcsim(void *sim)
int SaveState(void *sim)
The operation interface.
Definition Operations.h:358
size_t Sample(double v) const
Definition Alias.h:85
uint_fast64_t qubit_t
The type of a qubit.
Definition Types.h:21
std::vector< qubit_t > qubits_vector
The type of a vector of qubits.
Definition Types.h:23