Maestro 0.2.11
Unified interface for quantum circuit simulation
Loading...
Searching...
No Matches
MPSDummySimulator.h
Go to the documentation of this file.
1
15
16#pragma once
17
18#ifndef _MPSDUMMYSIMULATOR_H_
19#define _MPSDUMMYSIMULATOR_H_
20
21#include <algorithm>
22#include <deque>
23#include <random>
24
25#include <Eigen/Eigen>
26#include <unsupported/Eigen/CXX11/Tensor>
27
28#include "QuantumGate.h"
29#include "MPSSimulatorInterface.h"
30
31namespace Simulators {
32
34 public:
35 using IndexType = long long int;
36 using MatrixClass = QC::TensorNetworks::MPSSimulatorInterface::MatrixClass;
37 using GateClass = QC::TensorNetworks::MPSSimulatorInterface::GateClass;
38
39 MPSDummySimulator(size_t N) : nrQubits(N), maxVirtualExtent(0) {
40 InitQubitsMap();
41
43 }
44
45 std::unique_ptr<MPSDummySimulator> Clone() const {
46 auto clone = std::unique_ptr<MPSDummySimulator>(
47 new MPSDummySimulator(nrQubits, LightweightInitTag{}));
48 clone->qubitsMap = qubitsMap;
49 clone->qubitsMapInv = qubitsMapInv;
50 clone->maxVirtualExtent = maxVirtualExtent;
51 clone->bondCost = bondCost;
52 clone->maxBondDim = maxBondDim;
53 clone->currentBondDim = currentBondDim;
54
55 clone->totalSwappingCost = totalSwappingCost;
56
57 clone->growthFactorSwap = growthFactorSwap;
58 clone->growthFactorGate = growthFactorGate;
59
60 return clone;
61 }
62
63 size_t getNrQubits() const { return nrQubits; }
64
65 double getGrowthFactorSwap() const { return growthFactorSwap; }
66 double getGrowthFactorGate() const { return growthFactorGate; }
67
68 void setGrowthFactorSwap(double factor) { growthFactorSwap = factor; }
69 void setGrowthFactorGate(double factor) { growthFactorGate = factor; }
70
71 void Clear() { InitQubitsMap(); }
72
74 maxVirtualExtent = val;
75
76 if (nrQubits == 0) return;
77
78 const double untruncatedMaxExtent = std::pow(physExtent, nrQubits / 2);
79 IndexType maxVirtualExtentLimit =
80 static_cast<IndexType>(untruncatedMaxExtent);
81
82 if (untruncatedMaxExtent >= std::numeric_limits<IndexType>::max() ||
83 std::isnan(untruncatedMaxExtent) || std::isinf(untruncatedMaxExtent))
84 maxVirtualExtentLimit = std::numeric_limits<IndexType>::max() - 1;
85 else if (untruncatedMaxExtent < 2)
86 maxVirtualExtentLimit = 2;
87
88 maxVirtualExtent = maxVirtualExtent == 0
89 ? maxVirtualExtentLimit
90 : std::min(maxVirtualExtent, maxVirtualExtentLimit);
91
92 bondCost.resize(nrQubits - 1);
93 maxBondDim.resize(nrQubits - 1);
94 currentBondDim.assign(nrQubits - 1, 1.0);
95
96 // the checks here are overkill, but better safe than sorry
97 // we're dealing with large values here, overflows are to be expected
98 for (IndexType i = 0; i < static_cast<IndexType>(nrQubits); ++i) {
99 double maxExtent1 = std::pow((double)physExtent, (double)i + 1.);
100 double maxExtent2 =
101 std::pow((double)physExtent, (double)nrQubits - i - 1.);
102
103 if (maxExtent1 >= (double)std::numeric_limits<size_t>::max() ||
104 std::isnan(maxExtent1) || std::isinf(maxExtent1))
105 maxExtent1 = (double)std::numeric_limits<size_t>::max() - 1;
106 else if (maxExtent1 < 1)
107 maxExtent1 = 1;
108
109 if (maxExtent2 > (double)std::numeric_limits<size_t>::max() ||
110 std::isnan(maxExtent2) || std::isinf(maxExtent2))
111 maxExtent2 = (double)std::numeric_limits<size_t>::max() - 1;
112 else if (maxExtent2 < 1)
113 maxExtent2 = 1;
114
115 size_t maxRightExtent = (size_t)std::min<double>(
116 {maxExtent1, maxExtent2, (double)maxVirtualExtent});
117 if (maxRightExtent < 2) maxRightExtent = 2;
118
119 if (i < static_cast<IndexType>(nrQubits) - 1) {
120 maxBondDim[i] = static_cast<double>(maxRightExtent);
121 const double maxRightExtentD = static_cast<double>(maxRightExtent);
122 bondCost[i] = maxRightExtentD * maxRightExtentD * maxRightExtentD;
123 }
124 }
125 }
126
127 void print() const {
128 std::cout << "Qubits map: ";
129 for (int q = 0; q < static_cast<int>(qubitsMap.size()); ++q)
130 std::cout << q << "->" << qubitsMap[q] << " ";
131 std::cout << std::endl;
132 }
133
134 void ApplyGate(const QC::Gates::AppliedGate<MatrixClass>& gate) {
135 ApplyGate(gate, gate.getQubit1(), gate.getQubit2());
136 }
137
138 void ApplyGate(const std::shared_ptr<Circuits::IOperation<>>& gate) {
139 const auto qbits = gate->AffectedQubits();
140
141 if (qbits.size() == 1) {
142 // 1-qubit gates are no-ops in the dummy simulator (no swap logic)
143 return;
144 } else if (qbits.size() == 2) {
145 // Reuse a single static dummy 2-qubit gate to avoid heap allocation.
146 // Only getQubitsNumber() is checked; the matrix is never read.
147 static const QC::Gates::AppliedGate<MatrixClass> dummy2qGate(
148 MatrixClass::Identity(4, 4), 0, 1);
149 ApplyGate(dummy2qGate, qbits[0], qbits[1]);
150 } else if (qbits.size() == 3 && gate->GetType() == Circuits::OperationType::kGate) {
151 static const QC::Gates::AppliedGate<MatrixClass> dummy2qGate(
152 MatrixClass::Identity(4, 4), 0, 1);
153
154 const auto gateptr =
155 std::static_pointer_cast<Circuits::IQuantumGate<>>(gate);
156
157 const size_t q1 = qbits[0]; // control 1
158 const size_t q2 = qbits[1]; // control 2
159 const size_t q3 = qbits[2]; // target
160 if (gateptr->GetGateType() == Circuits::QuantumGateType::kCCXGateType) {
161 ApplyGate(dummy2qGate, q3, q2);
162 ApplyGate(dummy2qGate, q2, q1);
163 ApplyGate(dummy2qGate, q3, q2);
164 ApplyGate(dummy2qGate, q2, q1);
165 ApplyGate(dummy2qGate, q3, q1);
166 } else { // cswap
167 ApplyGate(dummy2qGate, q2, q3);
168 ApplyGate(dummy2qGate, q3, q2);
169 ApplyGate(dummy2qGate, q2, q1);
170 ApplyGate(dummy2qGate, q3, q2);
171 ApplyGate(dummy2qGate, q2, q1);
172 ApplyGate(dummy2qGate, q3, q1);
173 ApplyGate(dummy2qGate, q2, q3);
174 }
175 }
176 }
177
178 void ApplyGate(const GateClass& gate, IndexType qubit,
179 IndexType controllingQubit1 = 0) {
180 if (qubit < 0 || qubit >= static_cast<IndexType>(getNrQubits()))
181 throw std::invalid_argument("Qubit index out of bounds");
182 else if (controllingQubit1 < 0 ||
183 controllingQubit1 >= static_cast<IndexType>(getNrQubits()))
184 throw std::invalid_argument("Qubit index out of bounds");
185
186 // for two qubit gates:
187 // if the qubits are not adjacent, apply swap gates until they are
188 // don't forget to update the qubitsMap
189 if (gate.getQubitsNumber() > 1) {
190 IndexType qubit1 = qubitsMap[qubit];
191 IndexType qubit2 = qubitsMap[controllingQubit1];
192
193 if (abs(qubit1 - qubit2) > 1) {
194 SwapQubits(qubit, controllingQubit1);
195 qubit1 = qubitsMap[qubit];
196 qubit2 = qubitsMap[controllingQubit1];
197 assert(abs(qubit1 - qubit2) == 1);
198 }
199
200 // any 2-qubit gate (whether swaps were needed or not) grows the
201 // bond dimension at the bond between the two adjacent qubits.
202 // Charge the cost based on the bond dimension *before* the growth,
203 // consistent with SwapQubitsToPosition: the cost reflects the size of
204 // the SVD/contraction actually performed for this operation, and the
205 // grown bond dimension only affects subsequent operations.
206 const IndexType bond = std::min(qubit1, qubit2);
207 totalSwappingCost += bondCost[bond];
208 growBondDimension(bond, false);
209 }
210 }
211
213 const std::vector<QC::Gates::AppliedGate<MatrixClass>>& gates) {
214 for (const auto& gate : gates) ApplyGate(gate);
215 }
216
218 const std::vector<std::shared_ptr<Circuits::IOperation<>>>& gates) {
219 for (const auto& gate : gates) ApplyGate(gate);
220 }
221
222 void SetInitialQubitsMap(const std::vector<long long int>& initialMap) {
223 qubitsMap = initialMap;
224 for (size_t i = 0; i < initialMap.size(); ++i)
225 qubitsMapInv[initialMap[i]] = i;
226
227 totalSwappingCost = 0;
228 std::fill(currentBondDim.begin(), currentBondDim.end(), 1.0);
229 for (size_t i = 0; i < bondCost.size(); ++i) bondCost[i] = 1;
230 }
231
232 void SetCurrentBondDimensions(const std::vector<double>& dims) {
233 assert(dims.size() == currentBondDim.size());
234 for (size_t i = 0; i < dims.size(); ++i)
235 currentBondDim[i] = std::min(dims[i], maxBondDim[i]);
236
237 // update bondCost as well, since it depends on the current bond dimensions
238 for (size_t i = 0; i < currentBondDim.size(); ++i)
239 bondCost[i] = currentBondDim[i] * currentBondDim[i] * currentBondDim[i];
240 }
241
242 const std::vector<double>& getCurrentBondDimensions() const {
243 return currentBondDim;
244 }
245
246 const std::vector<double>& getMaxBondDimensions() const { return maxBondDim; }
247
248 const std::vector<IndexType>& getQubitsMap() const { return qubitsMap; }
249
250 const std::vector<IndexType>& getQubitsMapInv() const { return qubitsMapInv; }
251
252 void setTotalSwappingCost(double cost) { totalSwappingCost = cost; }
253 double getTotalSwappingCost() const { return totalSwappingCost; }
254
255 // Evaluate the total cost
256 // position, applying the current 2-qubit gate, and then simulating the
257 // next lookaheadDepth 2-qubit gates from upcomingGates.
259 IndexType meetPosition,
260 const std::vector<std::shared_ptr<Circuits::IOperation<>>>& upcomingGates,
261 long long int currentGateIndex, int lookaheadDepth,
262 int lookaheadDepthWithHeuristic, double currentCost, double& bestCost,
263 bool useSameDummy = false) {
264 if (currentGateIndex >= static_cast<long long int>(upcomingGates.size())) {
265 if (currentCost < bestCost) bestCost = currentCost;
266 return;
267 }
268
269 // skip the 1 qubit gates, advance to the next 2-qubit gate
270 while (currentGateIndex <
271 static_cast<long long int>(upcomingGates.size()) &&
272 upcomingGates[currentGateIndex]->AffectedQubits().size() < 2)
273 ++currentGateIndex;
274
275 if (currentGateIndex >= static_cast<long long int>(upcomingGates.size())) {
276 if (currentCost < bestCost) bestCost = currentCost;
277 return;
278 }
279
280 const auto& op = upcomingGates[currentGateIndex];
281 const auto qbits = op->AffectedQubits();
282
283 assert(qbits.size() >= 2);
284
285 const IndexType qubit1 = static_cast<IndexType>(qbits[0]);
286 const IndexType qubit2 = static_cast<IndexType>(qbits[1]);
287
288 // Perform the swap to the meeting position;
289 // totalSwappingCost starts at 0 and SwapQubitsToPosition accumulates
290 // per-bond costs, so the cost correctly depends on meetPosition.
291
292 IndexType realq1 = qubitsMap[qubit1];
293 IndexType realq2 = qubitsMap[qubit2];
294 if (realq1 > realq2) std::swap(realq1, realq2);
295
296 if (useSameDummy) {
297 const double ccost = getTotalSwappingCost();
298
299 if (realq2 - realq1 > 1)
300 SwapQubitsToPosition(qubit1, qubit2, meetPosition);
301
302 ApplyGate(op);
303
304 // const double importanceFactor =
305 // pow(dimFactor, (lookaheadDepthWithHeuristic - lookaheadDepth - 1));
306 currentCost += (getTotalSwappingCost() - ccost) ;
307
308 if (currentCost >= bestCost) return;
309 // No more lookahead depth: return without further recursion
310 if (lookaheadDepth <= 0) {
311 if (currentCost < bestCost) bestCost = currentCost;
312 return;
313 }
314
315 // skip the 1 qubit gates, advance over the current gate on the next
316 // 2-qubit gate
317 ++currentGateIndex;
318 while (currentGateIndex <
319 static_cast<long long int>(upcomingGates.size()) &&
320 upcomingGates[currentGateIndex]->AffectedQubits().size() < 2)
321 ++currentGateIndex;
322
323 if (currentGateIndex >=
324 static_cast<long long int>(upcomingGates.size())) {
325 if (currentCost < bestCost) bestCost = currentCost;
326 return;
327 }
328
329 FindBestMeetingPosition(upcomingGates, currentGateIndex,
330 lookaheadDepth - 1, lookaheadDepthWithHeuristic,
331 currentCost, bestCost);
332 } else {
333 MPSDummySimulator dummySim(nrQubits, LightweightInitTag{});
334 dummySim.maxVirtualExtent = maxVirtualExtent;
335 dummySim.maxBondDim = maxBondDim;
336 dummySim.currentBondDim = currentBondDim;
337 dummySim.bondCost = bondCost;
338 dummySim.qubitsMap = qubitsMap;
339 dummySim.qubitsMapInv.resize(qubitsMap.size());
340 for (size_t i = 0; i < qubitsMap.size(); ++i)
341 dummySim.qubitsMapInv[qubitsMap[i]] = static_cast<IndexType>(i);
342
343 dummySim.setTotalSwappingCost(0);
344
345 if (realq2 - realq1 > 1)
346 dummySim.SwapQubitsToPosition(qubit1, qubit2, meetPosition);
347
348 dummySim.ApplyGate(op); // this also updates the bond dimensions in
349 // dummySim according to the applied gate
350
351 currentCost += dummySim.getTotalSwappingCost();
352
353 // Pruning: if current accumulated cost already exceeds best known, prune
354 if (currentCost >= bestCost) return;
355
356 // No more lookahead depth: return without further recursion
357 if (lookaheadDepth <= 0) {
358 if (currentCost < bestCost) bestCost = currentCost;
359 return;
360 }
361
362 // skip the 1 qubit gates, advance over the current gate on the next
363 // 2-qubit gate
364 ++currentGateIndex;
365 while (currentGateIndex <
366 static_cast<long long int>(upcomingGates.size()) &&
367 upcomingGates[currentGateIndex]->AffectedQubits().size() < 2)
368 ++currentGateIndex;
369
370 if (currentGateIndex >=
371 static_cast<long long int>(upcomingGates.size())) {
372 if (currentCost < bestCost) bestCost = currentCost;
373 return;
374 }
375
376 dummySim.FindBestMeetingPosition(
377 upcomingGates, currentGateIndex, lookaheadDepth - 1,
378 lookaheadDepthWithHeuristic, currentCost, bestCost);
379 }
380 }
381
382 // Find the meeting position that minimizes the combined cost of the
383 // current swap + lookahead gates. Returns the optimal bond index.
384 // outBestCost receives the total accumulated cost for the best position.
386 const std::vector<std::shared_ptr<Circuits::IOperation<>>>& upcomingGates,
387 long long int currentGateIndex, int lookaheadDepth,
388 int lookaheadDepthWithHeuristic, double currentCost, double& bestCost) {
389 const auto& op = upcomingGates[currentGateIndex];
390 const auto qbits = op->AffectedQubits();
391
392 assert(qbits.size() >= 2);
393
394 const IndexType qubit1 = static_cast<IndexType>(qbits[0]);
395 const IndexType qubit2 = static_cast<IndexType>(qbits[1]);
396
397 IndexType realq1 = qubitsMap[qubit1];
398 IndexType realq2 = qubitsMap[qubit2];
399
400 if (realq1 > realq2) std::swap(realq1, realq2);
401
402 if (lookaheadDepth <= 0) {
403 IndexType pos = (realq2 - realq1 > 1)
404 ? ComputeHeuristicMeetPosition(realq1, realq2)
405 : realq1;
406 EvaluateMeetingPositionCost(pos, upcomingGates, currentGateIndex, 0,
407 lookaheadDepthWithHeuristic, currentCost,
408 bestCost, false);
409 return pos;
410 }
411
412 if (realq2 - realq1 <= 1) {
414 realq1, upcomingGates, currentGateIndex, lookaheadDepth,
415 lookaheadDepthWithHeuristic, currentCost, bestCost,
416 lookaheadDepth <= lookaheadDepthWithHeuristic);
417
418 return realq1;
419 }
420
421 IndexType bestPosition = realq1;
422
423 if (lookaheadDepth <= lookaheadDepthWithHeuristic) {
424 bestPosition = ComputeHeuristicMeetPosition(realq1, realq2);
425 EvaluateMeetingPositionCost(bestPosition, upcomingGates, currentGateIndex,
426 lookaheadDepth, lookaheadDepthWithHeuristic,
427 currentCost, bestCost, true);
428 } else {
429 for (IndexType m = realq1; m < realq2; ++m) {
430 const double oldCost = bestCost;
431 EvaluateMeetingPositionCost(m, upcomingGates, currentGateIndex,
432 lookaheadDepth, lookaheadDepthWithHeuristic,
433 currentCost, bestCost);
434
435 if (bestCost < oldCost) bestPosition = m;
436 }
437 }
438
439 return bestPosition;
440 }
441
442 std::vector<long long int> ComputeOptimalQubitsMap(
443 const std::vector<std::shared_ptr<Circuits::Circuit<>>>& layers,
444 int nrShuffles = 25, int nrSwaps = 10) {
445 const IndexType nrQubits = getNrQubits();
446
447 if (layers.empty() || nrQubits <= 2) return qubitsMap;
448
449
450 auto evaluateCost =
451 [&, this](const std::vector<IndexType>& candidateMap) -> double {
452 auto saveQubitsMap = qubitsMap;
453 auto saveQubitsMapInv = qubitsMapInv;
454 auto saveCurrentBondDim = currentBondDim;
455 auto saveTotalSwappingCost = totalSwappingCost;
456 auto saveBondCost = bondCost;
457
458 SetInitialQubitsMap(candidateMap);
459
460 for (const auto& layer : layers)
461 ApplyGates(layer->GetOperations());
462 auto cost = getTotalSwappingCost();
463 qubitsMap = std::move(saveQubitsMap);
464 qubitsMapInv = std::move(saveQubitsMapInv);
465 currentBondDim = std::move(saveCurrentBondDim);
466 totalSwappingCost = saveTotalSwappingCost;
467 bondCost = std::move(saveBondCost);
468
469 return cost;
470 };
471
472 auto evaluateCostBounded = [&, this](
473 const std::vector<IndexType>& candidateMap,
474 double bound) -> double {
475 auto saveQubitsMap = qubitsMap;
476 auto saveQubitsMapInv = qubitsMapInv;
477 auto saveCurrentBondDim = currentBondDim;
478 auto saveTotalSwappingCost = totalSwappingCost;
479 auto saveBondCost = bondCost;
480
481 SetInitialQubitsMap(candidateMap);
482
483 for (const auto& layer : layers) {
484 ApplyGates(layer->GetOperations());
485 auto cost = getTotalSwappingCost();
486 if (cost >= bound) {
487 // restore state
488 qubitsMap = std::move(saveQubitsMap);
489 qubitsMapInv = std::move(saveQubitsMapInv);
490 currentBondDim = std::move(saveCurrentBondDim);
491 totalSwappingCost = saveTotalSwappingCost;
492 bondCost = std::move(saveBondCost);
493 return cost;
494 }
495 }
496 auto cost = getTotalSwappingCost();
497 // restore state
498 qubitsMap = std::move(saveQubitsMap);
499 qubitsMapInv = std::move(saveQubitsMapInv);
500 currentBondDim = std::move(saveCurrentBondDim);
501 totalSwappingCost = saveTotalSwappingCost;
502 bondCost = std::move(saveBondCost);
503
504 return cost;
505 };
506
507
508 // Collect 2-qubit pairs from each layer, preserving layer boundaries
509 struct QubitPair {
510 IndexType q1, q2;
511 };
512 std::vector<std::vector<QubitPair>> layerPairs;
513
514 for (size_t li = 0; li < layers.size(); ++li) {
515 const auto& layer = layers[li];
516 std::vector<QubitPair> lp;
517 for (const auto& op : layer->GetOperations()) {
518 auto qbits = op->AffectedQubits();
519 if (qbits.size() >= 2) {
520 if (qbits[0] > qbits[1]) std::swap(qbits[0], qbits[1]);
521
522 lp.push_back({static_cast<IndexType>(qbits[0]),
523 static_cast<IndexType>(qbits[1])});
524 }
525 }
526 if (!lp.empty()) layerPairs.push_back(std::move(lp));
527 }
528 if (layerPairs.empty()) return qubitsMap;
529
530 // Build a linear qubit chain from ordered pairs using a group-merging
531 // strategy. Pairs are processed in order (earlier = more important).
532 // - If neither qubit is placed, create a new group [q1, q2].
533 // - If one qubit is already in a group, add the other at the group
534 // end (front or back) that minimizes distance to the existing one.
535 // - If both are in different groups, merge as [gA][gB] or [gB][gA],
536 // whichever places q1 and q2 closer together.
537 // - If both are already in the same group, do nothing.
538 auto buildChain = [&](const std::vector<std::vector<QubitPair>>& orderedLP)
539 -> std::vector<long long int> {
540 std::vector<std::deque<IndexType>> groups;
541 std::vector<int> qubitGroup(nrQubits, -1);
542 std::unordered_set<IndexType> placedQubits;
543
544 for (const auto& lp : orderedLP) {
545 for (const auto& p : lp) {
546 const int g1 = qubitGroup[p.q1];
547 const int g2 = qubitGroup[p.q2];
548
549 if (g1 < 0 && g2 < 0) {
550 // Neither placed: create a new group
551 const int gIdx = static_cast<int>(groups.size());
552 groups.push_back({p.q1, p.q2});
553 qubitGroup[p.q1] = gIdx;
554 qubitGroup[p.q2] = gIdx;
555 placedQubits.insert(p.q1);
556 placedQubits.insert(p.q2);
557 } else if (g1 >= 0 && g2 < 0) {
558 // q1 placed, q2 not: add q2 at the closer end of q1's group
559 auto& grp = groups[g1];
560 size_t idx = 0;
561 for (size_t i = 0; i < grp.size(); ++i)
562 if (grp[i] == p.q1) {
563 idx = i;
564 break;
565 }
566 // front distance: idx + 1, back distance: grp.size() - idx
567 if (idx + 1 <= grp.size() - idx)
568 grp.push_front(p.q2);
569 else
570 grp.push_back(p.q2);
571 qubitGroup[p.q2] = g1;
572 placedQubits.insert(p.q2);
573 } else if (g1 < 0 && g2 >= 0) {
574 // q2 placed, q1 not: add q1 at the closer end of q2's group
575 auto& grp = groups[g2];
576 size_t idx = 0;
577 for (size_t i = 0; i < grp.size(); ++i)
578 if (grp[i] == p.q2) {
579 idx = i;
580 break;
581 }
582 if (idx + 1 <= grp.size() - idx)
583 grp.push_front(p.q1);
584 else
585 grp.push_back(p.q1);
586 qubitGroup[p.q1] = g2;
587 placedQubits.insert(p.q1);
588 } else if (g1 != g2) {
589 // Both in different groups: merge to minimize distance
590 auto& grpA = groups[g1];
591 auto& grpB = groups[g2];
592 size_t idxA = 0, idxB = 0;
593 for (size_t i = 0; i < grpA.size(); ++i)
594 if (grpA[i] == p.q1) {
595 idxA = i;
596 break;
597 }
598 for (size_t i = 0; i < grpB.size(); ++i)
599 if (grpB[i] == p.q2) {
600 idxB = i;
601 break;
602 }
603 // [A][B]: q1 at idxA, q2 at |A|+idxB -> dist = |A|+idxB-idxA
604 // [B][A]: q2 at idxB, q1 at |B|+idxA -> dist = |B|+idxA-idxB
605 const size_t distAB = grpA.size() + idxB - idxA;
606 const size_t distBA = grpB.size() + idxA - idxB;
607
608 int mergedGroup;
609 if (distAB <= distBA) {
610 for (const auto q : grpB) grpA.push_back(q);
611 grpB.clear();
612 mergedGroup = g1;
613 } else {
614 for (const auto q : grpA) grpB.push_back(q);
615 grpA.clear();
616 mergedGroup = g2;
617 }
618 for (const auto q : groups[mergedGroup])
619 qubitGroup[q] = mergedGroup;
620 }
621 // Both in the same group: do nothing
622
623 if (placedQubits.size() == static_cast<size_t>(nrQubits))
624 break; // all qubits placed, can stop processing pairs
625 }
626
627 if (placedQubits.size() == static_cast<size_t>(nrQubits))
628 break; // all qubits placed, can stop processing pairs
629 }
630
631 // Concatenate all non-empty groups
632 std::vector<IndexType> chain;
633 chain.reserve(nrQubits);
634 for (const auto& grp : groups)
635 for (const auto q : grp) chain.push_back(q);
636
637 // Append any remaining unplaced qubits
638 for (IndexType q = 0; q < nrQubits; ++q)
639 if (qubitGroup[q] < 0) chain.push_back(q);
640
641 assert(chain.size() == static_cast<size_t>(nrQubits));
642
643 // Convert chain to qubitsMap: chain[physPos] = logicalQubit
644 std::vector<long long int> result(nrQubits);
645 for (size_t i = 0; i < chain.size(); ++i)
646 result[chain[i]] = static_cast<long long int>(i);
647
648 return result;
649 };
650
651 auto optMap = buildChain(layerPairs);
652 if (layers.size() <= 2) return optMap;
653 auto optCost = evaluateCost(optMap);
654
655 std::vector<long long int> qubitsMap(nrQubits);
656 std::iota(qubitsMap.begin(), qubitsMap.end(), 0);
657 double tryCost = evaluateCostBounded(qubitsMap, optCost);
658 if (tryCost < optCost) {
659 optCost = tryCost;
660 optMap = qubitsMap;
661 }
662
663 // try some random shuffles as well, in case the heuristic ordering is not
664 // optimal
665 std::mt19937 rng(42);
666 for (int i = 0; i < nrShuffles; ++i) {
667 std::shuffle(qubitsMap.begin(), qubitsMap.end(), rng);
668 tryCost = evaluateCostBounded(qubitsMap, optCost);
669 if (tryCost < optCost) {
670 optCost = tryCost;
671 optMap = qubitsMap;
672 }
673 }
674
675
676 std::uniform_int_distribution<IndexType> qubitDist(0, nrQubits - 1);
677 std::uniform_int_distribution<int> nrSwapsDist(
678 1, std::min<int>(3, static_cast<int>(nrQubits) / 2));
679
680 const int maxNoImprove = std::max(nrShuffles, static_cast<int>(nrQubits));
681 const int maxTotalShuffles = maxNoImprove * 3;
682 int noImproveCount = 0;
683
684 for (int s = 0; s < maxTotalShuffles && noImproveCount < maxNoImprove; ++s) {
685 auto tryMap = optMap;
686 const int nrSwaps = nrSwapsDist(rng);
687 for (int sw = 0; sw < nrSwaps; ++sw) {
688 const IndexType a = qubitDist(rng);
689 IndexType b = qubitDist(rng);
690 while (b == a) b = qubitDist(rng);
691 std::swap(tryMap[a], tryMap[b]);
692 }
693
694 auto cost = evaluateCostBounded(tryMap, optCost);
695 if (cost < optCost) {
696 optMap = tryMap;
697 optCost = cost;
698 noImproveCount = 0;
699 } else {
700 ++noImproveCount;
701 }
702 }
703
704 // 2-opt local search: iteratively swap pairs of positions in the best map
705 // and keep improvements, until no single swap can reduce the cost
706 {
707 auto candidate = optMap;
708 bool improved = true;
709 for (int improvementCount = 0; improved && improvementCount < nrSwaps;
710 ++improvementCount) {
711 improved = false;
712 for (IndexType i = 0; i < nrQubits; ++i) {
713 for (IndexType j = i + 1; j < nrQubits; ++j) {
714 // swap the mapped positions of qubits i and j
715 std::swap(candidate[i], candidate[j]);
716 auto cost = evaluateCostBounded(candidate, optCost);
717 if (cost < optCost) {
718 optMap = candidate;
719 optCost = cost;
720 improved = true;
721 } else {
722 // revert
723 std::swap(candidate[i], candidate[j]);
724 }
725 }
726 }
727 }
728 }
729
730 return optMap;
731 }
732
733 private:
734 void InitQubitsMap() {
735 qubitsMap.resize(getNrQubits());
736 qubitsMapInv.resize(getNrQubits());
737
738 for (IndexType i = 0; i < static_cast<IndexType>(getNrQubits()); ++i)
739 qubitsMapInv[i] = qubitsMap[i] = i;
740
741 totalSwappingCost = 0;
742 if (!currentBondDim.empty())
743 std::fill(currentBondDim.begin(), currentBondDim.end(), 1.0);
744 }
745
746 struct LightweightInitTag {};
747
748 // Lightweight constructor: sets up qubit maps but skips the expensive
749 // SetMaxBondDimension computation. Caller must populate bond arrays.
750 MPSDummySimulator(size_t N, LightweightInitTag)
751 : nrQubits(N) {
752 }
753
754 void SwapQubits(IndexType qubit1, IndexType qubit2) {
755 IndexType realq1 = qubitsMap[qubit1];
756 IndexType realq2 = qubitsMap[qubit2];
757 if (realq1 > realq2) {
758 std::swap(realq1, realq2);
759 std::swap(qubit1, qubit2);
760 }
761
762 if (realq2 - realq1 <= 1) return;
763
764 const IndexType meetPos = ComputeHeuristicMeetPosition(realq1, realq2);
765 SwapQubitsToPosition(qubit1, qubit2, meetPos);
766 }
767
768 // Pick the meeting position with the lowest bond cost (i.e., lowest
769 // current bond dimension), matching the real MPS simulator's
770 // FindBestMeetingPositionLocal behavior.
771 IndexType ComputeHeuristicMeetPosition(IndexType realq1,
772 IndexType realq2) const {
773 assert(realq1 < realq2);
774
775 IndexType bestPos = realq1;
776 double bestCost = bondCost[realq1];
777
778 for (IndexType m = realq1 + 1; m < realq2; ++m) {
779 if (bondCost[m] < bestCost) {
780 bestCost = bondCost[m];
781 bestPos = m;
782 }
783 }
784
785 return bestPos;
786 }
787
788 public:
789 // Swap two logical qubits so they meet at a specified bond position.
790 // meetPosition is the bond index (in real/chain coordinates) where
791 // the two qubits will end up adjacent: one at meetPosition, the other
792 // at meetPosition+1.
793 // meetPosition must be in [min(realq1,realq2), max(realq1,realq2)-1].
795 IndexType meetPosition) {
796 IndexType realq1 = qubitsMap[qubit1];
797 IndexType realq2 = qubitsMap[qubit2];
798 if (realq1 > realq2) {
799 std::swap(realq1, realq2);
800 std::swap(qubit1, qubit2);
801 }
802
803 if (realq2 - realq1 <= 1) return;
804
805 assert(meetPosition >= realq1 && meetPosition < realq2);
806
807 // Move lower qubit (qubit1) rightward from realq1 to meetPosition
808 {
809 IndexType movingReal = realq1;
810 while (movingReal < meetPosition) {
811 const IndexType toReal = movingReal + 1;
812 const IndexType toInv = qubitsMapInv[toReal];
813
814 qubitsMap[toInv] = movingReal;
815 qubitsMapInv[movingReal] = toInv;
816
817 qubitsMap[qubit1] = toReal;
818 qubitsMapInv[toReal] = qubit1;
819
820 totalSwappingCost += bondCost[movingReal];
821 growBondDimension(movingReal, true);
822 movingReal = toReal;
823 }
824 }
825
826 // Move upper qubit (qubit2) leftward from realq2 to meetPosition+1
827 {
828 IndexType movingReal = realq2;
829 while (movingReal > meetPosition + 1) {
830 const IndexType toReal = movingReal - 1;
831 const IndexType toInv = qubitsMapInv[toReal];
832
833 qubitsMap[toInv] = movingReal;
834 qubitsMapInv[movingReal] = toInv;
835
836 qubitsMap[qubit2] = toReal;
837 qubitsMapInv[toReal] = qubit2;
838
839 totalSwappingCost += bondCost[toReal];
840 growBondDimension(toReal, true);
841 movingReal = toReal;
842 }
843 }
844
845 assert(abs(qubitsMap[qubit1] - qubitsMap[qubit2]) == 1);
846 }
847
848 private:
849 size_t nrQubits;
850
851 std::vector<IndexType> qubitsMap;
852 std::vector<IndexType> qubitsMapInv;
853
854 static constexpr size_t physExtent = 2;
855 IndexType maxVirtualExtent = 0;
856 std::vector<double> bondCost;
857 std::vector<double> maxBondDim;
858 std::vector<double> currentBondDim;
859
860 double totalSwappingCost = 0;
861
862 double growthFactorSwap = 1.;
863 double growthFactorGate = 0.65;
864
865 void growBondDimension(IndexType bond, bool swap = true) {
866 // the left and right bond dimensions are relevant because:
867 // the initial configuration before applying the swap or other gate is:
868
869 // - O - O -
870 // | |
871
872 // The two physical legs have dimension 2
873 // this is contracted into:
874
875 // ---
876 // -| |-
877 // ---
878 // | |
879
880 // the left and right dimensions stay the same and also the physical legs have dimension 2
881
882 // then the swap or the other gate is applied, getting a result that looks graphically as above, but of course with different values inside the tensor
883 // swap is special, just swaps the values for (0, 1) and (1, 0) in
884 // the physical legs, while other gates can change all values in the tensor
885
886 // then the tensor is reshaped into a matrix, having dimensions 2 * leftDim x 2 * rightNeighborDim on this matrix SVD is applied, to separate out the
887 // qubits tensors again, and the bond dimension is the number of singular
888 // values kept after truncation (if done), or the number of non-zero
889 // singular values if no truncation is done. The bond dimension can be at
890 // most min(2 * min(leftDim, rightNeighborDim), maxBondDim[bond]) and the minimum is obviously 1
891
892
893 const IndexType leftBond = bond - 1;
894 const IndexType rightNeigborBond = bond + 1;
895 const double betweenDim = currentBondDim[bond];
896
897 const double leftDim = leftBond >= 0 ? currentBondDim[leftBond] : 1;
898 const double rightNeighborDim = rightNeigborBond < static_cast<IndexType>(currentBondDim.size()) ? currentBondDim[rightNeigborBond] : 1;
899
900 const double newMaxDim = (swap && leftDim == rightNeighborDim) ? betweenDim : 2. * std::min(leftDim, rightNeighborDim);
901
902 const double growthFactor = swap ? growthFactorSwap : growthFactorGate;
903
904 currentBondDim[bond] = std::min(newMaxDim * growthFactor, maxBondDim[bond]);
905
906 currentBondDim[bond] = std::max(currentBondDim[bond], 1.);
907
908 bondCost[bond] =
909 currentBondDim[bond] * currentBondDim[bond] * currentBondDim[bond];
910 }
911};
912
913} // namespace Simulators
914
915#endif
Circuit class for holding the sequence of operations.
Definition Circuit.h:46
The operation interface.
Definition Operations.h:358
void SetCurrentBondDimensions(const std::vector< double > &dims)
const std::vector< IndexType > & getQubitsMap() const
void setGrowthFactorSwap(double factor)
const std::vector< IndexType > & getQubitsMapInv() const
void ApplyGates(const std::vector< QC::Gates::AppliedGate< MatrixClass > > &gates)
const std::vector< double > & getCurrentBondDimensions() const
void ApplyGate(const QC::Gates::AppliedGate< MatrixClass > &gate)
std::unique_ptr< MPSDummySimulator > Clone() const
void SwapQubitsToPosition(IndexType qubit1, IndexType qubit2, IndexType meetPosition)
void EvaluateMeetingPositionCost(IndexType meetPosition, const std::vector< std::shared_ptr< Circuits::IOperation<> > > &upcomingGates, long long int currentGateIndex, int lookaheadDepth, int lookaheadDepthWithHeuristic, double currentCost, double &bestCost, bool useSameDummy=false)
QC::TensorNetworks::MPSSimulatorInterface::GateClass GateClass
std::vector< long long int > ComputeOptimalQubitsMap(const std::vector< std::shared_ptr< Circuits::Circuit<> > > &layers, int nrShuffles=25, int nrSwaps=10)
QC::TensorNetworks::MPSSimulatorInterface::MatrixClass MatrixClass
const std::vector< double > & getMaxBondDimensions() const
void ApplyGates(const std::vector< std::shared_ptr< Circuits::IOperation<> > > &gates)
IndexType FindBestMeetingPosition(const std::vector< std::shared_ptr< Circuits::IOperation<> > > &upcomingGates, long long int currentGateIndex, int lookaheadDepth, int lookaheadDepthWithHeuristic, double currentCost, double &bestCost)
void ApplyGate(const GateClass &gate, IndexType qubit, IndexType controllingQubit1=0)
void setGrowthFactorGate(double factor)
void SetInitialQubitsMap(const std::vector< long long int > &initialMap)
void SetMaxBondDimension(IndexType val)
void ApplyGate(const std::shared_ptr< Circuits::IOperation<> > &gate)
@ kGate
the usual quantum gate, result stays in simulator's state
Definition Operations.h:28