The Quantum Exact Simulation Toolkit v4.0.0
Loading...
Searching...
No Matches
operations.cpp
1/** @file
2 * Unit tests of the operations module.
3 *
4 * @author Tyson Jones
5 *
6 * @defgroup unitops Operations
7 * @ingroup unittests
8 */
9
10#include "quest/include/quest.h"
11
12#include <catch2/catch_test_macros.hpp>
13#include <catch2/catch_template_test_macros.hpp>
14#include <catch2/matchers/catch_matchers_string.hpp>
15#include <catch2/generators/catch_generators_range.hpp>
16
17#include "tests/utils/cache.hpp"
18#include "tests/utils/qvector.hpp"
19#include "tests/utils/qmatrix.hpp"
20#include "tests/utils/compare.hpp"
21#include "tests/utils/convert.hpp"
22#include "tests/utils/evolve.hpp"
23#include "tests/utils/linalg.hpp"
24#include "tests/utils/lists.hpp"
25#include "tests/utils/measure.hpp"
26#include "tests/utils/macros.hpp"
27#include "tests/utils/random.hpp"
28
29#include <tuple>
30
31using std::tuple;
32
33
34
35/*
36 * UTILITIES
37 */
38
39#define TEST_CATEGORY \
40 LABEL_UNIT_TAG "[operations]"
41
42
43/*
44 * reference operator matrices used by testing
45 */
46
47namespace FixedMatrices {
48
49 qmatrix H = {
50 {1/std::sqrt(2), 1/std::sqrt(2)},
51 {1/std::sqrt(2), -1/std::sqrt(2)}};
52
53 qmatrix X = getPauliMatrix(1);
54 qmatrix Y = getPauliMatrix(2);
55 qmatrix Z = getPauliMatrix(3);
56
57 qreal PI = 3.14159265358979323846;
58 qmatrix T = {
59 {1, 0},
60 {0, std::exp(1_i * PI/4)}};
61
62 qmatrix S = {
63 {1, 0},
64 {0, 1_i}};
65
66 qmatrix SWAP = {
67 {1, 0, 0, 0},
68 {0, 0, 1, 0},
69 {0, 1, 0, 0},
70 {0, 0, 0, 1}};
71
72 qmatrix sqrtSWAP = {
73 {1, 0, 0, 0},
74 {0, (1+1_i)/2, (1-1_i)/2, 0},
75 {0, (1-1_i)/2, (1+1_i)/2, 0},
76 {0, 0, 0, 1}};
77}
78
79namespace ParameterisedMatrices {
80
81 auto Rx = [](qreal p) { return getExponentialOfPauliMatrix(p, FixedMatrices::X); };
82 auto Ry = [](qreal p) { return getExponentialOfPauliMatrix(p, FixedMatrices::Y); };
83 auto Rz = [](qreal p) { return getExponentialOfPauliMatrix(p, FixedMatrices::Z); };
84
85 auto PS = [](qreal p) { return qmatrix{{1, 0}, {0, std::exp(p*1_i)}}; };
86 auto PS2 = [](qreal p) { return getControlledMatrix(PS(p), 1); };
87}
88
89namespace VariableSizeMatrices {
90
91 auto X = [](int n) { return getKroneckerProduct(FixedMatrices::X, n); };
92 auto PF = [](int n) { return getControlledMatrix(FixedMatrices::Z, n - 1); };
93}
94
95namespace VariableSizeParameterisedMatrices {
96
97 auto Z = [](qreal p, int n) {
98 qmatrix m = getKroneckerProduct(FixedMatrices::Z, n);
99 return getExponentialOfPauliMatrix(p, m);
100 };
101
102 auto PS = [](qreal p, int n) {
103 qmatrix m = ParameterisedMatrices::PS(p);
104 return getControlledMatrix(m, n - 1);
105 };
106}
107
108
109/*
110 * execute 'function' upon each kind of cached qureg
111 * (e.g. distributed, GPU-accelerated, etc) and a
112 * reference state (T1 = qvector or qmatrix), when
113 * both are initialised in the debug state, and
114 * thereafter confirm they approximately agree. Each
115 * qureg deployment is featured in a separate test
116 * section, so are accounted distinctly.
117 */
118
119void TEST_ON_CACHED_QUREGS(quregCache quregs, auto& reference, auto& function) {
120
121 for (auto& [label, qureg]: quregs) {
122
123 DYNAMIC_SECTION( label ) {
124
125 // no need to validate whether qureg successfully
126 // enters the debug state here, because the below
127 // serial setToDebugState() is guaranteed to succeed
128 initDebugState(qureg);
129 setToDebugState(reference);
130
131 function(qureg, reference);
132 REQUIRE_AGREE( qureg, reference );
133 }
134 }
135}
136
137void TEST_ON_CACHED_QUREG_AND_MATRIX(quregCache quregs, matrixCache matrices, auto apiFunc, auto refState, auto refMatr, auto refFunc) {
138
139 for (auto& [labelA, qureg]: quregs) {
140 for (auto& [labelB, matrix]: matrices) {
141
142 // skip illegal (distributed matrix, local qureg) combo
143 if (matrix.isDistributed && ! qureg.isDistributed)
144 continue;
145
146 DYNAMIC_SECTION( labelA + LABEL_DELIMITER + labelB ) {
147
148 // set qureg and reference to debug
149 initDebugState(qureg);
150 setToDebugState(refState);
151
152 // set API matrix to pre-initialised ref matrix
153 setFullStateDiagMatr(matrix, 0, getDiagonals(refMatr));
154
155 // API and reference functions should produce agreeing states
156 apiFunc(qureg, matrix);
157 refFunc(refState, refMatr);
158 REQUIRE_AGREE( qureg, refState );
159 }
160 }
161 }
162}
163
164
165/*
166 * simply avoids boilerplate
167 */
168
169#define PREPARE_TEST( numQubits, statevecQuregs, densmatrQuregs, statevecRef, densmatrRef ) \
170 int numQubits = getNumCachedQubits(); \
171 auto statevecQuregs = getCachedStatevecs(); \
172 auto densmatrQuregs = getCachedDensmatrs(); \
173 qvector statevecRef = getZeroVector(getPow2(numQubits)); \
174 qmatrix densmatrRef = getZeroMatrix(getPow2(numQubits));
175
176
177/*
178 * Template flags for specifying what kind of additional
179 * arguments (in addition to ctrls/states/targs below) are
180 * accepted by an API operation, when passing said operation
181 * to automated testing facilities. This is NOT consulted
182 * when generically invoking the API operation (we instead
183 * used variadic templates above for that), but IS used
184 * by the testing code to decide how to prepare inputs.
185 *
186 * For example:
187 * - applyHadamard: none
188 * - applyRotateX: scalar
189 * - applyRotateAroundAxis: axisrots
190 * - applyDiagMatr1: diagmatr
191 * - applyDiagMatrPower: diagpower
192 * - applyCompMatr: compmatr
193 * - applyPauliStr: paulistr
194 * - applyPauliGadgt: pauligad
195*/
196
197enum ArgsFlag { none, scalar, axisrots, diagmatr, diagpower, compmatr, paulistr, pauligad };
198
199
200/*
201 * Template flags for specifying how many control and
202 * target qubits are accepted by an API operation.
203 * Value 'anystates' is reserved for control qubits,
204 * indicating that ctrls must accompany ctrl-states in
205 * the API signature.
206 *
207 * For example:
208 * - applyHadamard: Ctrls=zero, Targs=one
209 * - applyControlledSwap: Ctrls=one, Targs=two
210 * - applyMultiControlledCompMatr: Ctrls=any, Targs=any
211 * - applyMultiStateControlledT: Ctrls=anystates, Targs=one
212 */
213
214enum NumQubitsFlag { zero, one, two, any, anystates };
215
216void assertNumQubitsFlagsAreValid(NumQubitsFlag ctrlsFlag, NumQubitsFlag targsFlag) {
217
218 DEMAND(
219 ctrlsFlag == zero ||
220 ctrlsFlag == one ||
221 ctrlsFlag == any ||
222 ctrlsFlag == anystates );
223
224 DEMAND(
225 targsFlag == zero ||
226 targsFlag == one ||
227 targsFlag == two ||
228 targsFlag == any);
229}
230
231void assertNumQubitsFlagsValid(
232 NumQubitsFlag ctrlsFlag, NumQubitsFlag targsFlag,
233 vector<int> ctrls, vector<int> states, vector<int> targs
234) {
235 assertNumQubitsFlagsAreValid(ctrlsFlag, targsFlag);
236
237 // we explicitly permit targsFlag=zero while
238 // targs.size() != 0, which occurs when targs
239 // are supplied to the API through an alternate
240 // argument (e.g. a PauliStr)
241
242 if (targsFlag == one)
243 DEMAND( targs.size() == 1 );
244
245 if (targsFlag == two)
246 DEMAND( targs.size() == 2 );
247
248 if (ctrlsFlag == zero)
249 DEMAND( ctrls.size() == 0 );
250
251 if (ctrlsFlag == one)
252 DEMAND( ctrls.size() == 1 );
253
254 if (ctrlsFlag == anystates)
255 DEMAND( states.size() == ctrls.size() );
256 else
257 DEMAND( states.size() == 0 );
258}
259
260
261/*
262 * extract the runtime values of the number of control
263 * and target qubtits from their compile-time templates.
264 * When their number is permitted to be 'any', we use
265 * a generator to successively test all possible numbers.
266 * As such, this function is named similar to a Catch2
267 * macro so the caller recognises it is a generator.
268 */
269
270template <NumQubitsFlag Ctrls, NumQubitsFlag Targs, ArgsFlag Args>
271int GENERATE_NUM_TARGS(int numQuregQubits) {
272
273 assertNumQubitsFlagsAreValid(zero, Targs);
274 DEMAND( Targs != one || numQuregQubits >= 1 );
275 DEMAND( Targs != two || numQuregQubits >= 2 );
276 DEMAND( numQuregQubits > 0 );
277
278 // single choice if #targs is compile-time set
279 if constexpr (Targs == one)
280 return 1;
281 if constexpr (Targs == two)
282 return 2;
283
284 if constexpr (Targs == any) {
285
286 // we can target all qubits...
287 int maxNumTargs = numQuregQubits;
288
289 // unless we're applying CompMatr in distributed
290 // mode. Technically we can support all targets
291 // if the Qureg is a density matrix or not
292 // distributed, but we safely choose the min here
293 if (Args == compmatr)
294 maxNumTargs = numQuregQubits - getLog2(getQuESTEnv().numNodes);
295
296 // we must also ensure there is space left for a forced ctrl
297 if (Ctrls == one && maxNumTargs == numQuregQubits)
298 maxNumTargs -= 1;
299
300 return GENERATE_COPY( range(1, maxNumTargs+1) );
301 }
302}
303
304template <NumQubitsFlag Ctrls>
305int GENERATE_NUM_CTRLS(int numFreeQubits) {
306
307 assertNumQubitsFlagsAreValid(Ctrls, one);
308 DEMAND( Ctrls != one || numFreeQubits >= 1 );
309 DEMAND( numFreeQubits >= 0 );
310
311 if constexpr (Ctrls == zero)
312 return 0;
313
314 if constexpr (Ctrls == one)
315 return 1;
316
317 if constexpr (Ctrls == any || Ctrls == anystates)
318 return GENERATE_COPY( range(0, numFreeQubits+1) );
319}
320
321
322/*
323 * invoke an API operation (e.g. applyHadamard), passing
324 * any elements of (ctrls,states,targs) that it accepts
325 * (as informed by the template values) along with any
326 * additional arguments. We explicitly accept numCtrls
327 * and numTargs, rather than inferring them from ctrls
328 * and targs, such that invalid numbers (like -1) can be
329 * passed by input validation testing.
330 *
331 * This big, ugly bespoke function is necessary, rather
332 * than a simple variadic template, because the QuEST
333 * API accepts fixed numbers of qubits as individual
334 * arguments, rather than as lists/vectors/pointers. Note
335 * that our use of variadic templates (args) means we do
336 * not need to include ArgsFlag as a template parameter.
337 */
338
339template <NumQubitsFlag Ctrls, NumQubitsFlag Targs>
340void invokeApiOperation(
341 auto operation, Qureg qureg,
342 vector<int> ctrls, vector<int> states, int numCtrls,
343 vector<int> targs, int numTargs,
344 auto&... args
345) {
346 assertNumQubitsFlagsValid(Ctrls, Targs, ctrls, states, targs);
347
348 if constexpr (Ctrls == zero) {
349 if constexpr (Targs == zero) operation(qureg, args...);
350 if constexpr (Targs == one) operation(qureg, targs[0], args...);
351 if constexpr (Targs == two) operation(qureg, targs[0], targs[1], args...);
352 if constexpr (Targs == any) operation(qureg, targs.data(), numTargs, args...);
353 }
354 if constexpr (Ctrls == one) {
355 if constexpr (Targs == zero) operation(qureg, ctrls[0], args...);
356 if constexpr (Targs == one) operation(qureg, ctrls[0], targs[0], args...);
357 if constexpr (Targs == two) operation(qureg, ctrls[0], targs[0], targs[1], args...);
358 if constexpr (Targs == any) operation(qureg, ctrls[0], targs.data(), numTargs, args...);
359 }
360 if constexpr (Ctrls == any) {
361 if constexpr (Targs == zero) operation(qureg, ctrls.data(), numCtrls, args...);
362 if constexpr (Targs == one) operation(qureg, ctrls.data(), numCtrls, targs[0], args...);
363 if constexpr (Targs == two) operation(qureg, ctrls.data(), numCtrls, targs[0], targs[1], args...);
364 if constexpr (Targs == any) operation(qureg, ctrls.data(), numCtrls, targs.data(), numTargs, args...);
365 }
366 if constexpr (Ctrls == anystates) {
367 if constexpr (Targs == zero) operation(qureg, ctrls.data(), states.data(), numCtrls, args...);
368 if constexpr (Targs == one) operation(qureg, ctrls.data(), states.data(), numCtrls, targs[0], args...);
369 if constexpr (Targs == two) operation(qureg, ctrls.data(), states.data(), numCtrls, targs[0], targs[1], args...);
370 if constexpr (Targs == any) operation(qureg, ctrls.data(), states.data(), numCtrls, targs.data(), numTargs, args...);
371 }
372}
373
374// overload to avoid passing numCtrls and numTargs
375
376template <NumQubitsFlag Ctrls, NumQubitsFlag Targs>
377void invokeApiOperation(auto operation, Qureg qureg, vector<int> ctrls, vector<int> states, vector<int> targs, auto&... args) {
378 invokeApiOperation<Ctrls,Targs>(operation, qureg, ctrls, states, ctrls.size(), targs, targs.size(), args...);
379}
380
381
382/*
383 * prepare an API matrix (e.g. CompMatr1), as per
384 * the given template parameters, with random
385 * elements. This is used for testing API functions
386 * which accept matrices.
387 */
388
389template <NumQubitsFlag Targs, ArgsFlag Args>
390auto getRandomApiMatrix(int numTargs) {
391
392 DEMAND(
393 Args == diagmatr ||
394 Args == diagpower ||
395 Args == compmatr );
396
397 qmatrix qm = (Args == compmatr)?
398 getRandomUnitary(numTargs) :
399 getRandomDiagonalUnitary(numTargs);
400
401 if constexpr (Args == compmatr && Targs == one)
402 return getCompMatr1(qm);
403
404 if constexpr (Args == compmatr && Targs == two)
405 return getCompMatr2(qm);
406
407 if constexpr (Args == compmatr && Targs == any) {
408 CompMatr cm = createCompMatr(numTargs); // must be freed
409 setCompMatr(cm, qm);
410 return cm;
411 }
412
413 qvector dv = getDiagonals(qm);
414 constexpr bool diag = (Args == diagmatr || Args == diagpower);
415
416 if constexpr (diag && Targs == one)
417 return getDiagMatr1(dv);
418
419 if constexpr (diag && Targs == two)
420 return getDiagMatr2(dv);
421
422 if constexpr (diag && Targs == any) {
423 DiagMatr dm = createDiagMatr(numTargs); // must be freed
424 setDiagMatr(dm, dv);
425 return dm;
426 }
427}
428
429
430/*
431 * chooses random values for the remaining arguments
432 * (after controls/states/targets) to API operations,
433 * with types informed by the template parameters.
434 * For example:
435 * - applyRotateX() accepts a scalar
436 * - applyCompMatr() accepts a CompMatr
437 * - applyDiagMatrPower() accepts a DiagMatr and qcomp
438 */
439
440template <NumQubitsFlag Targs, ArgsFlag Args>
441auto getRandomRemainingArgs(vector<int> targs) {
442
443 if constexpr (Args == none)
444 return tuple{ };
445
446 if constexpr (Args == scalar) {
447 qreal angle = getRandomPhase();
448 return tuple{ angle };
449 }
450
451 if constexpr (Args == axisrots) {
452 qreal angle = getRandomPhase();
453 qreal x = getRandomReal(-1, 1);
454 qreal y = getRandomReal(-1, 1);
455 qreal z = getRandomReal(-1, 1);
456 return tuple{ angle, x, y, z };
457 }
458
459 if constexpr (Args == compmatr || Args == diagmatr) {
460 auto matrix = getRandomApiMatrix<Targs,Args>(targs.size()); // allocates heap mem
461 return tuple{ matrix };
462 }
463
464 if constexpr (Args == diagpower) {
465 DiagMatr matrix = getRandomApiMatrix<Targs,Args>(targs.size()); // allocates heap mem
466 qcomp exponent = qcomp(getRandomReal(-3, 3), 0); // real for unitarity
467 return tuple{ matrix, exponent };
468 }
469
470 if constexpr (Args == paulistr) {
471 PauliStr str = getRandomPauliStr(targs);
472 return tuple{ str };
473 }
474
475 if constexpr (Args == pauligad) {
476 PauliStr str = getRandomPauliStr(targs);
477 qreal angle = getRandomPhase();
478 return tuple{ str, angle };
479 }
480}
481
482
483template <NumQubitsFlag Targs, ArgsFlag Args>
484void freeRemainingArgs(auto args) {
485
486 if constexpr (Targs == any && Args == compmatr)
487 destroyCompMatr(std::get<0>(args));
488
489 if constexpr (Targs == any && Args == diagmatr)
490 destroyDiagMatr(std::get<0>(args));
491
492 if constexpr (Targs == any && Args == diagpower)
493 destroyDiagMatr(std::get<0>(args));
494}
495
496
497/*
498 * unpack the given reference operator matrix (a qmatrix)
499 * for an API operation, which is passed to testOperation(),
500 * and which will be effected upon the reference state (a
501 * qvector or qmatrix). The type/form of matrixRefGen depends
502 * on the type of API operation, indicated by template parameter.
503 */
504
505template <NumQubitsFlag Targs, ArgsFlag Args>
506qmatrix getReferenceMatrix(auto matrixRefGen, vector<int> targs, auto additionalArgs) {
507
508 if constexpr (Args == none && Targs != any)
509 return matrixRefGen;
510
511 if constexpr (Args == none && Targs == any)
512 return matrixRefGen(targs.size());
513
514 if constexpr (Args == scalar && Targs != any) {
515 qreal angle = std::get<0>(additionalArgs);
516 return matrixRefGen(angle);
517 }
518
519 if constexpr (Args == scalar && Targs == any) {
520 qreal angle = std::get<0>(additionalArgs);
521 return matrixRefGen(angle, targs.size());
522 }
523
524 if constexpr (Args == axisrots) {
525 qreal angle = std::get<0>(additionalArgs);
526 qreal x = std::get<1>(additionalArgs);
527 qreal y = std::get<2>(additionalArgs);
528 qreal z = std::get<3>(additionalArgs);
529 return getExponentialOfNormalisedPauliVector(angle, x, y, z);
530 }
531
532 if constexpr (Args == compmatr || Args == diagmatr) {
533 auto apiMatrix = std::get<0>(additionalArgs);
534 return getMatrix(apiMatrix);
535 }
536
537 if constexpr (Args == diagpower) {
538 auto apiMatrix = std::get<0>(additionalArgs);
539 qmatrix diag = getMatrix(apiMatrix);
540 qcomp power = std::get<1>(additionalArgs);
541 return getPowerOfDiagonalMatrix(diag, power);
542 }
543
544 if constexpr (Args == paulistr) {
545 PauliStr str = std::get<0>(additionalArgs);
546 return getMatrix(str, targs);
547 }
548
549 if constexpr (Args == pauligad) {
550 PauliStr str = std::get<0>(additionalArgs);
551 qreal angle = std::get<1>(additionalArgs);
552 qmatrix matr = getMatrix(str, targs);
553 return getExponentialOfPauliMatrix(angle, matr);
554 }
555}
556
557
558/*
559 * display all/only relevant inputs given to an
560 * API operation when its subsequent test fails.
561 * This is like a customisation of CAPTURE, although
562 * we must use UNSCOPED_INFO (and ergo re-implement
563 * some printers) because our branching makes scopes
564 * which end CAPTURE's lifetime.
565 */
566
567
568// @todo surely this should live somewhere else,
569// and/or re-use printer_ functions as much as possible
570
571std::string toString(vector<int> list) {
572
573 std::string out = "{ ";
574 for (int& e : list)
575 out += std::to_string(e) + " ";
576 out += "}";
577 return out;
578}
579
580extern int paulis_getPauliAt(PauliStr str, int ind);
581
582std::string toString(PauliStr str, vector<int> targs) {
583
584 std::string labels = "IXYZ";
585
586 // ugly but adequate - like me (call me)
587 std::string out = "";
588 for (int i=targs.size()-1; i>=0; i--)
589 out += labels[paulis_getPauliAt(str, targs[i])];
590
591 return out;
592}
593
594template <NumQubitsFlag Ctrls, NumQubitsFlag Targs, ArgsFlag Args>
595void CAPTURE_RELEVANT( vector<int> ctrls, vector<int> states, vector<int> targs, auto& args ) {
596
597 // display ctrls
598 if (Ctrls == one)
599 UNSCOPED_INFO( "control := " << ctrls[0] );
600 if (Ctrls == any || Ctrls == anystates )
601 UNSCOPED_INFO( "controls := " << toString(ctrls) );
602
603 // display states
604 if (Ctrls == anystates)
605 UNSCOPED_INFO( "states := " << toString(states) );
606
607 // display targs
608 if (Targs == one)
609 UNSCOPED_INFO( "target := " << targs[0] );
610
611 if (Targs == two) {
612 UNSCOPED_INFO( "target A := " << targs[0] );
613 UNSCOPED_INFO( "target B := " << targs[1] );
614 }
615 if (Targs == any)
616 UNSCOPED_INFO( "targets := " << toString(targs) );
617
618 // display rotation angle
619 if constexpr (Args == scalar)
620 UNSCOPED_INFO( "angle := " << std::get<0>(args) );
621
622 // display rotation angle and axis
623 if constexpr (Args == axisrots) {
624 UNSCOPED_INFO( "angle := " << std::get<0>(args) );
625 UNSCOPED_INFO( "x := " << std::get<1>(args) );
626 UNSCOPED_INFO( "y := " << std::get<2>(args) );
627 UNSCOPED_INFO( "z := " << std::get<3>(args) );
628 }
629
630 // note but don't display API matrices
631 if constexpr (Args == compmatr || Args == diagmatr || Args == diagpower)
632 UNSCOPED_INFO( "matrix := (omitted)" );
633
634 // display exponent of diagonal matrix
635 if constexpr (Args == diagpower) {
636 qcomp p = std::get<1>(args);
637 UNSCOPED_INFO( "exponent := " << std::real(p) << " + (" << std::imag(p) << ")i" );
638 }
639
640 // display PauliStr
641 if constexpr (Args == paulistr || Args == pauligad)
642 UNSCOPED_INFO( "paulis := " << toString(std::get<0>(args), targs) );
643
644 // display PauliStr angle
645 if constexpr (Args == pauligad)
646 UNSCOPED_INFO( "angle := " << std::get<1>(args) );
647}
648
649
650/*
651 * perform a unit test on an API operation. The
652 * template parameters are compile-time clues
653 * about what inputs to prepare and pass to the
654 * operation, and how its reference matrix (arg
655 * matrixRefGen) is formatted.
656 */
657
658template <NumQubitsFlag Ctrls, NumQubitsFlag Targs, ArgsFlag Args>
659void testOperation(auto operation, auto matrixRefGen, bool multiplyOnly) {
660
661 assertNumQubitsFlagsAreValid(Ctrls, Targs);
662
663 PREPARE_TEST( numQubits, statevecQuregs, densmatrQuregs, statevecRef, densmatrRef );
664
665 SECTION( LABEL_CORRECTNESS ) {
666
667 // try all possible number of ctrls and targs
668 int numTargs = GENERATE_NUM_TARGS<Ctrls,Targs,Args>(numQubits);
669 int numCtrls = GENERATE_NUM_CTRLS<Ctrls>(numQubits - numTargs);
670
671 // either try all possible ctrls and targs, or randomise them
672 auto ctrlsAndTargs = GENERATE_CTRLS_AND_TARGS( numQubits, numCtrls, numTargs );
673 vector<int> ctrls = std::get<0>(ctrlsAndTargs);
674 vector<int> targs = std::get<1>(ctrlsAndTargs);
675
676 // randomise control states (if operation accepts them)
677 vector<int> states = getRandomOutcomes(numCtrls * (Ctrls == anystates));
678
679 // randomise remaining operation parameters
680 auto primaryArgs = tuple{ ctrls, states, targs };
681 auto furtherArgs = getRandomRemainingArgs<Targs,Args>(targs); // may allocate heap memory
682
683 // obtain the reference matrix for this operation
684 qmatrix matrixRef = getReferenceMatrix<Targs,Args>(matrixRefGen, targs, furtherArgs);
685
686 // PauliStr arg replaces target qubit list in API operations
687 constexpr NumQubitsFlag RevTargs = (Args==paulistr||Args==pauligad)? zero : Targs;
688
689 // disabling unitary-check validation for compmatr, since it's hard to
690 // generate numerically-precise random unitaries upon many qubits, or
691 // upon few qubits are single-precision. So we disable completely until
692 // we re-implement 'input validation' checks which force us to fix thresholds
693 (Args == compmatr)?
696
697 // prepare test function which will receive both statevectors and density matrices
698 auto testFunc = [&](Qureg qureg, auto& stateRef) -> void {
699
700 // invoke API operation, passing all args (unpacking variadic)
701 auto apiFunc = [](auto&&... args) { return invokeApiOperation<Ctrls,RevTargs>(args...); };
702 auto allArgs = std::tuple_cat(tuple{operation, qureg}, primaryArgs, furtherArgs);
703 std::apply(apiFunc, allArgs);
704
705 // update reference state
706 (multiplyOnly)?
707 multiplyReferenceOperator(stateRef, ctrls, states, targs, matrixRef):
708 applyReferenceOperator(stateRef, ctrls, states, targs, matrixRef);
709 };
710
711 // report operation's input parameters if any subsequent test fails
712 CAPTURE_RELEVANT<Ctrls,Targs,Args>( ctrls, states, targs, furtherArgs );
713
714 // test API operation on all available deployment combinations (e.g. OMP, MPI, MPI+GPU, etc)
715 SECTION( LABEL_STATEVEC ) { TEST_ON_CACHED_QUREGS(statevecQuregs, statevecRef, testFunc); }
716 SECTION( LABEL_DENSMATR ) { TEST_ON_CACHED_QUREGS(densmatrQuregs, densmatrRef, testFunc); }
717
718 // free any heap-alloated API matrices and restore epsilon
719 freeRemainingArgs<Targs,Args>(furtherArgs);
721 }
722
723 // @todo
724 // input validation!
725}
726
727template <NumQubitsFlag Ctrls, NumQubitsFlag Targs, ArgsFlag Args>
728void testOperation(auto operation, auto matrixRefGen) {
729
730 // default multiplyOnly=false
731 testOperation<Ctrls,Targs,Args>(operation, matrixRefGen, false);
732}
733
734
735/*
736 * perform unit tests for the four distinctly-controlled
737 * variants of the given API operation (with the specified
738 * function name suffix). 'numtargs' indicates the number
739 * of target qubits accepted by the operation, and 'argtype'
740 * indicates the types of remaining arguments (if any exist).
741 * 'matrixgen' is the matrix representation (of varying
742 * formats) of the operation, against which it will be compared.
743 */
744
745#define TEST_ANY_CTRL_OPERATION( namesuffix, numtargs, argtype, matrixgen ) \
746 TEST_CASE( "apply" #namesuffix, TEST_CATEGORY ) { testOperation<zero, numtargs,argtype>( apply ## namesuffix, matrixgen); } \
747 TEST_CASE( "applyControlled" #namesuffix, TEST_CATEGORY ) { testOperation<one, numtargs,argtype>( applyControlled ## namesuffix, matrixgen); } \
748 TEST_CASE( "applyMultiControlled" #namesuffix, TEST_CATEGORY ) { testOperation<any, numtargs,argtype>( applyMultiControlled ## namesuffix, matrixgen); } \
749 TEST_CASE( "applyMultiStateControlled" #namesuffix, TEST_CATEGORY ) { testOperation<anystates,numtargs,argtype>( applyMultiStateControlled ## namesuffix, matrixgen); }
750
751
752
753/**
754 * TESTS
755 *
756 * @ingroup unitops
757 * @{
758 */
759
760TEST_ANY_CTRL_OPERATION( PauliStr, any, paulistr, nullptr );
761TEST_ANY_CTRL_OPERATION( PauliGadget, any, pauligad, nullptr );
762TEST_ANY_CTRL_OPERATION( CompMatr1, one, compmatr, nullptr );
763TEST_ANY_CTRL_OPERATION( CompMatr2, two, compmatr, nullptr );
764TEST_ANY_CTRL_OPERATION( CompMatr, any, compmatr, nullptr );
765TEST_ANY_CTRL_OPERATION( DiagMatr1, one, diagmatr, nullptr );
766TEST_ANY_CTRL_OPERATION( DiagMatr2, two, diagmatr, nullptr );
767TEST_ANY_CTRL_OPERATION( DiagMatr, any, diagmatr, nullptr );
768TEST_ANY_CTRL_OPERATION( DiagMatrPower, any, diagpower, nullptr );
769
770TEST_ANY_CTRL_OPERATION( Hadamard, one, none, FixedMatrices::H );
771TEST_ANY_CTRL_OPERATION( PauliX, one, none, FixedMatrices::X );
772TEST_ANY_CTRL_OPERATION( PauliY, one, none, FixedMatrices::Y );
773TEST_ANY_CTRL_OPERATION( PauliZ, one, none, FixedMatrices::Z );
774TEST_ANY_CTRL_OPERATION( T, one, none, FixedMatrices::T );
775TEST_ANY_CTRL_OPERATION( S, one, none, FixedMatrices::S );
776TEST_ANY_CTRL_OPERATION( Swap, two, none, FixedMatrices::SWAP );
777TEST_ANY_CTRL_OPERATION( SqrtSwap, two, none, FixedMatrices::sqrtSWAP );
778
779TEST_ANY_CTRL_OPERATION( RotateX, one, scalar, ParameterisedMatrices::Rx );
780TEST_ANY_CTRL_OPERATION( RotateY, one, scalar, ParameterisedMatrices::Ry );
781TEST_ANY_CTRL_OPERATION( RotateZ, one, scalar, ParameterisedMatrices::Rz );
782
783TEST_ANY_CTRL_OPERATION( RotateAroundAxis, one, axisrots, nullptr );
784
785TEST_ANY_CTRL_OPERATION( MultiQubitNot, any, none, VariableSizeMatrices::X );
786
787TEST_ANY_CTRL_OPERATION( PhaseGadget, any, scalar, VariableSizeParameterisedMatrices::Z );
788
789TEST_CASE( "multiplyPauliStr", TEST_CATEGORY ) { testOperation<zero,any,paulistr>(multiplyPauliStr, nullptr, true); }
790TEST_CASE( "multiplyPauliGadget", TEST_CATEGORY ) { testOperation<zero,any,pauligad>(multiplyPauliGadget, nullptr, true); }
791
792TEST_CASE( "multiplyCompMatr1", TEST_CATEGORY ) { testOperation<zero,one,compmatr>(multiplyCompMatr1, nullptr, true); }
793TEST_CASE( "multiplyCompMatr2", TEST_CATEGORY ) { testOperation<zero,two,compmatr>(multiplyCompMatr2, nullptr, true); }
794TEST_CASE( "multiplyCompMatr", TEST_CATEGORY ) { testOperation<zero,any,compmatr>(multiplyCompMatr, nullptr, true); }
795
796TEST_CASE( "multiplyDiagMatr1", TEST_CATEGORY ) { testOperation<zero,one,diagmatr>(multiplyDiagMatr1, nullptr, true); }
797TEST_CASE( "multiplyDiagMatr2", TEST_CATEGORY ) { testOperation<zero,two,diagmatr>(multiplyDiagMatr2, nullptr, true); }
798TEST_CASE( "multiplyDiagMatr", TEST_CATEGORY ) { testOperation<zero,any,diagmatr>(multiplyDiagMatr, nullptr, true); }
799
800TEST_CASE( "multiplyDiagMatrPower", TEST_CATEGORY ) { testOperation<zero,any,diagpower>(multiplyDiagMatrPower, nullptr, true); }
801
802TEST_CASE( "multiplyMultiQubitNot", TEST_CATEGORY ) { testOperation<zero,any,none>(multiplyMultiQubitNot, VariableSizeMatrices::X, true); }
803TEST_CASE( "multiplyPhaseGadget", TEST_CATEGORY ) { testOperation<zero,any,scalar>(multiplyPhaseGadget, VariableSizeParameterisedMatrices::Z, true); }
804
805TEST_CASE( "applyPhaseFlip", TEST_CATEGORY ) { testOperation<zero,one,none>(applyPhaseFlip, VariableSizeMatrices::PF(1)); }
806TEST_CASE( "applyTwoQubitPhaseFlip", TEST_CATEGORY ) { testOperation<zero,two,none>(applyTwoQubitPhaseFlip, VariableSizeMatrices::PF(2)); }
807TEST_CASE( "applyMultiQubitPhaseFlip", TEST_CATEGORY ) { testOperation<zero,any,none>(applyMultiQubitPhaseFlip, VariableSizeMatrices::PF); }
808
809TEST_CASE( "applyPhaseShift", TEST_CATEGORY ) { testOperation<zero,one,scalar>(applyPhaseShift, ParameterisedMatrices::PS); }
810TEST_CASE( "applyTwoQubitPhaseShift", TEST_CATEGORY ) { testOperation<zero,two,scalar>(applyTwoQubitPhaseShift, ParameterisedMatrices::PS2); }
811TEST_CASE( "applyMultiQubitPhaseShift", TEST_CATEGORY ) { testOperation<zero,any,scalar>(applyMultiQubitPhaseShift, VariableSizeParameterisedMatrices::PS); }
812
813
814TEST_CASE( "applyQuantumFourierTransform", TEST_CATEGORY ) {
815
816 PREPARE_TEST( numQubits, statevecQuregs, densmatrQuregs, statevecRef, densmatrRef );
817
818 SECTION( LABEL_CORRECTNESS ) {
819
820 int numTargs = GENERATE_COPY( range(1,numQubits+1) );
821 auto targs = GENERATE_TARGS( numQubits, numTargs );
822
823 CAPTURE( targs );
824
825 SECTION( LABEL_STATEVEC ) {
826
827 auto testFunc = [&](Qureg qureg, qvector& ref) {
828 applyQuantumFourierTransform(qureg, targs.data(), targs.size());
829 ref = getDisceteFourierTransform(ref, targs);
830 };
831
832 TEST_ON_CACHED_QUREGS(statevecQuregs, statevecRef, testFunc);
833 }
834
835 SECTION( LABEL_DENSMATR ) {
836
837 // prepare a random mixture
838 auto states = getRandomOrthonormalStateVectors(numQubits, getRandomInt(1,10));
839 auto probs = getRandomProbabilities(states.size());
840
841 auto testFunc = [&](Qureg qureg, qmatrix& ref) {
842
843 // overwrite the Qureg debug state set by caller to above mixture
844 setQuregToReference(qureg, getMixture(states, probs));
845 applyQuantumFourierTransform(qureg, targs.data(), targs.size());
846
847 ref = getZeroMatrix(ref.size());
848 for (size_t i=0; i<states.size(); i++) {
849 qvector vec = getDisceteFourierTransform(states[i], targs);
850 ref += probs[i] * getOuterProduct(vec, vec);
851 }
852 };
853
854 TEST_ON_CACHED_QUREGS(densmatrQuregs, densmatrRef, testFunc);
855 }
856 }
857
858 /// @todo input validation
859}
860
861
862TEST_CASE( "applyFullQuantumFourierTransform", TEST_CATEGORY ) {
863
864 PREPARE_TEST( numQubits, statevecQuregs, densmatrQuregs, statevecRef, densmatrRef );
865
866 SECTION( LABEL_CORRECTNESS ) {
867
868 GENERATE( range(0,10) );
869
870 SECTION( LABEL_STATEVEC ) {
871
872 auto testFunc = [&](Qureg qureg, qvector& ref) {
874 ref = getDisceteFourierTransform(ref);
875 };
876
877 TEST_ON_CACHED_QUREGS(statevecQuregs, statevecRef, testFunc);
878 }
879
880 SECTION( LABEL_DENSMATR ) {
881
882 // prepare a random mixture
883 auto states = getRandomOrthonormalStateVectors(numQubits, getRandomInt(1,10));
884 auto probs = getRandomProbabilities(states.size());
885
886 auto testFunc = [&](Qureg qureg, qmatrix& ref) {
887
888 // overwrite the Qureg debug state set by caller to above mixture
889 setQuregToReference(qureg, getMixture(states, probs));
891
892 ref = getZeroMatrix(ref.size());
893 for (size_t i=0; i<states.size(); i++) {
894 qvector vec = getDisceteFourierTransform(states[i]);
895 ref += probs[i] * getOuterProduct(vec, vec);
896 }
897 };
898
899 TEST_ON_CACHED_QUREGS(densmatrQuregs, densmatrRef, testFunc);
900 }
901 }
902
903 /// @todo input validation
904}
905
906
907TEST_CASE( "applyQubitProjector", TEST_CATEGORY ) {
908
909 PREPARE_TEST( numQubits, statevecQuregs, densmatrQuregs, statevecRef, densmatrRef );
910
911 SECTION( LABEL_CORRECTNESS ) {
912
913 GENERATE( range(0,10) );
914 int target = GENERATE_COPY( range(0,numQubits) );
915 int outcome = GENERATE( 0, 1 );
916
917 qmatrix projector = getProjector(outcome);
918
919 auto testFunc = [&](Qureg qureg, auto& ref) {
920 applyQubitProjector(qureg, target, outcome);
921 applyReferenceOperator(ref, {target}, projector);
922 };
923
924 CAPTURE( target, outcome );
925 SECTION( LABEL_STATEVEC ) { TEST_ON_CACHED_QUREGS(statevecQuregs, statevecRef, testFunc); }
926 SECTION( LABEL_DENSMATR ) { TEST_ON_CACHED_QUREGS(densmatrQuregs, densmatrRef, testFunc); }
927 }
928
929 /// @todo input validation
930}
931
932
933TEST_CASE( "applyMultiQubitProjector", TEST_CATEGORY ) {
934
935 PREPARE_TEST( numQubits, statevecQuregs, densmatrQuregs, statevecRef, densmatrRef );
936
937 SECTION( LABEL_CORRECTNESS ) {
938
939 int numTargs = GENERATE_COPY( range(1,numQubits+1) );
940 auto targets = GENERATE_TARGS( numQubits, numTargs );
941 auto outcomes = getRandomOutcomes(numTargs);
942
943 qmatrix projector = getProjector(targets, outcomes, numQubits);
944
945 auto testFunc = [&](Qureg qureg, auto& ref) {
946 applyMultiQubitProjector(qureg, targets.data(), outcomes.data(), numTargs);
947 applyReferenceOperator(ref, projector);
948 };
949
950 CAPTURE( targets, outcomes );
951 SECTION( LABEL_STATEVEC ) { TEST_ON_CACHED_QUREGS(statevecQuregs, statevecRef, testFunc); }
952 SECTION( LABEL_DENSMATR ) { TEST_ON_CACHED_QUREGS(densmatrQuregs, densmatrRef, testFunc); }
953 }
954
955 /// @todo input validation
956}
957
958
959TEST_CASE( "applyForcedQubitMeasurement", TEST_CATEGORY ) {
960
961 PREPARE_TEST( numQubits, statevecQuregs, densmatrQuregs, statevecRef, densmatrRef );
962
963 SECTION( LABEL_CORRECTNESS ) {
964
965 GENERATE( range(0,10) );
966 int target = GENERATE_COPY( range(0,numQubits) );
967 int outcome = GENERATE( 0, 1 );
968
969 qmatrix projector = getProjector(outcome);
970
971 auto testFunc = [&](Qureg qureg, auto& ref) {
972
973 // overwrite caller's setting of initDebugState, since
974 // that precludes outcomes=|0><0| due to zero-probability
975 setToRandomState(ref);
976 setQuregToReference(qureg, ref);
977
978 // compare the probabilities...
979 qreal apiProb = applyForcedQubitMeasurement(qureg, target, outcome);
980 qreal refProb = getReferenceProbability(ref, {target}, {outcome});
981 REQUIRE_AGREE( apiProb, refProb );
982
983 // and the post-projection states (caller calls subsequent REQUIRE_AGREE)
984 applyReferenceOperator(ref, {target}, projector);
985 ref /= (qureg.isDensityMatrix)?
986 refProb : std::sqrt(refProb);
987 };
988
989 CAPTURE( target, outcome );
990 SECTION( LABEL_STATEVEC ) { TEST_ON_CACHED_QUREGS(statevecQuregs, statevecRef, testFunc); }
991 SECTION( LABEL_DENSMATR ) { TEST_ON_CACHED_QUREGS(densmatrQuregs, densmatrRef, testFunc); }
992 }
993
994 /// @todo input validation
995}
996
997
998TEST_CASE( "applyForcedMultiQubitMeasurement", TEST_CATEGORY ) {
999
1000 PREPARE_TEST( numQubits, statevecQuregs, densmatrQuregs, statevecRef, densmatrRef );
1001
1002 SECTION( LABEL_CORRECTNESS ) {
1003
1004 int numTargs = GENERATE_COPY( range(1,numQubits+1) );
1005 auto targets = GENERATE_TARGS( numQubits, numTargs );
1006 auto outcomes = getRandomOutcomes(numTargs);
1007
1008 qmatrix projector = getProjector(targets, outcomes, numQubits);
1009
1010 auto testFunc = [&](Qureg qureg, auto& ref) {
1011
1012 // overwrite caller's setting of initDebugState, since
1013 // that precludes outcomes=|0><0| due to zero-probability
1014 setToRandomState(ref);
1015 setQuregToReference(qureg, ref);
1016
1017 // compare the probabilities...
1018 qreal apiProb = applyForcedMultiQubitMeasurement(qureg, targets.data(), outcomes.data(), numTargs);
1019 qreal refProb = getReferenceProbability(ref, targets, outcomes);
1020 REQUIRE_AGREE( apiProb, refProb );
1021
1022 // and the post-measurement states (caller calls subsequent REQUIRE_AGREE)
1023 applyReferenceOperator(ref, projector);
1024 ref /= (qureg.isDensityMatrix)?
1025 refProb : std::sqrt(refProb);
1026 };
1027
1028 CAPTURE( targets, outcomes );
1029 SECTION( LABEL_STATEVEC ) { TEST_ON_CACHED_QUREGS(statevecQuregs, statevecRef, testFunc); }
1030 SECTION( LABEL_DENSMATR ) { TEST_ON_CACHED_QUREGS(densmatrQuregs, densmatrRef, testFunc); }
1031 }
1032
1033 /// @todo input validation
1034}
1035
1036
1037TEST_CASE( "applyMultiQubitMeasurement", TEST_CATEGORY ) {
1038
1039 PREPARE_TEST( numQubits, statevecQuregs, densmatrQuregs, statevecRef, densmatrRef );
1040
1041 SECTION( LABEL_CORRECTNESS ) {
1042
1043 int numTargs = GENERATE_COPY( range(1,numQubits+1) );
1044 auto targets = GENERATE_TARGS( numQubits, numTargs );
1045
1046 auto testFunc = [&](Qureg qureg, auto& ref) {
1047
1048 // overwrite caller's setting of initDebugState, since
1049 // sampling requires the outcome probs are normalised
1050 setToRandomState(ref);
1051 setQuregToReference(qureg, ref);
1052
1053 // the output API state...
1054 qindex apiOut = applyMultiQubitMeasurement(qureg, targets.data(), numTargs);
1055
1056 // informs the projector which determines the post-measurement reference
1057 auto apiOutBits = getBits(apiOut, numTargs);
1058 qmatrix projector = getProjector(targets, apiOutBits, numQubits);
1059 applyReferenceOperator(ref, projector);
1060 qreal refProb = getReferenceProbability(ref, targets, apiOutBits);
1061 ref /= (qureg.isDensityMatrix)?
1062 refProb : std::sqrt(refProb);
1063 };
1064
1065 CAPTURE( targets );
1066 SECTION( LABEL_STATEVEC ) { TEST_ON_CACHED_QUREGS(statevecQuregs, statevecRef, testFunc); }
1067 SECTION( LABEL_DENSMATR ) { TEST_ON_CACHED_QUREGS(densmatrQuregs, densmatrRef, testFunc); }
1068 }
1069
1070 /// @todo input validation
1071}
1072
1073
1074TEST_CASE( "applyMultiQubitMeasurementAndGetProb", TEST_CATEGORY ) {
1075
1076 PREPARE_TEST( numQubits, statevecQuregs, densmatrQuregs, statevecRef, densmatrRef );
1077
1078 SECTION( LABEL_CORRECTNESS ) {
1079
1080 int numTargs = GENERATE_COPY( range(1,numQubits+1) );
1081 auto targets = GENERATE_TARGS( numQubits, numTargs );
1082
1083 auto testFunc = [&](Qureg qureg, auto& ref) {
1084
1085 // overwrite caller's setting of initDebugState, since
1086 // sampling requires the outcome probs are normalised
1087 setToRandomState(ref);
1088 setQuregToReference(qureg, ref);
1089
1090 // compare the measurement probability...
1091 qreal apiProb = -1;
1092 qindex apiOut = applyMultiQubitMeasurementAndGetProb(qureg, targets.data(), numTargs, &apiProb);
1093 auto apiOutBits = getBits(apiOut, numTargs);
1094 qreal refProb = getReferenceProbability(ref, targets, apiOutBits);
1095 REQUIRE_AGREE( apiProb, refProb );
1096
1097 // and the post-measurement states (caller calls subsequent REQUIRE_AGREE)
1098 qmatrix projector = getProjector(targets, apiOutBits, numQubits);
1099 applyReferenceOperator(ref, projector);
1100 ref /= (qureg.isDensityMatrix)?
1101 refProb : std::sqrt(refProb);
1102 };
1103
1104 CAPTURE( targets );
1105 SECTION( LABEL_STATEVEC ) { TEST_ON_CACHED_QUREGS(statevecQuregs, statevecRef, testFunc); }
1106 SECTION( LABEL_DENSMATR ) { TEST_ON_CACHED_QUREGS(densmatrQuregs, densmatrRef, testFunc); }
1107 }
1108
1109 /// @todo input validation
1110}
1111
1112
1113TEST_CASE( "applyQubitMeasurement", TEST_CATEGORY ) {
1114
1115 PREPARE_TEST( numQubits, statevecQuregs, densmatrQuregs, statevecRef, densmatrRef );
1116
1117 SECTION( LABEL_CORRECTNESS ) {
1118
1119 GENERATE( range(0,10) );
1120 int target = GENERATE_COPY( range(0,numQubits) );
1121
1122 auto testFunc = [&](Qureg qureg, auto& ref) {
1123
1124 // overwrite caller's setting of initDebugState, since
1125 // sampling requires the outcome probs are normalised
1126 setToRandomState(ref);
1127 setQuregToReference(qureg, ref);
1128
1129 // the output API state...
1130 int apiOut = applyQubitMeasurement(qureg, target);
1131
1132 // informs the projector which determines the post-measurement reference
1133 qmatrix projector = getProjector(apiOut);
1134 applyReferenceOperator(ref, {target}, projector);
1135 qreal refProb = getReferenceProbability(ref, {target}, {apiOut});
1136 ref /= (qureg.isDensityMatrix)?
1137 refProb : std::sqrt(refProb);
1138 };
1139
1140 CAPTURE( target );
1141 SECTION( LABEL_STATEVEC ) { TEST_ON_CACHED_QUREGS(statevecQuregs, statevecRef, testFunc); }
1142 SECTION( LABEL_DENSMATR ) { TEST_ON_CACHED_QUREGS(densmatrQuregs, densmatrRef, testFunc); }
1143 }
1144
1145 /// @todo input validation
1146}
1147
1148
1149TEST_CASE( "applyQubitMeasurementAndGetProb", TEST_CATEGORY ) {
1150
1151 PREPARE_TEST( numQubits, statevecQuregs, densmatrQuregs, statevecRef, densmatrRef );
1152
1153 SECTION( LABEL_CORRECTNESS ) {
1154
1155 GENERATE( range(0,10) );
1156 int target = GENERATE_COPY( range(0,numQubits) );
1157
1158 auto testFunc = [&](Qureg qureg, auto& ref) {
1159
1160 // overwrite caller's setting of initDebugState, since
1161 // sampling requires the outcome probs are normalised
1162 setToRandomState(ref);
1163 setQuregToReference(qureg, ref);
1164
1165 // compare the measurement probability...
1166 qreal apiProb = -1;
1167 int apiOut = applyQubitMeasurementAndGetProb(qureg, target, &apiProb);
1168 qreal refProb = getReferenceProbability(ref, {target}, {apiOut});
1169 REQUIRE_AGREE( apiProb, refProb );
1170
1171 // and the post-measurement states (caller calls subsequent REQUIRE_AGREE)
1172 qmatrix projector = getProjector(apiOut);
1173 applyReferenceOperator(ref, {target}, projector);
1174 ref /= (qureg.isDensityMatrix)?
1175 refProb : std::sqrt(refProb);
1176 };
1177
1178 CAPTURE( target );
1179 SECTION( LABEL_STATEVEC ) { TEST_ON_CACHED_QUREGS(statevecQuregs, statevecRef, testFunc); }
1180 SECTION( LABEL_DENSMATR ) { TEST_ON_CACHED_QUREGS(densmatrQuregs, densmatrRef, testFunc); }
1181 }
1182
1183 /// @todo input validation
1184}
1185
1186
1187TEST_CASE( "multiplyFullStateDiagMatr", TEST_CATEGORY LABEL_MIXED_DEPLOY_TAG ) {
1188
1189 PREPARE_TEST( numQubits, cachedSV, cachedDM, refSV, refDM );
1190
1191 auto cachedMatrs = getCachedFullStateDiagMatrs();
1192
1193 SECTION( LABEL_CORRECTNESS ) {
1194
1195 qmatrix refMatr = getRandomDiagonalMatrix(getPow2(numQubits));
1196 auto apiFunc = multiplyFullStateDiagMatr;
1197
1198 GENERATE( range(0, TEST_NUM_MIXED_DEPLOYMENT_REPETITIONS) );
1199
1200 SECTION( LABEL_STATEVEC ) {
1201
1202 auto refFunc = [&] (qvector& state, qmatrix matr) { multiplyReferenceOperator(state, matr); };
1203
1204 TEST_ON_CACHED_QUREG_AND_MATRIX( cachedSV, cachedMatrs, apiFunc, refSV, refMatr, refFunc);
1205 }
1206
1207 SECTION( LABEL_DENSMATR ) {
1208
1209 auto refFunc = [&] (qmatrix& state, qmatrix matr) { multiplyReferenceOperator(state, matr); };
1210
1211 TEST_ON_CACHED_QUREG_AND_MATRIX( cachedDM, cachedMatrs, apiFunc, refDM, refMatr, refFunc);
1212 }
1213 }
1214
1215 /// @todo input validation
1216}
1217
1218
1219TEST_CASE( "multiplyFullStateDiagMatrPower", TEST_CATEGORY LABEL_MIXED_DEPLOY_TAG ) {
1220
1221 PREPARE_TEST( numQubits, cachedSV, cachedDM, refSV, refDM );
1222
1223 auto cachedMatrs = getCachedFullStateDiagMatrs();
1224
1225 SECTION( LABEL_CORRECTNESS ) {
1226
1227 qmatrix refMatr = getRandomDiagonalMatrix(getPow2(numQubits));
1228 qcomp exponent = getRandomComplex();
1229
1230 auto apiFunc = [&](Qureg qureg, FullStateDiagMatr matr) {
1231 return multiplyFullStateDiagMatrPower(qureg, matr, exponent);
1232 };
1233
1234 CAPTURE( exponent );
1235
1236 GENERATE( range(0, TEST_NUM_MIXED_DEPLOYMENT_REPETITIONS) );
1237
1238 SECTION( LABEL_STATEVEC ) {
1239
1240 auto refFunc = [&] (qvector& state, qmatrix matr) {
1241 matr = getPowerOfDiagonalMatrix(matr, exponent);
1242 multiplyReferenceOperator(state, matr);
1243 };
1244
1245 TEST_ON_CACHED_QUREG_AND_MATRIX( cachedSV, cachedMatrs, apiFunc, refSV, refMatr, refFunc);
1246 }
1247
1248 SECTION( LABEL_DENSMATR ) {
1249
1250 auto refFunc = [&] (qmatrix& state, qmatrix matr) {
1251 matr = getPowerOfDiagonalMatrix(matr, exponent);
1252 multiplyReferenceOperator(state, matr);
1253 };
1254
1255 TEST_ON_CACHED_QUREG_AND_MATRIX( cachedDM, cachedMatrs, apiFunc, refDM, refMatr, refFunc);
1256 }
1257 }
1258
1259 /// @todo input validation
1260}
1261
1262
1263TEST_CASE( "applyFullStateDiagMatr", TEST_CATEGORY LABEL_MIXED_DEPLOY_TAG ) {
1264
1265 PREPARE_TEST( numQubits, cachedSV, cachedDM, refSV, refDM );
1266
1267 auto cachedMatrs = getCachedFullStateDiagMatrs();
1268
1269 SECTION( LABEL_CORRECTNESS ) {
1270
1271 qmatrix refMatr = getRandomDiagonalUnitary(numQubits);
1272 auto apiFunc = applyFullStateDiagMatr;
1273
1274 GENERATE( range(0, TEST_NUM_MIXED_DEPLOYMENT_REPETITIONS) );
1275
1276 SECTION( LABEL_STATEVEC ) {
1277
1278 auto refFunc = [&] (qvector& state, qmatrix matr) { applyReferenceOperator(state, matr); };
1279
1280 TEST_ON_CACHED_QUREG_AND_MATRIX( cachedSV, cachedMatrs, apiFunc, refSV, refMatr, refFunc);
1281 }
1282
1283 SECTION( LABEL_DENSMATR ) {
1284
1285 auto refFunc = [&] (qmatrix& state, qmatrix matr) { applyReferenceOperator(state, matr); };
1286
1287 TEST_ON_CACHED_QUREG_AND_MATRIX( cachedDM, cachedMatrs, apiFunc, refDM, refMatr, refFunc);
1288 }
1289 }
1290
1291 /// @todo input validation
1292}
1293
1294
1295TEST_CASE( "applyFullStateDiagMatrPower", TEST_CATEGORY LABEL_MIXED_DEPLOY_TAG ) {
1296
1297 PREPARE_TEST( numQubits, cachedSV, cachedDM, refSV, refDM );
1298
1299 auto cachedMatrs = getCachedFullStateDiagMatrs();
1300
1301 SECTION( LABEL_CORRECTNESS ) {
1302
1303 qmatrix refMatr = getRandomDiagonalUnitary(numQubits);
1304
1305 // supplying a complex exponent requires disabling
1306 // numerical validation to relax unitarity
1307 bool testRealExp = GENERATE( true, false );
1308 qcomp exponent = (testRealExp)?
1309 qcomp(getRandomReal(-2, 2), 0):
1311
1312 auto apiFunc = [&](Qureg qureg, FullStateDiagMatr matr) {
1313 return applyFullStateDiagMatrPower(qureg, matr, exponent);
1314 };
1315
1316 CAPTURE( exponent );
1317
1318 GENERATE( range(0, TEST_NUM_MIXED_DEPLOYMENT_REPETITIONS) );
1319
1320 if (!testRealExp)
1322
1323 SECTION( LABEL_STATEVEC ) {
1324
1325 auto refFunc = [&] (qvector& state, qmatrix matr) {
1326 matr = getPowerOfDiagonalMatrix(matr, exponent);
1327 applyReferenceOperator(state, matr);
1328 };
1329
1330 TEST_ON_CACHED_QUREG_AND_MATRIX( cachedSV, cachedMatrs, apiFunc, refSV, refMatr, refFunc);
1331 }
1332
1333 SECTION( LABEL_DENSMATR ) {
1334
1335 auto refFunc = [&] (qmatrix& state, qmatrix matr) {
1336 matr = getPowerOfDiagonalMatrix(matr, exponent);
1337 applyReferenceOperator(state, matr);
1338 };
1339
1340 TEST_ON_CACHED_QUREG_AND_MATRIX( cachedDM, cachedMatrs, apiFunc, refDM, refMatr, refFunc);
1341 }
1342
1344 }
1345
1346 /// @todo input validation
1347}
1348
1349
1350TEST_CASE( "multiplyPauliStrSum", TEST_CATEGORY LABEL_MIXED_DEPLOY_TAG ) {
1351
1352 PREPARE_TEST( numQubits, statevecQuregs, densmatrQuregs, statevecRef, densmatrRef );
1353
1354 SECTION( LABEL_CORRECTNESS ) {
1355
1356 int numQubits = getNumCachedQubits();
1357 int numTerms = GENERATE_COPY( 1, 2, 10 );
1358
1359 PauliStrSum sum = createRandomPauliStrSum(numQubits, numTerms);
1360
1361 auto testFunc = [&](Qureg qureg, auto& ref) {
1362
1363 // must use (and ergo make) an identically-deployed workspace
1364 Qureg workspace = createCloneQureg(qureg);
1365 multiplyPauliStrSum(qureg, sum, workspace);
1366 destroyQureg(workspace);
1367
1368 ref = getMatrix(sum, numQubits) * ref;
1369 };
1370
1371 CAPTURE( numTerms );
1372 SECTION( LABEL_STATEVEC ) { TEST_ON_CACHED_QUREGS(statevecQuregs, statevecRef, testFunc); }
1373 SECTION( LABEL_DENSMATR ) { TEST_ON_CACHED_QUREGS(densmatrQuregs, densmatrRef, testFunc); }
1374 }
1375
1376 /// @todo input validation
1377}
1378
1379
1380/** @} (end defgroup) */
1381
1382
1383
1384/**
1385 * @todo
1386 * UNTESTED FUNCTIONS
1387 */
1388
1389void applyTrotterizedPauliStrSumGadget(Qureg qureg, PauliStrSum sum, qreal angle, int order, int reps);
void setValidationEpsilonToDefault()
Definition debug.cpp:96
void setValidationEpsilon(qreal eps)
Definition debug.cpp:88
QuESTEnv getQuESTEnv()
void initDebugState(Qureg qureg)
CompMatr createCompMatr(int numQubits)
Definition matrices.cpp:211
DiagMatr createDiagMatr(int numQubits)
Definition matrices.cpp:246
void destroyDiagMatr(DiagMatr matrix)
Definition matrices.cpp:403
void destroyCompMatr(CompMatr matrix)
Definition matrices.cpp:402
static CompMatr2 getCompMatr2(qcomp **in)
Definition matrices.h:320
static CompMatr1 getCompMatr1(qcomp **in)
Definition matrices.h:304
static DiagMatr2 getDiagMatr2(qcomp *in)
Definition matrices.h:352
static DiagMatr1 getDiagMatr1(qcomp *in)
Definition matrices.h:338
void setDiagMatr(DiagMatr out, qcomp *in)
Definition matrices.cpp:435
void setCompMatr(CompMatr matr, qcomp **vals)
Definition matrices.cpp:427
void setFullStateDiagMatr(FullStateDiagMatr out, qindex startInd, qcomp *in, qindex numElems)
Definition matrices.cpp:447
void multiplyCompMatr1(Qureg qureg, int target, CompMatr1 matrix)
void multiplyCompMatr2(Qureg qureg, int target1, int target2, CompMatr2 matr)
void multiplyCompMatr(Qureg qureg, int *targets, int numTargets, CompMatr matr)
void multiplyDiagMatr1(Qureg qureg, int target, DiagMatr1 matr)
void multiplyDiagMatr2(Qureg qureg, int target1, int target2, DiagMatr2 matr)
void multiplyDiagMatr(Qureg qureg, int *targets, int numTargets, DiagMatr matrix)
void multiplyDiagMatrPower(Qureg qureg, int *targets, int numTargets, DiagMatr matrix, qcomp exponent)
void applyFullStateDiagMatr(Qureg qureg, FullStateDiagMatr matrix)
void multiplyFullStateDiagMatrPower(Qureg qureg, FullStateDiagMatr matrix, qcomp exponent)
void applyFullStateDiagMatrPower(Qureg qureg, FullStateDiagMatr matrix, qcomp exponent)
void multiplyFullStateDiagMatr(Qureg qureg, FullStateDiagMatr matrix)
qreal applyForcedQubitMeasurement(Qureg qureg, int target, int outcome)
qindex applyMultiQubitMeasurement(Qureg qureg, int *qubits, int numQubits)
int applyQubitMeasurement(Qureg qureg, int target)
qreal applyForcedMultiQubitMeasurement(Qureg qureg, int *qubits, int *outcomes, int numQubits)
qindex applyMultiQubitMeasurementAndGetProb(Qureg qureg, int *qubits, int numQubits, qreal *probability)
int applyQubitMeasurementAndGetProb(Qureg qureg, int target, qreal *probability)
void multiplyMultiQubitNot(Qureg qureg, int *targets, int numTargets)
void multiplyPauliGadget(Qureg qureg, PauliStr str, qreal angle)
void multiplyPauliStr(Qureg qureg, PauliStr str)
void multiplyPauliStrSum(Qureg qureg, PauliStrSum sum, Qureg workspace)
void applyTrotterizedPauliStrSumGadget(Qureg qureg, PauliStrSum sum, qreal angle, int order, int reps)
void applyMultiQubitPhaseShift(Qureg qureg, int *targets, int numTargets, qreal angle)
void applyTwoQubitPhaseShift(Qureg qureg, int target1, int target2, qreal angle)
void applyTwoQubitPhaseFlip(Qureg qureg, int target1, int target2)
void applyPhaseShift(Qureg qureg, int target, qreal angle)
void applyPhaseFlip(Qureg qureg, int target)
void multiplyPhaseGadget(Qureg qureg, int *targets, int numTargets, qreal angle)
void applyMultiQubitPhaseFlip(Qureg qureg, int *targets, int numTargets)
void applyMultiQubitProjector(Qureg qureg, int *qubits, int *outcomes, int numQubits)
void applyQubitProjector(Qureg qureg, int target, int outcome)
void applyFullQuantumFourierTransform(Qureg qureg)
void applyQuantumFourierTransform(Qureg qureg, int *targets, int numTargets)
Qureg createCloneQureg(Qureg qureg)
Definition qureg.cpp:313
void destroyQureg(Qureg qureg)
Definition qureg.cpp:328
qmatrix getKroneckerProduct(qmatrix a, qmatrix b)
Definition linalg.cpp:523
qmatrix getExponentialOfPauliMatrix(qreal arg, qmatrix m)
Definition linalg.cpp:348
qmatrix getZeroMatrix(size_t dim)
Definition qmatrix.cpp:18
qcomp getRandomComplex()
Definition random.cpp:85
qmatrix getRandomUnitary(int numQb)
Definition random.cpp:326
qreal getRandomReal(qreal min, qreal maxIncl)
Definition random.cpp:60
vector< qreal > getRandomProbabilities(int numProbs)
Definition random.cpp:138
int getRandomInt(int min, int maxExcl)
Definition random.cpp:76
TEST_CASE("calcExpecPauliStr", TEST_CATEGORY)
TEST_ANY_CTRL_OPERATION(PauliStr, any, paulistr, nullptr)
Definition qureg.h:42