The Quantum Exact Simulation Toolkit v4.2.0
Loading...
Searching...
No Matches
multiplication.cpp
1/** @file
2 * API definitions for directly pre- and post-multiplying
3 * operators upon density matrices, likely constituting
4 * non-physical operations which break state normalisation.
5 *
6 * @author Tyson Jones
7 */
8
9#include "quest/include/qureg.h"
10#include "quest/include/paulis.h"
11#include "quest/include/matrices.h"
12#include "quest/include/multiplication.h"
13
14#include "quest/src/core/validation.hpp"
15#include "quest/src/core/utilities.hpp"
16#include "quest/src/core/localiser.hpp"
17#include "quest/src/core/paulilogic.hpp"
18
19#include <vector>
20
21using std::vector;
22
23
24
25/*
26 * CompMatr1
27 */
28
29extern "C" {
30
31void leftapplyCompMatr1(Qureg qureg, int target, CompMatr1 matrix) {
32 validate_quregFields(qureg, __func__);
33 validate_target(qureg, target, __func__);
34 validate_matrixFields(matrix, __func__);
35
36 bool conj = false;
37 bool transp = false;
38 localiser_statevec_anyCtrlOneTargDenseMatr(qureg, {}, {}, target, matrix, conj, transp);
39}
40
41void rightapplyCompMatr1(Qureg qureg, int target, CompMatr1 matrix) {
42 validate_quregFields(qureg, __func__);
43 validate_quregIsDensityMatrix(qureg, __func__);
44 validate_target(qureg, target, __func__);
45 validate_matrixFields(matrix, __func__);
46
47 // rho matrix ~ transpose(rho) (x) I ||rho>>
48 bool conj = false;
49 bool transp = true;
50 int qubit = util_getBraQubit(target, qureg);
51 localiser_statevec_anyCtrlOneTargDenseMatr(qureg, {}, {}, qubit, matrix, conj, transp);
52}
53
54} // end de-mangler
55
56
57
58/*
59 * CompMatr2
60 */
61
62extern "C" {
63
64void leftapplyCompMatr2(Qureg qureg, int target1, int target2, CompMatr2 matrix) {
65 validate_quregFields(qureg, __func__);
66 validate_twoTargets(qureg, target1, target2, __func__);
67 validate_matrixFields(matrix, __func__);
68 validate_mixedAmpsFitInNode(qureg, 2, __func__);
69
70 bool conj = false;
71 bool transp = false;
72 localiser_statevec_anyCtrlTwoTargDenseMatr(qureg, {}, {}, target1, target2, matrix, conj, transp);
73}
74
75void rightapplyCompMatr2(Qureg qureg, int target1, int target2, CompMatr2 matrix) {
76 validate_quregFields(qureg, __func__);
77 validate_quregIsDensityMatrix(qureg, __func__);
78 validate_twoTargets(qureg, target1, target2, __func__);
79 validate_matrixFields(matrix, __func__);
80 validate_mixedAmpsFitInNode(qureg, 2, __func__);
81
82 // rho matrix ~ transpose(rho) (x) I ||rho>>
83 bool conj = false;
84 bool transp = true;
85 int qubit1 = util_getBraQubit(target1, qureg);
86 int qubit2 = util_getBraQubit(target2, qureg);
87 localiser_statevec_anyCtrlTwoTargDenseMatr(qureg, {}, {}, qubit1, qubit2, matrix, conj, transp);
88}
89
90} // end de-mangler
91
92
93
94/*
95 * CompMatr
96 */
97
98extern "C" {
99
100void leftapplyCompMatr(Qureg qureg, int* targets, int numTargets, CompMatr matrix) {
101 validate_quregFields(qureg, __func__);
102 validate_targets(qureg, targets, numTargets, __func__);
103 validate_matrixDimMatchesTargets(matrix, numTargets, __func__); // also validates fields and is-sync
104 validate_mixedAmpsFitInNode(qureg, numTargets, __func__);
105
106 bool conj = false;
107 bool transp = false;
108 localiser_statevec_anyCtrlAnyTargDenseMatr(qureg, {}, {}, util_getVector(targets, numTargets), matrix, conj, transp);
109}
110
111void rightapplyCompMatr(Qureg qureg, int* targets, int numTargets, CompMatr matrix) {
112 validate_quregFields(qureg, __func__);
113 validate_quregIsDensityMatrix(qureg, __func__);
114 validate_targets(qureg, targets, numTargets, __func__);
115 validate_matrixDimMatchesTargets(matrix, numTargets, __func__); // also validates fields and is-sync
116 validate_mixedAmpsFitInNode(qureg, numTargets, __func__);
117
118 // rho matrix ~ transpose(rho) (x) I ||rho>>
119 bool conj = false;
120 bool transp = true;
121 auto qubits = util_getBraQubits(util_getVector(targets, numTargets), qureg);
122 localiser_statevec_anyCtrlAnyTargDenseMatr(qureg, {}, {}, qubits, matrix, conj, transp);
123}
124
125} // end de-mangler
126
127void leftapplyCompMatr(Qureg qureg, vector<int> targets, CompMatr matr) {
128
129 leftapplyCompMatr(qureg, targets.data(), targets.size(), matr);
130}
131
132void rightapplyCompMatr(Qureg qureg, vector<int> targets, CompMatr matr) {
133
134 rightapplyCompMatr(qureg, targets.data(), targets.size(), matr);
135}
136
137
138
139/*
140 * DiagMatr1
141 */
142
143extern "C" {
144
145void leftapplyDiagMatr1(Qureg qureg, int target, DiagMatr1 matrix) {
146 validate_quregFields(qureg, __func__);
147 validate_target(qureg, target, __func__);
148 validate_matrixFields(matrix, __func__);
149
150 bool conj = false;
151 localiser_statevec_anyCtrlOneTargDiagMatr(qureg, {}, {}, target, matrix, conj);
152}
153
154void rightapplyDiagMatr1(Qureg qureg, int target, DiagMatr1 matrix) {
155 validate_quregFields(qureg, __func__);
156 validate_quregIsDensityMatrix(qureg, __func__);
157 validate_target(qureg, target, __func__);
158 validate_matrixFields(matrix, __func__);
159
160 bool conj = false;
161 int qubit = util_getBraQubit(target, qureg);
162 localiser_statevec_anyCtrlOneTargDiagMatr(qureg, {}, {}, qubit, matrix, conj);
163}
164
165} // end de-mangler
166
167
168
169/*
170 * DiagMatr2
171 */
172
173extern "C" {
174
175void leftapplyDiagMatr2(Qureg qureg, int target1, int target2, DiagMatr2 matrix) {
176 validate_quregFields(qureg, __func__);
177 validate_twoTargets(qureg, target1, target2, __func__);
178 validate_matrixFields(matrix, __func__);
179
180 bool conj = false;
181 localiser_statevec_anyCtrlTwoTargDiagMatr(qureg, {}, {}, target1, target2, matrix, conj);
182}
183
184void rightapplyDiagMatr2(Qureg qureg, int target1, int target2, DiagMatr2 matrix) {
185 validate_quregFields(qureg, __func__);
186 validate_quregIsDensityMatrix(qureg, __func__);
187 validate_twoTargets(qureg, target1, target2, __func__);
188 validate_matrixFields(matrix, __func__);
189
190 bool conj = false;
191 int qubit1 = util_getBraQubit(target1, qureg);
192 int qubit2 = util_getBraQubit(target2, qureg);
193 localiser_statevec_anyCtrlTwoTargDiagMatr(qureg, {}, {}, qubit1, qubit2, matrix, conj);
194}
195
196} // end de-mangler
197
198
199
200/*
201 * DiagMatr
202 */
203
204extern "C" {
205
206void leftapplyDiagMatr(Qureg qureg, int* targets, int numTargets, DiagMatr matrix) {
207 validate_quregFields(qureg, __func__);
208 validate_targets(qureg, targets, numTargets, __func__);
209 validate_matrixDimMatchesTargets(matrix, numTargets, __func__); // also validates fields and is-sync
210
211 bool conj = false;
212 qcomp exponent = 1;
213 auto qubits = util_getVector(targets, numTargets);
214 localiser_statevec_anyCtrlAnyTargDiagMatr(qureg, {}, {}, qubits, matrix, exponent, conj);
215}
216
217void rightapplyDiagMatr(Qureg qureg, int* targets, int numTargets, DiagMatr matrix) {
218 validate_quregFields(qureg, __func__);
219 validate_quregIsDensityMatrix(qureg, __func__);
220 validate_targets(qureg, targets, numTargets, __func__);
221 validate_matrixDimMatchesTargets(matrix, numTargets, __func__); // also validates fields and is-sync
222
223 bool conj = false;
224 qcomp exponent = 1;
225 auto qubits = util_getBraQubits(util_getVector(targets, numTargets), qureg);
226 localiser_statevec_anyCtrlAnyTargDiagMatr(qureg, {}, {}, qubits, matrix, exponent, conj);
227}
228
229} // end de-mangler
230
231void leftapplyDiagMatr(Qureg qureg, vector<int> targets, DiagMatr matrix) {
232
233 leftapplyDiagMatr(qureg, targets.data(), targets.size(), matrix);
234}
235
236void rightapplyDiagMatr(Qureg qureg, vector<int> targets, DiagMatr matrix) {
237
238 rightapplyDiagMatr(qureg, targets.data(), targets.size(), matrix);
239}
240
241
242
243/*
244 * DiagMatrPower
245 */
246
247extern "C" {
248
249void leftapplyDiagMatrPower(Qureg qureg, int* targets, int numTargets, DiagMatr matrix, qcomp exponent) {
250 validate_quregFields(qureg, __func__);
251 validate_targets(qureg, targets, numTargets, __func__);
252 validate_matrixDimMatchesTargets(matrix, numTargets, __func__); // also validates fields and is-sync, but not unitarity
253 validate_matrixExpIsNonDiverging(matrix, exponent, __func__); // harmlessly re-validates fields and is-sync
254
255 bool conj = false;
256 auto qubits = util_getVector(targets, numTargets);
257 localiser_statevec_anyCtrlAnyTargDiagMatr(qureg, {}, {}, qubits, matrix, exponent, conj);
258}
259
260void rightapplyDiagMatrPower(Qureg qureg, int* targets, int numTargets, DiagMatr matrix, qcomp exponent) {
261 validate_quregFields(qureg, __func__);
262 validate_quregIsDensityMatrix(qureg, __func__);
263 validate_targets(qureg, targets, numTargets, __func__);
264 validate_matrixDimMatchesTargets(matrix, numTargets, __func__); // also validates fields and is-sync, but not unitarity
265 validate_matrixExpIsNonDiverging(matrix, exponent, __func__); // harmlessly re-validates fields and is-sync
266
267 bool conj = false;
268 auto qubits = util_getBraQubits(util_getVector(targets, numTargets), qureg);
269 localiser_statevec_anyCtrlAnyTargDiagMatr(qureg, {}, {}, qubits, matrix, exponent, conj);
270}
271
272} // end de-mangler
273
274void leftapplyDiagMatrPower(Qureg qureg, vector<int> targets, DiagMatr matrix, qcomp exponent) {
275
276 leftapplyDiagMatrPower(qureg, targets.data(), targets.size(), matrix, exponent);
277}
278
279void rightapplyDiagMatrPower(Qureg qureg, vector<int> targets, DiagMatr matrix, qcomp exponent) {
280
281 rightapplyDiagMatrPower(qureg, targets.data(), targets.size(), matrix, exponent);
282}
283
284
285
286/*
287 * FullStateDiagMatr (and power)
288 */
289
290extern "C" {
291
293 validate_quregFields(qureg, __func__);
294 validate_matrixFields(matrix, __func__);
295 validate_matrixAndQuregAreCompatible(matrix, qureg, false, __func__); // matrix can be non-unitary
296
297 leftapplyFullStateDiagMatrPower(qureg, matrix, 1); // harmlessly re-validates
298}
299
300void leftapplyFullStateDiagMatrPower(Qureg qureg, FullStateDiagMatr matrix, qcomp exponent) {
301 validate_quregFields(qureg, __func__);
302 validate_matrixFields(matrix, __func__);
303 validate_matrixAndQuregAreCompatible(matrix, qureg, false, __func__); // matrix can be non-unitary
304 validate_matrixExpIsNonDiverging(matrix, exponent, __func__);
305
306 // rho -> matrix^exponent rho
307 bool leftMultiply = true;
308 bool rightMultiply = false;
309 bool rightConj = false;
310
311 (qureg.isDensityMatrix)?
312 localiser_densmatr_allTargDiagMatr(qureg, matrix, exponent, leftMultiply, rightMultiply, rightConj):
313 localiser_statevec_allTargDiagMatr(qureg, matrix, exponent);
314}
315
317 validate_quregFields(qureg, __func__);
318 validate_quregIsDensityMatrix(qureg, __func__);
319 validate_matrixFields(matrix, __func__);
320 validate_matrixAndQuregAreCompatible(matrix, qureg, false, __func__); // matrix can be non-unitary
321
322 rightapplyFullStateDiagMatrPower(qureg, matrix, 1); // harmlessly re-validates
323}
324
325void rightapplyFullStateDiagMatrPower(Qureg qureg, FullStateDiagMatr matrix, qcomp exponent) {
326 validate_quregFields(qureg, __func__);
327 validate_quregIsDensityMatrix(qureg, __func__);
328 validate_matrixFields(matrix, __func__);
329 validate_matrixAndQuregAreCompatible(matrix, qureg, false, __func__); // matrix can be non-unitary
330 validate_matrixExpIsNonDiverging(matrix, exponent, __func__);
331
332 // rho -> rho matrix^exponent
333 bool leftMultiply = false;
334 bool rightMultiply = true;
335 bool rightConj = false;
336 localiser_densmatr_allTargDiagMatr(qureg, matrix, exponent, leftMultiply, rightMultiply, rightConj);
337}
338
339} // end de-mangler
340
341
342
343/*
344 * swap
345 */
346
347extern "C" {
348
349void leftapplySwap(Qureg qureg, int qubit1, int qubit2) {
350 validate_quregFields(qureg, __func__);
351 validate_twoTargets(qureg, qubit1, qubit2, __func__);
352
353 localiser_statevec_anyCtrlSwap(qureg, {}, {}, qubit1, qubit2);
354}
355
356void rightapplySwap(Qureg qureg, int qubit1, int qubit2) {
357 validate_quregFields(qureg, __func__);
358 validate_quregIsDensityMatrix(qureg, __func__);
359 validate_twoTargets(qureg, qubit1, qubit2, __func__);
360
361 qubit1 = util_getBraQubit(qubit1, qureg);
362 qubit2 = util_getBraQubit(qubit2, qureg);
363 localiser_statevec_anyCtrlSwap(qureg, {}, {}, qubit1, qubit2);
364}
365
366} // end de-mangler
367
368
369
370/*
371 * individual Paulis
372 */
373
374extern "C" {
375
376void leftapplyPauliX(Qureg qureg, int target) {
377 validate_quregFields(qureg, __func__);
378 validate_target(qureg, target, __func__);
379
380 PauliStr str = getPauliStr("X", {target});
381 localiser_statevec_anyCtrlPauliTensor(qureg, {}, {}, str);
382}
383
384void leftapplyPauliY(Qureg qureg, int target) {
385 validate_quregFields(qureg, __func__);
386 validate_target(qureg, target, __func__);
387
388 PauliStr str = getPauliStr("Y", {target});
389 localiser_statevec_anyCtrlPauliTensor(qureg, {}, {}, str);
390}
391
392void leftapplyPauliZ(Qureg qureg, int target) {
393 validate_quregFields(qureg, __func__);
394 validate_target(qureg, target, __func__);
395
396 PauliStr str = getPauliStr("Z", {target});
397 localiser_statevec_anyCtrlPauliTensor(qureg, {}, {}, str);
398}
399
400void rightapplyPauliX(Qureg qureg, int target) {
401 validate_quregFields(qureg, __func__);
402 validate_quregIsDensityMatrix(qureg, __func__);
403 validate_target(qureg, target, __func__);
404
405 PauliStr str = getPauliStr("X", {target});
406 str = paulis_getShiftedPauliStr(str, qureg.numQubits);
407 localiser_statevec_anyCtrlPauliTensor(qureg, {}, {}, str);
408}
409
410void rightapplyPauliY(Qureg qureg, int target) {
411 validate_quregFields(qureg, __func__);
412 validate_quregIsDensityMatrix(qureg, __func__);
413 validate_target(qureg, target, __func__);
414
415 qcomp factor = -1; // undo transpose
416 PauliStr str = getPauliStr("Y", {target});
417 str = paulis_getShiftedPauliStr(str, qureg.numQubits);
418 localiser_statevec_anyCtrlPauliTensor(qureg, {}, {}, str, factor);
419}
420
421void rightapplyPauliZ(Qureg qureg, int target) {
422 validate_quregFields(qureg, __func__);
423 validate_quregIsDensityMatrix(qureg, __func__);
424 validate_target(qureg, target, __func__);
425
426 PauliStr str = getPauliStr("Z", {target});
427 str = paulis_getShiftedPauliStr(str, qureg.numQubits);
428 localiser_statevec_anyCtrlPauliTensor(qureg, {}, {}, str);
429}
430
431} // end de-mangler
432
433
434
435/*
436 * Pauli strings
437 */
438
439extern "C" {
440
442 validate_quregFields(qureg, __func__);
443 validate_pauliStrTargets(qureg, str, __func__);
444
445 localiser_statevec_anyCtrlPauliTensor(qureg, {}, {}, str);
446}
447
449 validate_quregFields(qureg, __func__);
450 validate_quregIsDensityMatrix(qureg, __func__);
451 validate_pauliStrTargets(qureg, str, __func__);
452
453 qcomp factor = paulis_getSignOfPauliStrConj(str); // undo transpose
454 str = paulis_getShiftedPauliStr(str, qureg.numQubits);
455 localiser_statevec_anyCtrlPauliTensor(qureg, {}, {}, str, factor);
456}
457
458} // end de-mangler
459
460
461
462/*
463 * Pauli gadgets
464 */
465
466extern "C" {
467
468void leftapplyPauliGadget(Qureg qureg, PauliStr str, qreal angle) {
469 validate_quregFields(qureg, __func__);
470 validate_pauliStrTargets(qureg, str, __func__);
471
472 qreal phase = util_getPhaseFromGateAngle(angle);
473 localiser_statevec_anyCtrlPauliGadget(qureg, {}, {}, str, phase);
474}
475
476void rightapplyPauliGadget(Qureg qureg, PauliStr str, qreal angle) {
477 validate_quregFields(qureg, __func__);
478 validate_quregIsDensityMatrix(qureg, __func__);
479 validate_pauliStrTargets(qureg, str, __func__);
480
481 qreal factor = paulis_getSignOfPauliStrConj(str);
482 qreal phase = factor * util_getPhaseFromGateAngle(angle);
483 str = paulis_getShiftedPauliStr(str, qureg.numQubits);
484 localiser_statevec_anyCtrlPauliGadget(qureg, {}, {}, str, phase);
485}
486
487} // end de-mangler
488
489
490
491/*
492 * phase gadgets
493 */
494
495extern "C" {
496
497void leftapplyPhaseGadget(Qureg qureg, int* targets, int numTargets, qreal angle) {
498 validate_quregFields(qureg, __func__);
499 validate_targets(qureg, targets, numTargets, __func__);
500
501 qreal phase = util_getPhaseFromGateAngle(angle);
502 auto qubits = util_getVector(targets, numTargets);
503 localiser_statevec_anyCtrlPhaseGadget(qureg, {}, {}, qubits, phase);
504}
505
506void rightapplyPhaseGadget(Qureg qureg, int* targets, int numTargets, qreal angle) {
507 validate_quregFields(qureg, __func__);
508 validate_quregIsDensityMatrix(qureg, __func__);
509 validate_targets(qureg, targets, numTargets, __func__);
510
511 qreal phase = util_getPhaseFromGateAngle(angle);
512 auto qubits = util_getBraQubits(util_getVector(targets, numTargets), qureg);
513 localiser_statevec_anyCtrlPhaseGadget(qureg, {}, {}, qubits, phase);
514}
515
516} // end de-mangler
517
518void leftapplyPhaseGadget(Qureg qureg, vector<int> targets, qreal angle) {
519
520 leftapplyPhaseGadget(qureg, targets.data(), targets.size(), angle);
521}
522
523void rightapplyPhaseGadget(Qureg qureg, vector<int> targets, qreal angle) {
524
525 rightapplyPhaseGadget(qureg, targets.data(), targets.size(), angle);
526}
527
528
529
530/*
531 * many-qubit NOTs
532 */
533
534extern "C" {
535
536void leftapplyMultiQubitNot(Qureg qureg, int* targets, int numTargets) {
537 validate_quregFields(qureg, __func__);
538 validate_targets(qureg, targets, numTargets, __func__);
539
540 // harmlessly re-validates
541 PauliStr str = getPauliStr(std::string(numTargets, 'X'), targets, numTargets);
542 leftapplyPauliStr(qureg, str);
543}
544
545void rightapplyMultiQubitNot(Qureg qureg, int* targets, int numTargets) {
546 validate_quregFields(qureg, __func__);
547 validate_quregIsDensityMatrix(qureg, __func__);
548 validate_targets(qureg, targets, numTargets, __func__);
549
550 // harmlessly re-validates
551 PauliStr str = getPauliStr(std::string(numTargets, 'X'), targets, numTargets);
552 rightapplyPauliStr(qureg, str);
553}
554
555} // end de-mangler
556
557void leftapplyMultiQubitNot(Qureg qureg, vector<int> targets) {
558
559 leftapplyMultiQubitNot(qureg, targets.data(), targets.size());
560}
561
562void rightapplyMultiQubitNot(Qureg qureg, vector<int> targets) {
563
564 rightapplyMultiQubitNot(qureg, targets.data(), targets.size());
565}
566
567
568
569/*
570 * projectors
571 */
572
573extern "C" {
574
575void leftapplyQubitProjector(Qureg qureg, int qubit, int outcome) {
576 validate_quregFields(qureg, __func__);
577 validate_target(qureg, qubit, __func__);
578 validate_measurementOutcomeIsValid(outcome, __func__);
579
580 qreal prob = 1;
581 localiser_statevec_multiQubitProjector(qureg, {qubit}, {outcome}, prob);
582}
583
584void leftapplyMultiQubitProjector(Qureg qureg, int* qubits, int* outcomes, int numQubits) {
585 validate_quregFields(qureg, __func__);
586 validate_targets(qureg, qubits, numQubits, __func__);
587 validate_measurementOutcomesAreValid(outcomes, numQubits, __func__);
588
589 qreal prob = 1;
590 auto qubitVec = util_getVector(qubits, numQubits);
591 auto outcomeVec = util_getVector(outcomes, numQubits);
592 localiser_statevec_multiQubitProjector(qureg, qubitVec, outcomeVec, prob);
593}
594
595void rightapplyQubitProjector(Qureg qureg, int qubit, int outcome) {
596 validate_quregFields(qureg, __func__);
597 validate_quregIsDensityMatrix(qureg, __func__);
598 validate_target(qureg, qubit, __func__);
599 validate_measurementOutcomeIsValid(outcome, __func__);
600
601 qreal prob = 1;
602 localiser_statevec_multiQubitProjector(qureg, {util_getBraQubit(qubit,qureg)}, {outcome}, prob);
603}
604
605void rightapplyMultiQubitProjector(Qureg qureg, int* qubits, int* outcomes, int numQubits) {
606 validate_quregFields(qureg, __func__);
607 validate_quregIsDensityMatrix(qureg, __func__);
608 validate_targets(qureg, qubits, numQubits, __func__);
609 validate_measurementOutcomesAreValid(outcomes, numQubits, __func__);
610
611 qreal prob = 1;
612 auto qubitVec = util_getBraQubits(util_getVector(qubits, numQubits), qureg);
613 auto outcomeVec = util_getVector(outcomes, numQubits);
614 localiser_statevec_multiQubitProjector(qureg, qubitVec, outcomeVec, prob);
615}
616
617} // end de-mangler
618
619void leftapplyMultiQubitProjector(Qureg qureg, vector<int> qubits, vector<int> outcomes) {
620 validate_measurementOutcomesMatchTargets(qubits.size(), outcomes.size(), __func__);
621
622 leftapplyMultiQubitProjector(qureg, qubits.data(), outcomes.data(), outcomes.size());
623}
624
625void rightapplyMultiQubitProjector(Qureg qureg, vector<int> qubits, vector<int> outcomes) {
626 validate_measurementOutcomesMatchTargets(qubits.size(), outcomes.size(), __func__);
627
628 rightapplyMultiQubitProjector(qureg, qubits.data(), outcomes.data(), outcomes.size());
629}
630
631
632
633/*
634 * Pauli string sums
635 */
636
637extern "C" {
638
639void leftapplyPauliStrSum(Qureg qureg, PauliStrSum sum, Qureg workspace) {
640 validate_quregFields(qureg, __func__);
641 validate_quregFields(workspace, __func__);
642 validate_quregCanBeWorkspace(qureg, workspace, __func__);
643 validate_pauliStrSumFields(sum, __func__);
644 validate_pauliStrSumTargets(sum, qureg, __func__);
645
646 // clone qureg to workspace, set qureg to blank
647 localiser_statevec_setQuregToClone(workspace, qureg);
648 localiser_statevec_initUniformState(qureg, 0);
649
650 // left-multiply each term in-turn, mixing into output qureg, then undo using idempotency
651 for (qindex i=0; i<sum.numTerms; i++) {
652 localiser_statevec_anyCtrlPauliTensor(workspace, {}, {}, sum.strings[i]);
653 localiser_statevec_setQuregToWeightedSum(qureg, {1, sum.coeffs[i]}, {qureg, workspace});
654 localiser_statevec_anyCtrlPauliTensor(workspace, {}, {}, sum.strings[i]);
655 }
656
657 // workspace -> qureg, and qureg -> sum * qureg
658}
659
660void rightapplyPauliStrSum(Qureg qureg, PauliStrSum sum, Qureg workspace) {
661 validate_quregFields(qureg, __func__);
662 validate_quregFields(workspace, __func__);
663 validate_quregIsDensityMatrix(qureg, __func__);
664 validate_quregCanBeWorkspace(qureg, workspace, __func__);
665 validate_pauliStrSumFields(sum, __func__);
666 validate_pauliStrSumTargets(sum, qureg, __func__);
667
668 // clone qureg to workspace, set qureg to blank
669 localiser_statevec_setQuregToClone(workspace, qureg);
670 localiser_statevec_initUniformState(qureg, 0);
671
672 // post-multiply each term in-turn, mixing into output qureg, then undo using idempotency
673 for (qindex i=0; i<sum.numTerms; i++) {
674 PauliStr str = paulis_getShiftedPauliStr(sum.strings[i], qureg.numQubits);
675 qcomp factor = paulis_getSignOfPauliStrConj(str); // undoes transpose
676
677 localiser_statevec_anyCtrlPauliTensor(workspace, {}, {}, str, factor);
678 localiser_statevec_setQuregToWeightedSum(qureg, {1, sum.coeffs[i]}, {qureg, workspace});
679 localiser_statevec_anyCtrlPauliTensor(workspace, {}, {}, str, factor);
680 }
681
682 // workspace -> qureg, and qureg -> sum * qureg
683}
684
685} // end de-mangler
void rightapplyCompMatr1(Qureg qureg, int target, CompMatr1 matrix)
void leftapplyCompMatr1(Qureg qureg, int target, CompMatr1 matrix)
void rightapplyCompMatr2(Qureg qureg, int target1, int target2, CompMatr2 matrix)
void leftapplyCompMatr2(Qureg qureg, int target1, int target2, CompMatr2 matrix)
void leftapplyCompMatr(Qureg qureg, int *targets, int numTargets, CompMatr matrix)
void rightapplyCompMatr(Qureg qureg, int *targets, int numTargets, CompMatr matrix)
void leftapplyDiagMatr1(Qureg qureg, int target, DiagMatr1 matrix)
void rightapplyDiagMatr1(Qureg qureg, int target, DiagMatr1 matrix)
void leftapplyDiagMatr2(Qureg qureg, int target1, int target2, DiagMatr2 matrix)
void rightapplyDiagMatr2(Qureg qureg, int target1, int target2, DiagMatr2 matrix)
void rightapplyDiagMatr(Qureg qureg, int *targets, int numTargets, DiagMatr matrix)
void leftapplyDiagMatrPower(Qureg qureg, int *targets, int numTargets, DiagMatr matrix, qcomp exponent)
void rightapplyDiagMatrPower(Qureg qureg, int *targets, int numTargets, DiagMatr matrix, qcomp exponent)
void leftapplyDiagMatr(Qureg qureg, int *targets, int numTargets, DiagMatr matrix)
void rightapplyFullStateDiagMatr(Qureg qureg, FullStateDiagMatr matrix)
void leftapplyFullStateDiagMatr(Qureg qureg, FullStateDiagMatr matrix)
void rightapplyFullStateDiagMatrPower(Qureg qureg, FullStateDiagMatr matrix, qcomp exponent)
void leftapplyFullStateDiagMatrPower(Qureg qureg, FullStateDiagMatr matrix, qcomp exponent)
void leftapplyMultiQubitNot(Qureg qureg, int *targets, int numTargets)
void rightapplyMultiQubitNot(Qureg qureg, int *targets, int numTargets)
void leftapplyPauliX(Qureg qureg, int target)
void rightapplyPauliY(Qureg qureg, int target)
void leftapplyPauliY(Qureg qureg, int target)
void leftapplyPauliZ(Qureg qureg, int target)
void rightapplyPauliX(Qureg qureg, int target)
void rightapplyPauliZ(Qureg qureg, int target)
void leftapplyPauliGadget(Qureg qureg, PauliStr str, qreal angle)
void rightapplyPauliGadget(Qureg qureg, PauliStr str, qreal angle)
void rightapplyPauliStr(Qureg qureg, PauliStr str)
void leftapplyPauliStr(Qureg qureg, PauliStr str)
void rightapplyPauliStrSum(Qureg qureg, PauliStrSum sum, Qureg workspace)
void leftapplyPauliStrSum(Qureg qureg, PauliStrSum sum, Qureg workspace)
void rightapplyPhaseGadget(Qureg qureg, int *targets, int numTargets, qreal angle)
void leftapplyPhaseGadget(Qureg qureg, int *targets, int numTargets, qreal angle)
void leftapplyMultiQubitProjector(Qureg qureg, int *qubits, int *outcomes, int numQubits)
void rightapplyQubitProjector(Qureg qureg, int qubit, int outcome)
void rightapplyMultiQubitProjector(Qureg qureg, int *qubits, int *outcomes, int numQubits)
void leftapplyQubitProjector(Qureg qureg, int qubit, int outcome)
void leftapplySwap(Qureg qureg, int qubit1, int qubit2)
void rightapplySwap(Qureg qureg, int qubit1, int qubit2)
PauliStr getPauliStr(const char *paulis, int *indices, int numPaulis)
Definition paulis.cpp:76
Definition qureg.h:49