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