test_operators.cpp
Go to the documentation of this file.
1 
2 #include "catch.hpp"
3 #include "QuEST.h"
4 #include "utilities.hpp"
5 
10 #define PREPARE_TEST(quregVec, quregMatr, refVec, refMatr) \
11  Qureg quregVec = createQureg(NUM_QUBITS, QUEST_ENV); \
12  Qureg quregMatr = createDensityQureg(NUM_QUBITS, QUEST_ENV); \
13  initDebugState(quregVec); \
14  initDebugState(quregMatr); \
15  QVector refVec = toQVector(quregVec); \
16  QMatrix refMatr = toQMatrix(quregMatr);
17 
19 #define CLEANUP_TEST(quregVec, quregMatr) \
20  destroyQureg(quregVec, QUEST_ENV); \
21  destroyQureg(quregMatr, QUEST_ENV);
22 
23 /* allows concise use of Contains in catch's REQUIRE_THROWS_WITH */
24 using Catch::Matchers::Contains;
25 
26 
27 
32 TEST_CASE( "applyDiagonalOp", "[operators]" ) {
33 
34  PREPARE_TEST( quregVec, quregMatr, refVec, refMatr );
35 
36  SECTION( "correctness" ) {
37 
38  // try 10 random operators
39  GENERATE( range(0,10) );
40 
41  // make a totally random (non-Hermitian) diagonal oeprator
43  for (long long int i=0; i<op.numElemsPerChunk; i++) {
44  op.real[i] = getRandomReal(-5, 5);
45  op.imag[i] = getRandomReal(-5, 5);
46  }
47  syncDiagonalOp(op);
48 
49  SECTION( "state-vector" ) {
50 
51  QVector ref = toQMatrix(op) * refVec;
52  applyDiagonalOp(quregVec, op);
53  REQUIRE( areEqual(quregVec, ref) );
54  }
55  SECTION( "density-matrix" ) {
56 
57  QMatrix ref = toQMatrix(op) * refMatr;
58  applyDiagonalOp(quregMatr, op);
59  REQUIRE( areEqual(quregMatr, ref, 100*REAL_EPS) );
60  }
61 
63  }
64  SECTION( "input validation" ) {
65 
66  SECTION( "mismatching size" ) {
67 
69 
70  REQUIRE_THROWS_WITH( applyDiagonalOp(quregVec, op), Contains("equal number of qubits"));
71  REQUIRE_THROWS_WITH( applyDiagonalOp(quregMatr, op), Contains("equal number of qubits"));
72 
74  }
75  }
76  CLEANUP_TEST( quregVec, quregMatr );
77 }
78 
79 
80 
85 TEST_CASE( "applyMatrix2", "[operators]" ) {
86 
87  PREPARE_TEST( quregVec, quregMatr, refVec, refMatr );
88 
89  // every test will use a unique random matrix
90  QMatrix op = getRandomQMatrix(2); // 2-by-2
91  ComplexMatrix2 matr = toComplexMatrix2(op);
92 
93  SECTION( "correctness" ) {
94 
95  int target = GENERATE( range(0,NUM_QUBITS) );
96 
97  // reference boilerplate
98  int* ctrls = NULL;
99  int numCtrls = 0;
100  int targs[] = {target};
101  int numTargs = 1;
102 
103  SECTION( "state-vector" ) {
104 
105  applyMatrix2(quregVec, target, matr);
106  applyReferenceMatrix(refVec, ctrls, numCtrls, targs, numTargs, op);
107 
108  REQUIRE( areEqual(quregVec, refVec) );
109  }
110  SECTION( "density-matrix" ) {
111 
112  applyMatrix2(quregMatr, target, matr);
113  applyReferenceMatrix(refMatr, ctrls, numCtrls, targs, numTargs, op);
114 
115  REQUIRE( areEqual(quregMatr, refMatr, 10*REAL_EPS) );
116  }
117  }
118  SECTION( "input validation" ) {
119 
120  SECTION( "qubit indices" ) {
121 
122  int target = GENERATE( -1, NUM_QUBITS );
123  REQUIRE_THROWS_WITH( applyMatrix2(quregVec, target, matr), Contains("Invalid target") );
124  }
125  }
126  CLEANUP_TEST( quregVec, quregMatr );
127 }
128 
129 
130 
135 TEST_CASE( "applyMatrix4", "[operators]" ) {
136 
137  PREPARE_TEST( quregVec, quregMatr, refVec, refMatr );
138 
139  // in distributed mode, each node must be able to fit all amps modified by matrix
140  REQUIRE( quregVec.numAmpsPerChunk >= 4 );
141 
142  // every test will use a unique random matrix
143  QMatrix op = getRandomQMatrix(4); // 4-by-4
144  ComplexMatrix4 matr = toComplexMatrix4(op);
145 
146  SECTION( "correctness" ) {
147 
148  int targ1 = GENERATE( range(0,NUM_QUBITS) );
149  int targ2 = GENERATE_COPY( filter([=](int t){ return t!=targ1; }, range(0,NUM_QUBITS)) );
150 
151  // reference boilerplate
152  int* ctrls = NULL;
153  int numCtrls = 0;
154  int targs[] = {targ1, targ2};
155  int numTargs = 2;
156 
157  SECTION( "state-vector" ) {
158 
159  applyMatrix4(quregVec, targ1, targ2, matr);
160  applyReferenceMatrix(refVec, ctrls, numCtrls, targs, numTargs, op);
161  REQUIRE( areEqual(quregVec, refVec) );
162  }
163  SECTION( "density-matrix" ) {
164 
165  applyMatrix4(quregMatr, targ1, targ2, matr);
166  applyReferenceMatrix(refMatr, ctrls, numCtrls, targs, numTargs, op);
167  REQUIRE( areEqual(quregMatr, refMatr, 10*REAL_EPS) );
168  }
169  }
170  SECTION( "input validation" ) {
171 
172  SECTION( "qubit indices" ) {
173 
174  int targ1 = GENERATE( -1, NUM_QUBITS );
175  int targ2 = 0;
176  REQUIRE_THROWS_WITH( applyMatrix4(quregVec, targ1, targ2, matr), Contains("Invalid target") );
177  REQUIRE_THROWS_WITH( applyMatrix4(quregVec, targ2, targ1, matr), Contains("Invalid target") );
178  }
179  SECTION( "repetition of targets" ) {
180 
181  int qb = 0;
182  REQUIRE_THROWS_WITH( applyMatrix4(quregVec, qb, qb, matr), Contains("target") && Contains("unique") );
183  }
184  SECTION( "matrix fits in node" ) {
185 
186  // pretend we have a very limited distributed memory
187  quregVec.numAmpsPerChunk = 1;
188  REQUIRE_THROWS_WITH( applyMatrix4(quregVec, 0, 1, matr), Contains("targets too many qubits"));
189  }
190  }
191  CLEANUP_TEST( quregVec, quregMatr );
192 }
193 
194 
195 
200 TEST_CASE( "applyMatrixN", "[operators]" ) {
201 
202  PREPARE_TEST( quregVec, quregMatr, refVec, refMatr );
203 
204  // figure out max-num (inclusive) targs allowed by hardware backend
205  int maxNumTargs = calcLog2(quregVec.numAmpsPerChunk);
206 
207  SECTION( "correctness" ) {
208 
209  // generate all possible qubit arrangements
210  int numTargs = GENERATE_COPY( range(1,maxNumTargs+1) ); // inclusive upper bound
211  int* targs = GENERATE_COPY( sublists(range(0,NUM_QUBITS), numTargs) );
212 
213  // for each qubit arrangement, use a new random matrix
214  QMatrix op = getRandomQMatrix(1 << numTargs);
215  ComplexMatrixN matr = createComplexMatrixN(numTargs);
216  toComplexMatrixN(op, matr);
217 
218  // reference boilerplate
219  int* ctrls = NULL;
220  int numCtrls = 0;
221 
222  SECTION( "state-vector" ) {
223 
224  applyMatrixN(quregVec, targs, numTargs, matr);
225  applyReferenceMatrix(refVec, ctrls, numCtrls, targs, numTargs, op);
226  REQUIRE( areEqual(quregVec, refVec) );
227  }
228  SECTION( "density-matrix" ) {
229 
230  applyMatrixN(quregMatr, targs, numTargs, matr);
231  applyReferenceMatrix(refMatr, ctrls, numCtrls, targs, numTargs, op);
232  REQUIRE( areEqual(quregMatr, refMatr, 100*REAL_EPS) );
233  }
234  destroyComplexMatrixN(matr);
235  }
236  SECTION( "input validation" ) {
237 
238  SECTION( "number of targets" ) {
239 
240  // there cannot be more targets than qubits in register
241  int numTargs = GENERATE( -1, 0, NUM_QUBITS+1 );
242  int targs[NUM_QUBITS+1]; // prevents seg-fault if validation doesn't trigger
243  ComplexMatrixN matr = createComplexMatrixN(NUM_QUBITS+1); // prevent seg-fault
244 
245  REQUIRE_THROWS_WITH( applyMatrixN(quregVec, targs, numTargs, matr), Contains("Invalid number of target"));
246  destroyComplexMatrixN(matr);
247  }
248  SECTION( "repetition in targets" ) {
249 
250  int numTargs = 3;
251  int targs[] = {1,2,2};
252  ComplexMatrixN matr = createComplexMatrixN(numTargs); // prevents seg-fault if validation doesn't trigger
253 
254  REQUIRE_THROWS_WITH( applyMatrixN(quregVec, targs, numTargs, matr), Contains("target") && Contains("unique"));
255  destroyComplexMatrixN(matr);
256  }
257  SECTION( "qubit indices" ) {
258 
259  int numTargs = 3;
260  int targs[] = {1,2,3};
261  ComplexMatrixN matr = createComplexMatrixN(numTargs); // prevents seg-fault if validation doesn't trigger
262 
263  int inv = GENERATE( -1, NUM_QUBITS );
264  targs[GENERATE_COPY( range(0,numTargs) )] = inv; // make invalid target
265  REQUIRE_THROWS_WITH( applyMatrixN(quregVec, targs, numTargs, matr), Contains("Invalid target") );
266 
267  destroyComplexMatrixN(matr);
268  }
269  SECTION( "matrix creation" ) {
270 
271  int numTargs = 3;
272  int targs[] = {1,2,3};
273 
274  /* compilers don't auto-initialise to NULL; the below circumstance
275  * only really occurs when 'malloc' returns NULL in createComplexMatrixN,
276  * which actually triggers its own validation. Hence this test is useless
277  * currently.
278  */
279  ComplexMatrixN matr;
280  matr.real = NULL;
281  matr.imag = NULL;
282  REQUIRE_THROWS_WITH( applyMatrixN(quregVec, targs, numTargs, matr), Contains("created") );
283  }
284  SECTION( "matrix dimensions" ) {
285 
286  int targs[2] = {1,2};
287  ComplexMatrixN matr = createComplexMatrixN(3); // intentionally wrong size
288 
289  REQUIRE_THROWS_WITH( applyMatrixN(quregVec, targs, 2, matr), Contains("matrix size"));
290  destroyComplexMatrixN(matr);
291  }
292  SECTION( "matrix fits in node" ) {
293 
294  // pretend we have a very limited distributed memory (judged by matr size)
295  quregVec.numAmpsPerChunk = 1;
296  int qb[] = {1,2};
297  ComplexMatrixN matr = createComplexMatrixN(2); // prevents seg-fault if validation doesn't trigger
298  REQUIRE_THROWS_WITH( applyMatrixN(quregVec, qb, 2, matr), Contains("targets too many qubits"));
299  destroyComplexMatrixN(matr);
300  }
301  }
302  CLEANUP_TEST( quregVec, quregMatr );
303 }
304 
305 
306 
311 TEST_CASE( "applyMultiControlledMatrixN", "[operators]" ) {
312 
313  PREPARE_TEST( quregVec, quregMatr, refVec, refMatr );
314 
315  // figure out max-num targs (inclusive) allowed by hardware backend
316  int maxNumTargs = calcLog2(quregVec.numAmpsPerChunk);
317  if (maxNumTargs >= NUM_QUBITS)
318  maxNumTargs = NUM_QUBITS - 1; // leave room for min-number of control qubits
319 
320  SECTION( "correctness" ) {
321 
322  // try all possible numbers of targets and controls
323  int numTargs = GENERATE_COPY( range(1,maxNumTargs+1) );
324  int maxNumCtrls = NUM_QUBITS - numTargs;
325  int numCtrls = GENERATE_COPY( range(1,maxNumCtrls+1) );
326 
327  // generate all possible valid qubit arrangements
328  int* targs = GENERATE_COPY( sublists(range(0,NUM_QUBITS), numTargs) );
329  int* ctrls = GENERATE_COPY( sublists(range(0,NUM_QUBITS), numCtrls, targs, numTargs) );
330 
331  // for each qubit arrangement, use a new random unitary
332  QMatrix op = getRandomQMatrix(1 << numTargs);
333  ComplexMatrixN matr = createComplexMatrixN(numTargs);
334  toComplexMatrixN(op, matr);
335 
336  SECTION( "state-vector" ) {
337 
338  applyMultiControlledMatrixN(quregVec, ctrls, numCtrls, targs, numTargs, matr);
339  applyReferenceMatrix(refVec, ctrls, numCtrls, targs, numTargs, op);
340  REQUIRE( areEqual(quregVec, refVec) );
341  }
342  SECTION( "density-matrix" ) {
343 
344  applyMultiControlledMatrixN(quregMatr, ctrls, numCtrls, targs, numTargs, matr);
345  applyReferenceMatrix(refMatr, ctrls, numCtrls, targs, numTargs, op);
346  REQUIRE( areEqual(quregMatr, refMatr, 100*REAL_EPS) );
347  }
348  destroyComplexMatrixN(matr);
349  }
350  SECTION( "input validation" ) {
351 
352  SECTION( "number of targets" ) {
353 
354  // there cannot be more targets than qubits in register
355  // (numTargs=NUM_QUBITS is caught elsewhere, because that implies ctrls are invalid)
356  int numTargs = GENERATE( -1, 0, NUM_QUBITS+1 );
357  int targs[NUM_QUBITS+1]; // prevents seg-fault if validation doesn't trigger
358  int ctrls[] = {0};
359  ComplexMatrixN matr = createComplexMatrixN(NUM_QUBITS+1); // prevent seg-fault
360  toComplexMatrixN(getRandomQMatrix( 1 << (NUM_QUBITS+1)), matr);
361 
362  REQUIRE_THROWS_WITH( applyMultiControlledMatrixN(quregVec, ctrls, 1, targs, numTargs, matr), Contains("Invalid number of target"));
363  destroyComplexMatrixN(matr);
364  }
365  SECTION( "repetition in targets" ) {
366 
367  int ctrls[] = {0};
368  int numTargs = 3;
369  int targs[] = {1,2,2};
370  ComplexMatrixN matr = createComplexMatrixN(numTargs); // prevents seg-fault if validation doesn't trigger
371  toComplexMatrixN(getRandomQMatrix(1 << numTargs), matr);
372 
373  REQUIRE_THROWS_WITH( applyMultiControlledMatrixN(quregVec, ctrls, 1, targs, numTargs, matr), Contains("target") && Contains("unique"));
374  destroyComplexMatrixN(matr);
375  }
376  SECTION( "number of controls" ) {
377 
378  int numCtrls = GENERATE( -1, 0, NUM_QUBITS, NUM_QUBITS+1 );
379  int ctrls[NUM_QUBITS+1]; // avoids seg-fault if validation not triggered
380  int targs[1] = {0};
382  toComplexMatrixN(getRandomQMatrix(1 << 1), matr);
383 
384  REQUIRE_THROWS_WITH( applyMultiControlledMatrixN(quregVec, ctrls, numCtrls, targs, 1, matr), Contains("Invalid number of control"));
385  destroyComplexMatrixN(matr);
386  }
387  SECTION( "repetition in controls" ) {
388 
389  int ctrls[] = {0,1,1};
390  int targs[] = {3};
392  toComplexMatrixN(getRandomQMatrix(1 << 1), matr);
393 
394  REQUIRE_THROWS_WITH( applyMultiControlledMatrixN(quregVec, ctrls, 3, targs, 1, matr), Contains("control") && Contains("unique"));
395  destroyComplexMatrixN(matr);
396  }
397  SECTION( "control and target collision" ) {
398 
399  int ctrls[] = {0,1,2};
400  int targs[] = {3,1,4};
402  toComplexMatrixN(getRandomQMatrix(1 << 3), matr);
403 
404  REQUIRE_THROWS_WITH( applyMultiControlledMatrixN(quregVec, ctrls, 3, targs, 3, matr), Contains("Control") && Contains("target") && Contains("disjoint"));
405  destroyComplexMatrixN(matr);
406  }
407  SECTION( "qubit indices" ) {
408 
409  // valid inds
410  int numQb = 2;
411  int qb1[2] = {0,1};
412  int qb2[2] = {2,3};
413  ComplexMatrixN matr = createComplexMatrixN(numQb);
414  toComplexMatrixN(getRandomQMatrix(1 << numQb), matr);
415 
416  // make qb1 invalid
417  int inv = GENERATE( -1, NUM_QUBITS );
418  qb1[GENERATE_COPY(range(0,numQb))] = inv;
419 
420  REQUIRE_THROWS_WITH( applyMultiControlledMatrixN(quregVec, qb1, numQb, qb2, numQb, matr), Contains("Invalid control") );
421  REQUIRE_THROWS_WITH( applyMultiControlledMatrixN(quregVec, qb2, numQb, qb1, numQb, matr), Contains("Invalid target") );
422  destroyComplexMatrixN(matr);
423  }
424  SECTION( "matrix creation" ) {
425 
426  int ctrls[1] = {0};
427  int targs[3] = {1,2,3};
428 
429  /* compilers don't auto-initialise to NULL; the below circumstance
430  * only really occurs when 'malloc' returns NULL in createComplexMatrixN,
431  * which actually triggers its own validation. Hence this test is useless
432  * currently.
433  */
434  ComplexMatrixN matr;
435  matr.real = NULL;
436  matr.imag = NULL;
437  REQUIRE_THROWS_WITH( applyMultiControlledMatrixN(quregVec, ctrls, 1, targs, 3, matr), Contains("created") );
438  }
439  SECTION( "matrix dimensions" ) {
440 
441  int ctrls[1] = {0};
442  int targs[2] = {1,2};
443  ComplexMatrixN matr = createComplexMatrixN(3); // intentionally wrong size
444  toComplexMatrixN(getRandomQMatrix(1 << 3), matr);
445 
446  REQUIRE_THROWS_WITH( applyMultiControlledMatrixN(quregVec, ctrls, 1, targs, 2, matr), Contains("matrix size"));
447  destroyComplexMatrixN(matr);
448  }
449  SECTION( "matrix fits in node" ) {
450 
451  // pretend we have a very limited distributed memory (judged by matr size)
452  quregVec.numAmpsPerChunk = 1;
453  int ctrls[1] = {0};
454  int targs[2] = {1,2};
456  toComplexMatrixN(getRandomQMatrix(1 << 2), matr);
457 
458  REQUIRE_THROWS_WITH( applyMultiControlledMatrixN(quregVec, ctrls, 1, targs, 2, matr), Contains("targets too many qubits"));
459  destroyComplexMatrixN(matr);
460  }
461  }
462  CLEANUP_TEST( quregVec, quregMatr );
463 }
464 
465 
466 
471 TEST_CASE( "applyPauliHamil", "[operators]" ) {
472 
477 
478  initDebugState(vecIn);
479  initDebugState(matIn);
480 
481  SECTION( "correctness" ) {
482 
483  /* it's too expensive to try every possible Pauli configuration, so
484  * we'll try 10 random codes, and for each, random coefficients
485  */
486  GENERATE( range(0,10) );
487 
488  int numTerms = GENERATE( 1, 2, 10, 15 );
489  PauliHamil hamil = createPauliHamil(NUM_QUBITS, numTerms);
490  setRandomPauliSum(hamil);
491  QMatrix refHamil = toQMatrix(hamil);
492 
493  SECTION( "state-vector" ) {
494 
495  QVector vecRef = toQVector(vecIn);
496  applyPauliHamil(vecIn, hamil, vecOut);
497 
498  // ensure vecIn barely changes under precision
499  REQUIRE( areEqual(vecIn, vecRef) );
500 
501  // ensure vecOut changed correctly
502  REQUIRE( areEqual(vecOut, refHamil * vecRef) );
503  }
504  SECTION( "density-matrix" ) {
505 
506  QMatrix matRef = toQMatrix(matIn);
507  applyPauliHamil(matIn, hamil, matOut);
508 
509  // ensure matIn barely changes under precision
510  REQUIRE( areEqual(matIn, matRef) );
511 
512  // ensure matOut changed correctly
513  REQUIRE( areEqual(matOut, refHamil * matRef, 1E2*REAL_EPS) );
514  }
515 
516  destroyPauliHamil(hamil);
517  }
518  SECTION( "input validation" ) {
519 
520  SECTION( "pauli codes" ) {
521 
522  int numTerms = 3;
523  PauliHamil hamil = createPauliHamil(NUM_QUBITS, numTerms);
524 
525  // make one pauli code wrong
526  hamil.pauliCodes[GENERATE_COPY( range(0,numTerms*NUM_QUBITS) )] = (pauliOpType) GENERATE( -1, 4 );
527  REQUIRE_THROWS_WITH( applyPauliHamil(vecIn, hamil, vecOut), Contains("Invalid Pauli code") );
528 
529  destroyPauliHamil(hamil);
530  }
531  SECTION( "qureg dimensions" ) {
532 
533  Qureg badVec = createQureg(NUM_QUBITS+1, QUEST_ENV);
535 
536  REQUIRE_THROWS_WITH( applyPauliHamil(vecIn, hamil, badVec), Contains("Dimensions of the qubit registers don't match") );
537 
538  destroyQureg(badVec, QUEST_ENV);
539  destroyPauliHamil(hamil);
540  }
541  SECTION( "qureg types" ) {
542 
544 
545  REQUIRE_THROWS_WITH( applyPauliHamil(vecIn, hamil, matOut), Contains("Registers must both be state-vectors or both be density matrices") );
546 
547  destroyPauliHamil(hamil);
548  }
549  SECTION( "matching hamiltonian qubits" ) {
550 
551  PauliHamil hamil = createPauliHamil(NUM_QUBITS + 1, 1);
552 
553  REQUIRE_THROWS_WITH( applyPauliHamil(vecIn, hamil, vecOut), Contains("same number of qubits") );
554  REQUIRE_THROWS_WITH( applyPauliHamil(matIn, hamil, matOut), Contains("same number of qubits") );
555 
556  destroyPauliHamil(hamil);
557  }
558  }
559  destroyQureg(vecIn, QUEST_ENV);
560  destroyQureg(vecOut, QUEST_ENV);
561  destroyQureg(matIn, QUEST_ENV);
562  destroyQureg(matOut, QUEST_ENV);
563 }
564 
565 
566 
571 TEST_CASE( "applyPauliSum", "[operators]" ) {
572 
577 
578  initDebugState(vecIn);
579  initDebugState(matIn);
580 
581  SECTION( "correctness" ) {
582 
583  /* it's too expensive to try ALL Pauli sequences, via
584  * pauliOpType* paulis = GENERATE_COPY( pauliseqs(numPaulis) );.
585  * Furthermore, take(10, pauliseqs(numTargs)) will try the same pauli codes.
586  * Hence, we instead opt to repeatedly randomly generate pauliseqs
587  */
588  GENERATE( range(0,10) ); // gen 10 random pauli-codes
589 
590  int numTerms = GENERATE( 1, 2, 10, 15);
591  int numPaulis = numTerms * NUM_QUBITS;
592 
593  // each test will use random coefficients
594  qreal coeffs[numTerms];
595  pauliOpType paulis[numPaulis];
596  setRandomPauliSum(coeffs, paulis, NUM_QUBITS, numTerms);
597  QMatrix pauliSum = toQMatrix(coeffs, paulis, NUM_QUBITS, numTerms);
598 
599  SECTION( "state-vector" ) {
600 
601  QVector vecRef = toQVector(vecIn);
602  applyPauliSum(vecIn, paulis, coeffs, numTerms, vecOut);
603 
604  // ensure vecIn barely changes under precision
605  REQUIRE( areEqual(vecIn, vecRef) );
606 
607  // ensure vecOut changed correctly
608  REQUIRE( areEqual(vecOut, pauliSum * vecRef) );
609  }
610  SECTION( "density-matrix" ) {
611 
612  QMatrix matRef = toQMatrix(matIn);
613  applyPauliSum(matIn, paulis, coeffs, numTerms, matOut);
614 
615  // ensure matIn barely changes under precision
616  REQUIRE( areEqual(matIn, matRef) );
617 
618  // ensure matOut changed correctly
619  REQUIRE( areEqual(matOut, pauliSum * matRef, 1E2*REAL_EPS) );
620  }
621  }
622  SECTION( "input validation" ) {
623 
624  SECTION( "number of terms" ) {
625 
626  int numTerms = GENERATE( -1, 0 );
627  REQUIRE_THROWS_WITH( applyPauliSum(vecIn, NULL, NULL, numTerms, vecOut), Contains("Invalid number of terms") );
628  }
629  SECTION( "pauli codes" ) {
630 
631  // valid codes
632  int numTerms = 3;
633  int numPaulis = numTerms*NUM_QUBITS;
634  pauliOpType paulis[numPaulis];
635  for (int i=0; i<numPaulis; i++)
636  paulis[i] = PAULI_I;
637 
638  // make one invalid
639  paulis[GENERATE_COPY( range(0,numPaulis) )] = (pauliOpType) GENERATE( -1, 4 );
640  REQUIRE_THROWS_WITH( applyPauliSum(vecIn, paulis, NULL, numTerms, vecOut), Contains("Invalid Pauli code") );
641  }
642  SECTION( "qureg dimensions" ) {
643 
644  Qureg badVec = createQureg(NUM_QUBITS+1, QUEST_ENV);
645  pauliOpType paulis[NUM_QUBITS];
646  REQUIRE_THROWS_WITH( applyPauliSum(vecIn, paulis, NULL, 1, badVec), Contains("Dimensions of the qubit registers don't match") );
647  destroyQureg(badVec, QUEST_ENV);
648  }
649  SECTION( "qureg types" ) {
650 
651  pauliOpType paulis[NUM_QUBITS];
652  REQUIRE_THROWS_WITH( applyPauliSum(vecIn, paulis, NULL, 1, matOut), Contains("Registers must both be state-vectors or both be density matrices") );
653  }
654  }
655  destroyQureg(vecIn, QUEST_ENV);
656  destroyQureg(vecOut, QUEST_ENV);
657  destroyQureg(matIn, QUEST_ENV);
658  destroyQureg(matOut, QUEST_ENV);
659 }
660 
661 
662 
667 TEST_CASE( "applyTrotterCircuit", "[operators]" ) {
668 
671  initDebugState(vec);
672  initDebugState(mat);
673 
674  Qureg vecRef = createCloneQureg(vec, QUEST_ENV);
675  Qureg matRef = createCloneQureg(mat, QUEST_ENV);
676 
677  SECTION( "correctness" ) {
678 
679  SECTION( "one term" ) {
680 
681  // a Hamiltonian with one term has an exact (trivial) Trotterisation
683 
684  // H = coeff X Y Z (on qubits 0,1,2)
685  qreal coeff = getRandomReal(-5, 5);
686  int numTargs = 3;
687  int targs[] = {0, 1, 2};
688  pauliOpType codes[] = {PAULI_X, PAULI_Y, PAULI_Z};
689  hamil.termCoeffs[0] = coeff;
690  for (int i=0; i<numTargs; i++)
691  hamil.pauliCodes[targs[i]] = codes[i];
692 
693  // time can be negative
694  qreal time = getRandomReal(-2,2);
695 
696  // by commutation, all reps & orders yield the same total unitary
697  int reps = GENERATE( range(1,5) );
698 
699  // applyTrotter(t, 1, 1) = exp(-i t coeff paulis)
700  // multiRotatePauli(a) = exp(- i a / 2 paulis)
701 
702  SECTION( "state-vector" ) {
703 
704  int order = GENERATE( 1, 2, 4 );
705 
706  applyTrotterCircuit(vec, hamil, time, order, reps);
707  multiRotatePauli(vecRef, targs, codes, numTargs, 2*time*coeff);
708  REQUIRE( areEqual(vec, vecRef, 10*REAL_EPS) );
709  }
710  SECTION( "density-matrix" ) {
711 
712  int order = GENERATE( 1, 2 ); // precision bites density-matrices quickly
713 
714  applyTrotterCircuit(mat, hamil, time, order, reps);
715  multiRotatePauli(matRef, targs, codes, numTargs, 2*time*coeff);
716  REQUIRE( areEqual(mat, matRef, 1E2*REAL_EPS) );
717  }
718 
719  destroyPauliHamil(hamil);
720  }
721  SECTION( "commuting terms" ) {
722 
723  // a Hamiltonian of commuting terms, Trotterises exactly
725 
726  // H = c0 X Y I + c1 X I Z (on qubits 0,1,2)
727  int targs[] = {0, 1, 2};
728  hamil.pauliCodes[0] = PAULI_X;
729  hamil.pauliCodes[0 + NUM_QUBITS] = PAULI_X;
730  hamil.pauliCodes[1] = PAULI_Y;
731  hamil.pauliCodes[1 + NUM_QUBITS] = PAULI_I;
732  hamil.pauliCodes[2] = PAULI_I;
733  hamil.pauliCodes[2 + NUM_QUBITS] = PAULI_Z;
734  for (int i=0; i<hamil.numSumTerms; i++)
735  hamil.termCoeffs[i] = getRandomReal(-5,5);
736 
737  // time can be negative
738  qreal time = getRandomReal(-2,2);
739 
740  // applyTrotter(t, 1, 1) = exp(-i t c0 paulis0) exp(-i t c1 paulis1)
741  // multiRotatePauli(a) = exp(- i a / 2 paulis)
742 
743  SECTION( "state-vector" ) {
744 
745  int reps = GENERATE( range(1,5) );
746  int order = GENERATE( 1, 2, 4 );
747 
748  applyTrotterCircuit(vec, hamil, time, order, reps);
749  multiRotatePauli(vecRef, targs, hamil.pauliCodes, 3, 2*time*hamil.termCoeffs[0]);
750  multiRotatePauli(vecRef, targs, &(hamil.pauliCodes[NUM_QUBITS]), 3, 2*time*hamil.termCoeffs[1]);
751  REQUIRE( areEqual(vec, vecRef, 10*REAL_EPS) );
752  }
753  SECTION( "density-matrix" ) {
754 
755  int reps = GENERATE( range(1,5) );
756  int order = GENERATE( 1, 2 ); // precision hurts density matrices quickly
757 
758  applyTrotterCircuit(mat, hamil, time, order, reps);
759  multiRotatePauli(matRef, targs, hamil.pauliCodes, 3, 2*time*hamil.termCoeffs[0]);
760  multiRotatePauli(matRef, targs, &(hamil.pauliCodes[NUM_QUBITS]), 3, 2*time*hamil.termCoeffs[1]);
761  REQUIRE( areEqual(mat, matRef, 1E2*REAL_EPS) );
762  }
763 
764  destroyPauliHamil(hamil);
765  }
766  SECTION( "general" ) {
767 
768  /* We'll consider an analytic time-evolved state, so that we can avoid
769  * comparing applyTrotterCircuit to other numerical approximations.
770  * We can construct such a state, by using a Hamiltonian with known
771  * analytic eigenvalues, and hence a known period. Time evolution of the
772  * period will just yield the input state.
773  *
774  * E.g. H = pi sqrt(2) X Y Z X Y + pi Y Z X Y Z + pi Z X Y Z X
775  * has (degenerate) eigenvalues +- 2 pi, so the period
776  * of the Hamiltonian is t=1.
777  */
778 
779  // hardcoded 5 qubits here in the Pauli codes
780  REQUIRE( NUM_QUBITS == 5 );
781 
782  // (verbose for C++ compatibility)
784  qreal coeffs[] = {M_PI * sqrt(2.0), M_PI, M_PI};
785  enum pauliOpType codes[] = {
789  initPauliHamil(hamil, coeffs, codes);
790 
791  // evolving to t=1 should leave the input state unchanged
792  qreal time = 1;
793 
794  // since unnormalised (initDebugState), max fid is 728359.8336
795  qreal fidNorm = 728359.8336;
796 
797  SECTION( "absolute" ) {
798 
799  // such a high order and reps should yield precise solution
800  int order = 4;
801  int reps = 20;
802  applyTrotterCircuit(vec, hamil, time, 4, 20);
803  qreal fid = calcFidelity(vec, vecRef) / fidNorm;
804 
805  REQUIRE( fid == Approx(1).epsilon(1E-8) );
806  }
807  SECTION( "repetitions scaling" ) {
808 
809  // exclude order 1; too few reps for monotonic increase of accuracy
810  int order = GENERATE( 2, 4, 6 );
811 
812  // accuracy should increase with increasing repetitions
813  int reps[] = {1, 5, 10};
814 
815  qreal prevFid = 0;
816  for (int i=0; i<3; i++) {
817  initDebugState(vec);
818  applyTrotterCircuit(vec, hamil, time, order, reps[i]);
819  qreal fid = calcFidelity(vec, vecRef) / fidNorm;
820 
821  REQUIRE( fid >= prevFid );
822  prevFid = fid;
823  }
824  }
825  SECTION( "order scaling" ) {
826 
827  // exclude order 1; too few reps for monotonic increase of accuracy
828  int reps = GENERATE( 5, 10 );
829 
830  // accuracy should increase with increasing repetitions
831  int orders[] = {1, 2, 4, 6};
832 
833  qreal prevFid = 0;
834  for (int i=0; i<4; i++) {
835  initDebugState(vec);
836  applyTrotterCircuit(vec, hamil, time, orders[i], reps);
837  qreal fid = calcFidelity(vec, vecRef) / fidNorm;
838 
839  REQUIRE( fid >= prevFid );
840  prevFid = fid;
841  }
842  }
843 
844  destroyPauliHamil(hamil);
845  }
846  }
847  SECTION( "input validation" ) {
848 
849  SECTION( "repetitions" ) {
850 
852  int reps = GENERATE( -1, 0 );
853 
854  REQUIRE_THROWS_WITH( applyTrotterCircuit(vec, hamil, 1, 1, reps), Contains("repetitions must be >=1") );
855 
856  destroyPauliHamil(hamil);
857  }
858  SECTION( "order" ) {
859 
861  int order = GENERATE( -1, 0, 3, 5, 7 );
862 
863  REQUIRE_THROWS_WITH( applyTrotterCircuit(vec, hamil, 1, order, 1), Contains("order must be 1, or an even number") );
864 
865  destroyPauliHamil(hamil);
866  }
867  SECTION( "pauli codes" ) {
868 
869  int numTerms = 3;
870  PauliHamil hamil = createPauliHamil(NUM_QUBITS, numTerms);
871 
872  // make one pauli code wrong
873  hamil.pauliCodes[GENERATE_COPY( range(0,numTerms*NUM_QUBITS) )] = (pauliOpType) GENERATE( -1, 4 );
874  REQUIRE_THROWS_WITH( applyTrotterCircuit(vec, hamil, 1, 1, 1), Contains("Invalid Pauli code") );
875 
876  destroyPauliHamil(hamil);
877  }
878  SECTION( "matching hamiltonian qubits" ) {
879 
880  PauliHamil hamil = createPauliHamil(NUM_QUBITS + 1, 1);
881 
882  REQUIRE_THROWS_WITH( applyTrotterCircuit(vec, hamil, 1, 1, 1), Contains("same number of qubits") );
883  REQUIRE_THROWS_WITH( applyTrotterCircuit(mat, hamil, 1, 1, 1), Contains("same number of qubits") );
884 
885  destroyPauliHamil(hamil);
886  }
887  }
888 
889  destroyQureg(vec, QUEST_ENV);
890  destroyQureg(mat, QUEST_ENV);
891  destroyQureg(vecRef, QUEST_ENV);
892  destroyQureg(matRef, QUEST_ENV);
893 }
894 
pauliOpType
Codes for specifying Pauli operators.
Definition: QuEST.h:96
void applyMultiControlledMatrixN(Qureg qureg, int *ctrls, int numCtrls, int *targs, int numTargs, ComplexMatrixN u)
Apply a general N-by-N matrix, which may be non-unitary, with additional controlled qubits.
Definition: QuEST.c:874
void initPauliHamil(PauliHamil hamil, qreal *coeffs, enum pauliOpType *codes)
Initialise a PauliHamil instance with the given term coefficients and Pauli codes (one for every qubi...
Definition: QuEST.c:1253
QuESTEnv QUEST_ENV
The global QuESTEnv instance, to be created and destroyed once in this main(), so that the MPI enviro...
Definition: main.cpp:20
@ PAULI_Z
Definition: QuEST.h:96
void destroyComplexMatrixN(ComplexMatrixN m)
Destroy a ComplexMatrixN instance created with createComplexMatrixN()
Definition: QuEST.c:1120
@ PAULI_I
Definition: QuEST.h:96
ComplexMatrixN createComplexMatrixN(int numQubits)
Create (dynamically) a square complex matrix which can be passed to the multi-qubit general unitary f...
Definition: QuEST.c:1099
void syncDiagonalOp(DiagonalOp op)
Copy the elements in DiagonalOp op.real and op.imag to the persisent GPU memory.
Definition: QuEST.c:1280
void destroyDiagonalOp(DiagonalOp op, QuESTEnv env)
Destroys a DiagonalOp created with createDiagonalOp(), freeing its memory.
Definition: QuEST.c:1273
#define NUM_QUBITS
The default number of qubits in the registers created for unit testing (both statevectors and density...
Definition: utilities.hpp:36
#define CLEANUP_TEST(quregVec, quregMatr)
Destroys the data structures made by PREPARE_TEST.
qreal calcFidelity(Qureg qureg, Qureg pureState)
Calculates the fidelity of qureg (a statevector or density matrix) against a reference pure state (ne...
Definition: QuEST.c:942
qreal getRandomReal(qreal min, qreal max)
Returns a random real between min (inclusive) and max (exclusive), from the uniform distribution.
Definition: utilities.cpp:410
Represents a 4x4 matrix of complex numbers.
Definition: QuEST.h:125
void applyPauliHamil(Qureg inQureg, PauliHamil hamil, Qureg outQureg)
Modifies outQureg to be the result of applying PauliHamil (a Hermitian but not necessarily unitary op...
Definition: QuEST.c:819
Represents a general 2^N by 2^N matrix of complex numbers.
Definition: QuEST.h:136
#define qreal
void multiRotatePauli(Qureg qureg, int *targetQubits, enum pauliOpType *targetPaulis, int numTargets, qreal angle)
Apply a multi-qubit multi-Pauli rotation on a selected number of qubits.
Definition: QuEST.c:642
QMatrix toQMatrix(ComplexMatrix2 src)
Returns a copy of the given 2-by-2 matrix.
Definition: utilities.cpp:869
unsigned int calcLog2(long unsigned int num)
returns log2 of numbers which must be gauranteed to be 2^n
@ PAULI_X
Definition: QuEST.h:96
TEST_CASE("applyDiagonalOp", "[operators]")
void setRandomPauliSum(qreal *coeffs, pauliOpType *codes, int numQubits, int numTerms)
Populates the coeffs array with random qreals in (-5, 5), and populates codes with random Pauli codes...
Definition: utilities.cpp:1054
void toComplexMatrixN(QMatrix qm, ComplexMatrixN cm)
Initialises cm with the values of qm.
Definition: utilities.cpp:858
std::vector< qcomp > QVector
A complex vector, which can be zero-initialised with QVector(numAmps).
Definition: utilities.hpp:60
qreal * imag
The imaginary values of the 2^numQubits complex elements.
Definition: QuEST.h:191
QVector toQVector(Qureg qureg)
Returns an equal-size copy of the given state-vector qureg.
Definition: utilities.cpp:938
void applyDiagonalOp(Qureg qureg, DiagonalOp op)
Apply a diagonal complex operator, which is possibly non-unitary and non-Hermitian,...
Definition: QuEST.c:887
qreal * termCoeffs
The coefficient of each Pauli product. This is a length numSumTerms array.
Definition: QuEST.h:164
void applyMatrix4(Qureg qureg, int targetQubit1, int targetQubit2, ComplexMatrix4 u)
Apply a general 4-by-4 matrix, which may be non-unitary.
Definition: QuEST.c:853
enum pauliOpType * pauliCodes
The Pauli operators acting on each qubit, flattened over every operator.
Definition: QuEST.h:162
ComplexMatrix4 toComplexMatrix4(QMatrix qm)
Returns a ComplexMatrix4 copy of QMatix qm.
Definition: utilities.cpp:852
int numSumTerms
The number of terms in the weighted sum, or the number of Pauli products.
Definition: QuEST.h:166
Represents a diagonal complex operator on the full Hilbert state of a Qureg.
Definition: QuEST.h:178
@ PAULI_Y
Definition: QuEST.h:96
Represents a weighted sum of pauli products.
Definition: QuEST.h:158
void destroyQureg(Qureg qureg, QuESTEnv env)
Deallocate a Qureg object representing a set of qubits.
Definition: QuEST.c:77
qreal ** real
Definition: QuEST.h:139
void initDebugState(Qureg qureg)
Initialises qureg to be in the un-normalised, non-physical state with with n-th complex amplitude (2n...
Definition: QuEST.c:1308
void applyPauliSum(Qureg inQureg, enum pauliOpType *allPauliCodes, qreal *termCoeffs, int numSumTerms, Qureg outQureg)
Modifies outQureg to be the result of applying the weighted sum of Pauli products (a Hermitian but no...
Definition: QuEST.c:808
QMatrix getRandomQMatrix(int dim)
Returns a dim-by-dim complex matrix, where the real and imaginary value of each element are independe...
Definition: utilities.cpp:368
void applyMatrixN(Qureg qureg, int *targs, int numTargs, ComplexMatrixN u)
Apply a general N-by-N matrix, which may be non-unitary, on any number of target qubits.
Definition: QuEST.c:863
Represents a system of qubits.
Definition: QuEST.h:203
qreal ** imag
Definition: QuEST.h:140
long long int numElemsPerChunk
The number of the 2^numQubits amplitudes stored on each distributed node.
Definition: QuEST.h:183
Qureg createCloneQureg(Qureg qureg, QuESTEnv env)
Create a new Qureg which is an exact clone of the passed qureg, which can be either a statevector or ...
Definition: QuEST.c:64
std::vector< std::vector< qcomp > > QMatrix
A complex square matrix.
Definition: utilities.hpp:49
ComplexMatrix2 toComplexMatrix2(QMatrix qm)
Returns a ComplexMatrix2 copy of QMatix qm.
Definition: utilities.cpp:846
Catch::Generators::GeneratorWrapper< int * > sublists(int *list, int len, int sublen)
Returns a Catch2 generator of every length-sublen sublist of length-len list, in increasing lexograph...
Definition: utilities.cpp:1199
qreal * real
The real values of the 2^numQubits complex elements.
Definition: QuEST.h:189
Qureg createQureg(int numQubits, QuESTEnv env)
Create a Qureg object representing a set of qubits which will remain in a pure state.
Definition: QuEST.c:36
void destroyPauliHamil(PauliHamil h)
Destroy a PauliHamil instance, created with either createPauliHamil() or createPauliHamilFromFile().
Definition: QuEST.c:1163
void applyMatrix2(Qureg qureg, int targetQubit, ComplexMatrix2 u)
Apply a general 2-by-2 matrix, which may be non-unitary.
Definition: QuEST.c:844
#define PREPARE_TEST(quregVec, quregMatr, refVec, refMatr)
Prepares the needed data structures for unit testing some operators.
void applyReferenceMatrix(QVector &state, int *ctrls, int numCtrls, int *targs, int numTargs, QMatrix op)
Modifies the state-vector state to be the result of left-multiplying the multi-target operator matrix...
Definition: utilities.cpp:686
bool areEqual(QVector a, QVector b)
Returns true if the absolute value of the difference between every amplitude in vectors a and b is le...
Definition: utilities.cpp:387
Qureg createDensityQureg(int numQubits, QuESTEnv env)
Create a Qureg for qubits which are represented by a density matrix, and can be in mixed states.
Definition: QuEST.c:50
void applyTrotterCircuit(Qureg qureg, PauliHamil hamil, qreal time, int order, int reps)
Applies a trotterisation of unitary evolution to qureg.
Definition: QuEST.c:830
PauliHamil createPauliHamil(int numQubits, int numSumTerms)
Create a PauliHamil instance, which is a Hamiltonian expressed as a real-weighted sum of products of ...
Definition: QuEST.c:1147
DiagonalOp createDiagonalOp(int numQubits, QuESTEnv env)
Creates a DiagonalOp representing a diagonal operator on the full Hilbert space of a Qureg.
Definition: QuEST.c:1267
Represents a 2x2 matrix of complex numbers.
Definition: QuEST.h:114