test_data_structures.cpp
Go to the documentation of this file.
1 
2 #include "catch.hpp"
3 #include "QuEST.h"
4 #include "utilities.hpp"
5 
6 /* allows concise use of Contains in catch's REQUIRE_THROWS_WITH */
7 using Catch::Matchers::Contains;
8 
9 
10 
15 TEST_CASE( "fromComplex", "[data_structures]" ) {
16 
17  Complex a = {.real=.5, .imag=-.2};
18  qcomp b = fromComplex(a);
19 
20  REQUIRE( a.real == real(b) );
21  REQUIRE( a.imag == imag(b) );
22 }
23 
24 
25 
30 TEST_CASE( "getStaticComplexMatrixN", "[data_structures]" ) {
31 
32  /* use of this function is illegal in C++ */
33  SUCCEED( );
34 }
35 
36 
37 
42 TEST_CASE( "toComplex", "[data_structures]" ) {
43 
44  qcomp a = qcomp(.5,-.2);
45  Complex b = toComplex(a);
46 
47  REQUIRE( real(a) == b.real );
48  REQUIRE( imag(a) == b.imag );
49 }
50 
51 
52 
57 TEST_CASE( "createCloneQureg", "[data_structures]" ) {
58 
59  SECTION( "state-vector" ) {
60 
63 
64  // check properties are the same
65  REQUIRE( b.isDensityMatrix == a.isDensityMatrix );
67  REQUIRE( b.numQubitsInStateVec == a.numQubitsInStateVec );
68  REQUIRE( b.numAmpsPerChunk == a.numAmpsPerChunk );
69  REQUIRE( b.numAmpsTotal == a.numAmpsTotal );
70 
71  // check state-vector is the same (works for GPU and distributed)
72  REQUIRE( areEqual(a, b) );
73 
76  }
77  SECTION( "density-matrix" ) {
78 
81 
82  // check properties are the same
83  REQUIRE( b.isDensityMatrix == a.isDensityMatrix );
85  REQUIRE( b.numQubitsInStateVec == a.numQubitsInStateVec );
86  REQUIRE( b.numAmpsPerChunk == a.numAmpsPerChunk );
87  REQUIRE( b.numAmpsTotal == a.numAmpsTotal );
88 
89  // check state-vector is the same (works for GPU and distributed)
90  REQUIRE( areEqual(a, b) );
91 
94  }
95 }
96 
97 
98 
103 TEST_CASE( "createComplexMatrixN", "[data_structures]" ) {
104 
105  SECTION( "correctness" ) {
106 
107  int numQb = GENERATE( range(1,10+1) );
109 
110  // ensure elems are created and initialised to 0
111  REQUIRE( areEqual(toQMatrix(m), getZeroMatrix(1<<numQb)) );
112 
114  }
115  SECTION( "input validation" ) {
116 
117  SECTION( "number of qubits" ) {
118 
119  int numQb = GENERATE( -1, 0 );
120  REQUIRE_THROWS_WITH( createComplexMatrixN(numQb), Contains("Invalid number of qubits") );
121  }
122  }
123 }
124 
125 
126 
131 TEST_CASE( "createDensityQureg", "[data_structures]" ) {
132 
133  // must be at least one amplitude per node
134  int minNumQb = calcLog2(QUEST_ENV.numRanks) - 1; // density matrix has 2*numQb in state-vec
135  if (minNumQb <= 0)
136  minNumQb = 1;
137 
138  SECTION( "correctness" ) {
139 
140  // try 10 valid number of qubits
141  int numQb = GENERATE_COPY( range(minNumQb, minNumQb+10) );
142  Qureg reg = createDensityQureg(numQb, QUEST_ENV);
143 
144  // ensure elems (CPU and/or GPU) are created, and reg begins in |0><0|
145  QMatrix ref = getZeroMatrix(1<<numQb);
146  ref[0][0] = 1; // |0><0|
147  REQUIRE( areEqual(reg, ref) );
148 
149  destroyQureg(reg, QUEST_ENV);
150  }
151  SECTION( "input validation") {
152 
153  SECTION( "number of qubits" ) {
154 
155  int numQb = GENERATE( -1, 0 );
156  REQUIRE_THROWS_WITH( createDensityQureg(numQb, QUEST_ENV), Contains("Invalid number of qubits") );
157  }
158  SECTION( "number of amplitudes" ) {
159 
160  // use local QuESTEnv to safely modify
161  QuESTEnv env = QUEST_ENV;
162 
163  // too many amplitudes to store in type
164  int maxQb = (int) calcLog2(SIZE_MAX) / 2;
165  REQUIRE_THROWS_WITH( createDensityQureg(maxQb+1, env), Contains("Too many qubits") && Contains("size_t type") );
166 
167  /* n-qubit density matrix contains 2^(2n) amplitudes
168  * so can be spread between at most 2^(2n) ranks
169  */
170  int minQb = GENERATE_COPY( range(3,maxQb) );
171  env.numRanks = (int) pow(2, 2*minQb);
172  int numQb = GENERATE_COPY( range(1,minQb) );
173  REQUIRE_THROWS_WITH( createDensityQureg(numQb, env), Contains("Too few qubits") );
174  }
175  SECTION( "available memory" ) {
176 
177  /* there is no reliable way to force the malloc statements to
178  * fail, and hence trigger the matrixInit validation */
179  SUCCEED( );
180  }
181  }
182 }
183 
184 
185 
190 TEST_CASE( "createDiagonalOp", "[data_structures]" ) {
191 
192  // must be at least one amplitude per node
193  int minNumQb = calcLog2(QUEST_ENV.numRanks);
194  if (minNumQb == 0)
195  minNumQb = 1;
196 
197  SECTION( "correctness" ) {
198 
199  // try 10 valid number of qubits
200  int numQb = GENERATE_COPY( range(minNumQb, minNumQb+10) );
202 
203  // check properties are correct
204  REQUIRE( op.numQubits == numQb );
205  REQUIRE( op.chunkId == QUEST_ENV.rank );
206  REQUIRE( op.numChunks == QUEST_ENV.numRanks );
207  REQUIRE( op.numElemsPerChunk == (1LL << numQb) / QUEST_ENV.numRanks );
208  REQUIRE( op.real != NULL );
209  REQUIRE( op.imag != NULL );
210 
211  // check all elements in CPU are zero
212  REQUIRE( areEqual(toQVector(op), QVector(1LL << numQb)) );
213 
214  // (no concise way to check this for GPU)
215 
217  }
218  SECTION( "input validation" ) {
219 
220  SECTION( "number of qubits" ) {
221 
222  int numQb = GENERATE( -1, 0 );
223  REQUIRE_THROWS_WITH( createDiagonalOp(numQb, QUEST_ENV), Contains("Invalid number of qubits") );
224  }
225  SECTION( "number of elements" ) {
226 
227  // use local QuESTEnv to safely modify
228  QuESTEnv env = QUEST_ENV;
229 
230  // too many amplitudes to store in type
231  int maxQb = (int) calcLog2(SIZE_MAX);
232  REQUIRE_THROWS_WITH( createDiagonalOp(maxQb+1, env), Contains("Too many qubits") && Contains("size_t type") );
233 
234  // too few amplitudes to distribute
235  int minQb = GENERATE_COPY( range(2,maxQb) );
236  env.numRanks = (int) pow(2, minQb);
237  int numQb = GENERATE_COPY( range(1,minQb) );
238  REQUIRE_THROWS_WITH( createDiagonalOp(numQb, env), Contains("Too few qubits") && Contains("distributed"));
239  }
240  SECTION( "available memory" ) {
241 
242  /* there is no reliable way to force the malloc statements to
243  * fail, and hence trigger the diagonalOpInit validation */
244  SUCCEED( );
245  }
246  }
247 }
248 
249 
250 
255 TEST_CASE( "createDiagonalOpFromPauliHamilFile", "[data_structures]" ) {
256 
257  // files created & populated during the test, and deleted afterward
258  char fnPrefix[] = "temp_createDiagonalOpFromPauliHamilFile";
259  char fn[100];
260 
261  // each test uses a unique filename (managed by master node), to avoid file IO locks
262  // (it's safe for sub-test to overwrite this, since after each test, all files
263  // with prefix fnPrefix are deleted)
264  setUniqueFilename(fn, fnPrefix);
265 
266  // diagonal op must have at least one amplitude per node
267  int minNumQb = calcLog2(QUEST_ENV.numRanks);
268  if (minNumQb == 0)
269  minNumQb = 1;
270 
271  SECTION( "correctness" ) {
272 
273  SECTION( "general" ) {
274 
275  // try several Pauli Hamiltonian sizes
276  int numQb = GENERATE_COPY( range(minNumQb, 6+minNumQb) );
277  int numTerms = GENERATE_COPY( 1, minNumQb, 10*minNumQb );
278 
279  // create a PauliHamil with random elements
280  PauliHamil hamil = createPauliHamil(numQb, numTerms);
282 
283  // write the Hamiltonian to file (with trailing whitespace, and trailing newline)
284  if (QUEST_ENV.rank == 0) {
285  FILE* file = fopen(fn, "w");
286  int i=0;
287  for (int n=0; n<numTerms; n++) {
288  fprintf(file, REAL_STRING_FORMAT, hamil.termCoeffs[n]);
289  fprintf(file, " ");
290  for (int q=0; q<numQb; q++)
291  fprintf(file, "%d ", (int) hamil.pauliCodes[i++]);
292  fprintf(file, "\n");
293  }
294  fprintf(file, "\n");
295  fclose(file);
296  }
298 
299  // load the file as a diagonal operator, and compare
301  REQUIRE( areEqual(toQMatrix(op), toQMatrix(hamil)) );
302 
303  destroyPauliHamil(hamil);
305  }
306  SECTION( "edge cases" ) {
307 
308  // prepare a valid single-term diagonal Pauli Hamiltonian
309  qreal coeffs[] = {.1};
310  pauliOpType codes[minNumQb];
311  for (int q=0; q<minNumQb; q++)
312  codes[q] = (q%2)? PAULI_I : PAULI_Z;
313 
314  QMatrix ref = toQMatrix(coeffs, codes, minNumQb, 1);
315 
316  // prepare basic encoding string
317  string line = to_string(coeffs[0]) + " ";
318  for (int q=0; q<minNumQb; q++)
319  line += to_string(codes[q]) + ((q<minNumQb-1)? " ":"");
320 
321  SECTION( "no trailing newline or space" ) {
322 
323  writeToFileSynch(fn, line);
324 
326  REQUIRE( areEqual(ref, toQMatrix(op)) );
327 
329  }
330  SECTION( "trailing newlines" ) {
331 
332  writeToFileSynch(fn, line + "\n\n\n");
333 
335  REQUIRE( areEqual(ref, toQMatrix(op)) );
336 
338  }
339  SECTION( "trailing spaces" ) {
340 
341  writeToFileSynch(fn, line + " ");
342 
344  REQUIRE( areEqual(ref, toQMatrix(op)) );
345 
347  }
348  }
349  }
350  SECTION( "input validation") {
351 
352  SECTION( "number of qubits" ) {
353 
354  writeToFileSynch(fn, ".1 "); // 0 qubits
355  REQUIRE_THROWS_WITH( createDiagonalOpFromPauliHamilFile(fn, QUEST_ENV), Contains("The number of qubits") && Contains("strictly positive"));
356  }
357  SECTION( "number of elements" ) {
358 
359  // too many amplitudes to store in type
360  int maxQb = (int) calcLog2(SIZE_MAX);
361 
362  // encode one more qubit than legal to file
363  string line = ".1 ";
364  for (int q=0; q<(maxQb+1); q++)
365  line += "3 "; // trailing space ok
366  writeToFileSynch(fn, line);
367 
368  REQUIRE_THROWS_WITH( createDiagonalOpFromPauliHamilFile(fn, QUEST_ENV), Contains("Too many qubits") && Contains("size_t type") );
369 
370  // use local QuESTEnv to safely modify
371  QuESTEnv env = QUEST_ENV;
372 
373  // too few elements to distribute
374  int minQb = GENERATE_COPY( range(2,maxQb) );
375  env.numRanks = (int) pow(2, minQb);
376  int numQb = GENERATE_COPY( range(1,minQb) );
377 
378  line = ".1 ";
379  for (int q=0; q<numQb; q++)
380  line += "3 "; // trailing space ok
381  setUniqueFilename(fn, fnPrefix);
382  writeToFileSynch(fn, line);
383 
384  REQUIRE_THROWS_WITH( createDiagonalOpFromPauliHamilFile(fn, env), Contains("Too few qubits") && Contains("distributed") );
385  }
386  SECTION( "coefficient type" ) {
387 
388  writeToFileSynch(fn, "notanumber 1 2 3");
389  REQUIRE_THROWS_WITH( createDiagonalOpFromPauliHamilFile(fn, QUEST_ENV), Contains("Failed to parse") && Contains("coefficient"));
390  }
391  SECTION( "pauli code" ) {
392 
393  writeToFileSynch(fn, ".1 0 3 2"); // final is invalid Y
394  REQUIRE_THROWS_WITH( createDiagonalOpFromPauliHamilFile(fn, QUEST_ENV), Contains("contained operators other than PAULI_Z and PAULI_I"));
395 
396  setUniqueFilename(fn, fnPrefix);
397  writeToFileSynch(fn, ".1 0 1 3"); // second is invalid X
398  REQUIRE_THROWS_WITH( createDiagonalOpFromPauliHamilFile(fn, QUEST_ENV), Contains("contained operators other than PAULI_Z and PAULI_I"));
399 
400  setUniqueFilename(fn, fnPrefix);
401  writeToFileSynch(fn, ".1 0 1 4"); // final is invalid Pauli code
402  REQUIRE_THROWS_WITH( createDiagonalOpFromPauliHamilFile(fn, QUEST_ENV), Contains("invalid pauli code"));
403 
404  setUniqueFilename(fn, fnPrefix);
405  writeToFileSynch(fn, ".1 3 0 notanumber"); // final is invalid type
406  REQUIRE_THROWS_WITH( createDiagonalOpFromPauliHamilFile(fn, QUEST_ENV), Contains("Failed to parse the next expected Pauli code"));
407  }
408  }
409 
410  // delete all files created above
411  deleteFilesWithPrefixSynch(fnPrefix);
412 }
413 
414 
415 
420 TEST_CASE( "createPauliHamil", "[data_structures]" ) {
421 
422  SECTION( "correctness" ) {
423 
424  int numQb = GENERATE( range(1,5) );
425  int numTerms = GENERATE( range(1,5) );
426  PauliHamil hamil = createPauliHamil(numQb, numTerms);
427 
428  // check fields are correct
429  REQUIRE( hamil.numQubits == numQb );
430  REQUIRE( hamil.numSumTerms == numTerms );
431 
432  // check all Pauli codes are identity
433  int numPaulis = numQb * numTerms;
434  for (int i=0; i<numPaulis; i++) {
435  REQUIRE( hamil.pauliCodes[i] == PAULI_I );
436  }
437 
438  // check all term coefficients can be written to (no seg fault)
439  for (int j=0; j<numTerms; j++) {
440  hamil.termCoeffs[j] = 1;
441  REQUIRE( hamil.termCoeffs[j] == 1 );
442  }
443 
444  destroyPauliHamil(hamil);
445  }
446  SECTION( "input validation") {
447 
448  SECTION( "number of qubits" ) {
449 
450  int numQb = GENERATE( -1, 0 );
451  REQUIRE_THROWS_WITH( createPauliHamil(numQb, 1), Contains("The number of qubits and terms in the PauliHamil must be strictly positive.") );
452  }
453  SECTION( "number of terms" ) {
454 
455  int numTerms = GENERATE( -1, 0 );
456  REQUIRE_THROWS_WITH( createPauliHamil(1, numTerms), Contains("The number of qubits and terms in the PauliHamil must be strictly positive.") );
457  }
458  }
459 }
460 
461 
462 
467 TEST_CASE( "createPauliHamilFromFile", "[data_structures]" ) {
468 
469  // files created & populated during the test, and deleted afterward
470  char fnPrefix[] = "temp_createPauliHamilFromFile";
471  char fn[100];
472 
473  // each test uses a unique filename (managed by master node), to avoid file IO locks
474  // (it's safe for sub-test to overwrite this, since after each test, all files
475  // with prefix fnPrefix are deleted)
476  setUniqueFilename(fn, fnPrefix);
477 
478  SECTION( "correctness" ) {
479 
480  SECTION( "general" ) {
481 
482  // for several sizes...
483  int numQb = GENERATE( 1, 5, 10, 15 );
484  int numTerms = GENERATE( 1, 10, 30 );
485  int numPaulis = numQb*numTerms;
486 
487  // create a PauliHamil with random elements
488  qreal coeffs[numTerms];
489  enum pauliOpType paulis[numPaulis];
490  setRandomPauliSum(coeffs, paulis, numQb, numTerms);
491 
492  // write the Hamiltonian to file (with trailing whitespace, and trailing newline)
493  if (QUEST_ENV.rank == 0) {
494  FILE* file = fopen(fn, "w");
495  int i=0;
496  for (int n=0; n<numTerms; n++) {
497  fprintf(file, REAL_STRING_FORMAT, coeffs[n]);
498  fprintf(file, " ");
499  for (int q=0; q<numQb; q++)
500  fprintf(file, "%d ", (int) paulis[i++]);
501  fprintf(file, "\n");
502  }
503  fprintf(file, "\n");
504  fclose(file);
505  }
507 
508  // load the file as a PauliHamil
510 
511  // check fields agree
512  REQUIRE( hamil.numQubits == numQb );
513  REQUIRE( hamil.numSumTerms == numTerms );
514 
515  // check elements agree
516  int j=0;
517  for (int n=0; n<numTerms; n++) {
518  REQUIRE( absReal(hamil.termCoeffs[n] - coeffs[n]) <= REAL_EPS );
519  for (int q=0; q<numQb; q++) {
520  REQUIRE( hamil.pauliCodes[j] == paulis[j] );
521  j++;
522  }
523  }
524 
525  destroyPauliHamil(hamil);
526  }
527  SECTION( "edge cases" ) {
528 
529  SECTION( "no trailing newline or space" ) {
530 
531  writeToFileSynch(fn, ".1 1 0 1");
532 
534  REQUIRE( hamil.numSumTerms == 1 );
535  REQUIRE( hamil.numQubits == 3 );
536 
537  destroyPauliHamil(hamil);
538  }
539  SECTION( "trailing newlines" ) {
540 
541  writeToFileSynch(fn, ".1 1 0 1\n\n\n");
542 
544  REQUIRE( hamil.numSumTerms == 1 );
545  REQUIRE( hamil.numQubits == 3 );
546 
547  destroyPauliHamil(hamil);
548  }
549  SECTION( "trailing spaces" ) {
550 
551  writeToFileSynch(fn, ".1 1 0 1 ");
552 
554  REQUIRE( hamil.numSumTerms == 1 );
555  REQUIRE( hamil.numQubits == 3 );
556 
557  destroyPauliHamil(hamil);
558  }
559  }
560  }
561  SECTION( "input validation") {
562 
563  SECTION( "number of qubits" ) {
564 
565  writeToFileSynch(fn, ".1 ");
566 
567  REQUIRE_THROWS_WITH( createPauliHamilFromFile(fn), Contains("The number of qubits") && Contains("strictly positive"));
568  }
569  SECTION( "coefficient type" ) {
570 
571  writeToFileSynch(fn, "notanumber 1 2 3");
572 
573  REQUIRE_THROWS_WITH( createPauliHamilFromFile(fn), Contains("Failed to parse") && Contains("coefficient"));
574  }
575  SECTION( "pauli code" ) {
576 
577  writeToFileSynch(fn, ".1 1 2 4"); // invalid int
578 
579  REQUIRE_THROWS_WITH( createPauliHamilFromFile(fn), Contains("invalid pauli code"));
580 
581  setUniqueFilename(fn, fnPrefix);
582  writeToFileSynch(fn, ".1 1 2 notanumber"); // invalid type
583 
584  REQUIRE_THROWS_WITH( createPauliHamilFromFile(fn), Contains("Failed to parse the next expected Pauli code"));
585  }
586  }
587 
588  // cleanup temp files
590 }
591 
592 
593 
598 TEST_CASE( "createQuESTEnv", "[data_structures]" ) {
599 
600  /* there is no meaningful way to test this */
601  SUCCEED( );
602 }
603 
604 
605 
610 TEST_CASE( "createQureg", "[data_structures]" ) {
611 
612  // must be at least one amplitude per node
613  int minNumQb = calcLog2(QUEST_ENV.numRanks);
614  if (minNumQb == 0)
615  minNumQb = 1;
616 
617  SECTION( "correctness" ) {
618 
619  // try 10 valid number of qubits
620  int numQb = GENERATE_COPY( range(minNumQb, minNumQb+10) );
621  Qureg reg = createQureg(numQb, QUEST_ENV);
622 
623  // ensure elems (CPU and/or GPU) are created, and reg begins in |0>
624  QVector ref = QVector(1<<numQb);
625  ref[0] = 1; // |0>
626  REQUIRE( areEqual(reg, ref) );
627 
628  destroyQureg(reg, QUEST_ENV);
629  }
630  SECTION( "input validation") {
631 
632  SECTION( "number of qubits" ) {
633 
634  int numQb = GENERATE( -1, 0 );
635  REQUIRE_THROWS_WITH( createQureg(numQb, QUEST_ENV), Contains("Invalid number of qubits") );
636  }
637  SECTION( "number of amplitudes" ) {
638 
639  // use local QuESTEnv to safely modify
640  QuESTEnv env = QUEST_ENV;
641 
642  // too many amplitudes to store in type
643  int maxQb = (int) calcLog2(SIZE_MAX);
644  REQUIRE_THROWS_WITH( createQureg(maxQb+1, env), Contains("Too many qubits") && Contains("size_t type") );
645 
646  // too few amplitudes to distribute
647  int minQb = GENERATE_COPY( range(2,maxQb) );
648  env.numRanks = (int) pow(2, minQb);
649  int numQb = GENERATE_COPY( range(1,minQb) );
650  REQUIRE_THROWS_WITH( createQureg(numQb, env), Contains("Too few qubits") );
651  }
652  SECTION( "available memory" ) {
653 
654  /* there is no reliable way to force the malloc statements to
655  * fail, and hence trigger the matrixInit validation */
656  SUCCEED( );
657  }
658  }
659 }
660 
661 
662 
667 TEST_CASE( "destroyComplexMatrixN", "[data_structures]" ) {
668 
669  SECTION( "correctness" ) {
670 
671  /* there is no meaningful way to test this */
672  SUCCEED( );
673  }
674  SECTION( "input validation" ) {
675 
676  SECTION( "matrix not created" ) {
677 
678  /* this is an artificial test case since nothing in the QuEST API
679  * automatically sets un-initialised ComplexMatrixN fields to
680  * the NULL pointer.
681  */
682  ComplexMatrixN m;
683  m.real = NULL;
684 
685  /* the error message is also somewhat unrelated, but oh well
686  */
687  REQUIRE_THROWS_WITH( destroyComplexMatrixN(m), Contains("The ComplexMatrixN was not successfully created") );
688  }
689  }
690 }
691 
692 
693 
698 TEST_CASE( "destroyDiagonalOp", "[data_structures]" ) {
699 
700  /* there is no meaningful way to test this */
701  SUCCEED( );
702 }
703 
704 
705 
710 TEST_CASE( "destroyPauliHamil", "[data_structures]" ) {
711 
712  /* there is no meaningful way to test this.
713  * We e.g. cannot check that the pointers are NULL because
714  * they are not updated; this function passes the struct by value,
715  * not by reference. We also cannot reliably monitor the
716  * memory used in the heap at runtime.
717  */
718  SUCCEED( );
719 }
720 
721 
722 
727 TEST_CASE( "destroyQuESTEnv", "[data_structures]" ) {
728 
729  /* there is no meaningful way to test this */
730  SUCCEED( );
731 }
732 
733 
734 
739 TEST_CASE( "destroyQureg", "[data_structures]" ) {
740 
741  /* there is no meaningful way to test this.
742  * We e.g. cannot check that the pointers are NULL because
743  * they are not updated; this function passes the struct by value,
744  * not by reference. We also cannot reliably monitor the
745  * memory used in the heap at runtime.
746  */
747  SUCCEED( );
748 }
749 
750 
751 
756 TEST_CASE( "initComplexMatrixN", "[data_structures]" ) {
757 
758  /* use of this function is illegal in C++ */
759  SUCCEED( );
760 }
761 
762 
763 
768 TEST_CASE( "initDiagonalOp", "[data_structures]" ) {
769 
770  // must be at least one amplitude per node
771  int minNumQb = calcLog2(QUEST_ENV.numRanks);
772  if (minNumQb == 0)
773  minNumQb = 1;
774 
775  SECTION( "correctness" ) {
776 
777  // try 10 valid number of qubits
778  int numQb = GENERATE_COPY( range(minNumQb, minNumQb+10) );
780 
781  long long int len = (1LL << numQb);
782  qreal reals[len];
783  qreal imags[len];
784  long long int n;
785  for (n=0; n<len; n++) {
786  reals[n] = (qreal) n;
787  imags[n] = (qreal) -2*n; // (n - 2n i)
788  }
789  initDiagonalOp(op, reals, imags);
790 
791  // check that op.real and op.imag were modified...
792  REQUIRE( areEqual(toQVector(op), reals, imags) );
793 
794  // and also that GPU real and imag were modified
795  // via if it modifies an all-unity state-vec correctly
796  Qureg qureg = createQureg(numQb, QUEST_ENV);
797  for (long long int i=0; i<qureg.numAmpsPerChunk; i++) {
798  qureg.stateVec.real[i] = 1;
799  qureg.stateVec.imag[i] = 1;
800  }
801  copyStateToGPU(qureg);
802 
803  QVector prodRef = toQMatrix(op) * toQVector(qureg);
804 
805  // (n - 2n i) * (1 + 1i) = 3n - n*i
806  applyDiagonalOp(qureg, op);
807  copyStateFromGPU(qureg);
808  QVector result = toQVector(qureg);
809  REQUIRE( areEqual(prodRef, result) );
810 
811  destroyQureg(qureg, QUEST_ENV);
813  }
814 }
815 
816 
817 
822 TEST_CASE( "initDiagonalOpFromPauliHamil", "[data_structures]" ) {
823 
824  // distributed diagonal op must contain at least one amplitude per node
825  int minNumQb = calcLog2(QUEST_ENV.numRanks);
826  if (minNumQb == 0)
827  minNumQb = 1;
828 
829  SECTION( "correctness" ) {
830 
831  // try (at most) 10 valid number of qubits (even for validation)
832  int numQb = GENERATE_COPY( range(minNumQb, min(10,minNumQb+10)) );
834 
835  // try several sized random all-Z Hamiltonians
836  int numTerms = GENERATE_COPY( 1, numQb, 5*numQb );
837  PauliHamil hamil = createPauliHamil(numQb, numTerms);
839 
840  initDiagonalOpFromPauliHamil(op, hamil);
841  REQUIRE( areEqual(toQMatrix(op), toQMatrix(hamil)) );
842 
843  destroyPauliHamil(hamil);
845  }
846  SECTION( "input validation" ) {
847 
848  SECTION( "hamiltonian parameters" ) {
849 
850  DiagonalOp op = createDiagonalOp(minNumQb, QUEST_ENV);
851  PauliHamil hamil;
852 
853  hamil.numQubits = GENERATE( -1, 0 );
854  hamil.numSumTerms = 1;
855  REQUIRE_THROWS_WITH( initDiagonalOpFromPauliHamil(op, hamil), Contains("number of qubits") && Contains("strictly positive") );
856 
857  hamil.numQubits = minNumQb;
858  hamil.numSumTerms = GENERATE( -1, 0 );
859  REQUIRE_THROWS_WITH( initDiagonalOpFromPauliHamil(op, hamil), Contains("terms") && Contains("strictly positive") );
860 
862  }
863  SECTION( "mismatching dimensions" ) {
864 
865  int numQbA = minNumQb+1;
866  int numQbB = GENERATE_COPY( numQbA-1, numQbA+1 );
867 
868  DiagonalOp op = createDiagonalOp(numQbA, QUEST_ENV);
869  PauliHamil hamil = createPauliHamil(numQbB, 1);
870 
871  REQUIRE_THROWS_WITH( initDiagonalOpFromPauliHamil(op, hamil), Contains("Pauli Hamiltonian and diagonal operator have different, incompatible dimensions") );
872 
874  destroyPauliHamil(hamil);
875  }
876  SECTION( "pauli codes" ) {
877 
878  DiagonalOp op = createDiagonalOp(minNumQb, QUEST_ENV);
879  PauliHamil hamil = createPauliHamil(minNumQb, 5); // all I
880 
881  // make only one code invalid
882  int numCodes = minNumQb * hamil.numSumTerms;
883  int ind = GENERATE_COPY( range(0,numCodes) );
884  hamil.pauliCodes[ind] = GENERATE( PAULI_X, PAULI_Y );
885 
886  REQUIRE_THROWS_WITH( initDiagonalOpFromPauliHamil(op, hamil), Contains("contained operators other than PAULI_Z and PAULI_I") );
887 
889  destroyPauliHamil(hamil);
890  }
891  }
892 }
893 
894 
895 
900 TEST_CASE( "initPauliHamil", "[data_structures]" ) {
901 
902  SECTION( "correctness" ) {
903 
904  PauliHamil hamil = createPauliHamil(3, 2);
905 
906  qreal coeffs[] = {-5, 5};
907  enum pauliOpType codes[] = {
910  initPauliHamil(hamil, coeffs, codes);
911 
912  // check everything written correctly
913  for (int t=0; t<2; t++) {
914  REQUIRE( coeffs[t] == hamil.termCoeffs[t] );
915  for (int q=0; q<3; q++) {
916  int ind = 3*t+q;
917  REQUIRE( codes[ind] == hamil.pauliCodes[ind] );
918  }
919  }
920 
921  destroyPauliHamil(hamil);
922  }
923  SECTION( "input validation" ) {
924 
925  SECTION( "parameters" ) {
926 
927  // parameters checked before codes, so safe to leave un-initialised
928  qreal coeffs[1];
929  enum pauliOpType codes[1];
930  PauliHamil hamil;
931 
932  hamil.numQubits = GENERATE( -1, 0 );
933  hamil.numSumTerms = 1;
934  REQUIRE_THROWS_WITH( initPauliHamil(hamil, coeffs, codes), Contains("number of qubits") && Contains("strictly positive") );
935 
936  hamil.numQubits = 1;
937  hamil.numSumTerms = GENERATE( -1, 0 );
938  REQUIRE_THROWS_WITH( initPauliHamil(hamil, coeffs, codes), Contains("terms") && Contains("strictly positive") );
939  }
940  SECTION( "Pauli codes" ) {
941 
942  int numQb = 3;
943  int numTerms = 2;
944  int numCodes = numQb * numTerms;
945  qreal coeffs[numTerms];
946  enum pauliOpType codes[numCodes];
947 
948  // make only one code invalid
949  for (int i=0; i<numCodes; i++)
950  codes[i] = PAULI_I;
951  codes[GENERATE_COPY( range(0,numCodes) )] = (pauliOpType) GENERATE( -1, 4 );
952 
953  PauliHamil hamil = createPauliHamil(numQb, numTerms);
954  REQUIRE_THROWS_WITH( initPauliHamil(hamil, coeffs, codes), Contains("Invalid Pauli code") );
955  destroyPauliHamil(hamil);
956  }
957  }
958 }
959 
960 
961 
966 TEST_CASE( "setDiagonalOpElems", "[data_structures]" ) {
967 
968  // must be at least one amplitude per node
969  int minNumQb = calcLog2(QUEST_ENV.numRanks);
970  if (minNumQb == 0)
971  minNumQb = 1;
972 
973  // try 10 valid number of qubits (even for validation)
974  int numQb = GENERATE_COPY( range(minNumQb, minNumQb+10) );
976 
977  SECTION( "correctness" ) {
978 
979  // make entire array on every node
980  long long int len = (1LL << numQb);
981  qreal reals[len];
982  qreal imags[len];
983  long long int n;
984  for (n=0; n<len; n++) {
985  reals[n] = (qreal) n;
986  imags[n] = (qreal) -2*n; // (n - 2n i)
987  }
988 
989  // set one value at a time (only relevant nodes will update)
990  for (n=0; n<len; n++)
991  setDiagonalOpElems(op, n, &reals[n], &imags[n], 1);
992 
993  // check op.real and op.imag updated correctly
994  REQUIRE( areEqual(toQVector(op), reals, imags) );
995 
996  // no check that GPU values updated (occurs in initDiagonalOp)
997  }
998  SECTION( "input validation" ) {
999 
1000  long long int maxInd = (1LL << numQb);
1001  qreal *reals;
1002  qreal *imags;
1003 
1004  SECTION( "start index" ) {
1005 
1006  int startInd = GENERATE_COPY( -1, maxInd );
1007  int numAmps = 1;
1008  REQUIRE_THROWS_WITH( setDiagonalOpElems(op, startInd, reals, imags, numAmps), Contains("Invalid element index") );
1009  }
1010 
1011  SECTION( "number of elements" ) {
1012 
1013  // independent
1014  int startInd = 0;
1015  int numAmps = GENERATE_COPY( -1, maxInd+1 );
1016  REQUIRE_THROWS_WITH( setDiagonalOpElems(op, startInd, reals, imags, numAmps), Contains("Invalid number of elements") );
1017 
1018  // invalid considering start-index
1019  startInd = maxInd - 1;
1020  numAmps = 2;
1021  REQUIRE_THROWS_WITH( setDiagonalOpElems(op, startInd, reals, imags, numAmps), Contains("More elements given than exist") );
1022  }
1023  }
1024 
1026 }
1027 
1028 
1029 
1034 TEST_CASE( "syncDiagonalOp", "[data_structures]" ) {
1035 
1036  // must be at least one amplitude per node
1037  int minNumQb = calcLog2(QUEST_ENV.numRanks);
1038  if (minNumQb == 0)
1039  minNumQb = 1;
1040 
1041  SECTION( "correctness" ) {
1042 
1043  // try 10 valid number of qubits
1044  int numQb = GENERATE_COPY( range(minNumQb, minNumQb+10) );
1045  DiagonalOp op = createDiagonalOp(numQb, QUEST_ENV);
1046 
1047  // check that changes get sync'd to the GPU...
1048  long long int n;
1049  for (n=0; n<op.numElemsPerChunk; n++) {
1050  op.real[n] = (qreal) n;
1051  op.imag[n] = (qreal) -2*n; // (n - 2n i)
1052  }
1053  syncDiagonalOp(op);
1054 
1055  // via if it modifies an all-unity state-vec correctly
1056  Qureg qureg = createQureg(numQb, QUEST_ENV);
1057  for (long long int i=0; i<qureg.numAmpsPerChunk; i++) {
1058  qureg.stateVec.real[i] = 1;
1059  qureg.stateVec.imag[i] = 1;
1060  }
1061  copyStateToGPU(qureg);
1062 
1063  // (n - 2n i) * (1 + 1i) = 3n - n*i
1064  applyDiagonalOp(qureg, op);
1065  copyStateFromGPU(qureg);
1066  for (n=0; n<qureg.numAmpsPerChunk; n++) {
1067  REQUIRE( qureg.stateVec.real[n] == 3*n );
1068  REQUIRE( qureg.stateVec.imag[n] == -n );
1069  }
1070 
1071  destroyQureg(qureg, QUEST_ENV);
1073  }
1074 }
pauliOpType
Codes for specifying Pauli operators.
Definition: QuEST.h:96
void copyStateFromGPU(Qureg qureg)
In GPU mode, this copies the state-vector (or density matrix) from GPU memory (qureg....
Definition: QuEST_cpu.c:45
void initPauliHamil(PauliHamil hamil, qreal *coeffs, enum pauliOpType *codes)
Initialise PauliHamil instance hamil with the given term coefficients and Pauli codes (one for every ...
Definition: QuEST.c:1504
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
void syncQuESTEnv(QuESTEnv env)
Guarantees that all code up to the given point has been executed on all nodes (if running in distribu...
#define fromComplex(comp)
@ PAULI_Z
Definition: QuEST.h:96
int rank
Definition: QuEST.h:364
void destroyComplexMatrixN(ComplexMatrixN m)
Destroy a ComplexMatrixN instance created with createComplexMatrixN()
Definition: QuEST.c:1369
int numChunks
The number of nodes between which the elements of this operator are split.
Definition: QuEST.h:304
@ PAULI_I
Definition: QuEST.h:96
TEST_CASE("fromComplex", "[data_structures]")
DiagonalOp createDiagonalOpFromPauliHamilFile(char *fn, QuESTEnv env)
Creates and initialiases a diagonal operator from the Z Pauli Hamiltonian encoded in file with filena...
Definition: QuEST.c:1558
ComplexMatrixN createComplexMatrixN(int numQubits)
Allocate dynamic memory for a square complex matrix of any size, which can be passed to functions lik...
Definition: QuEST.c:1348
void syncDiagonalOp(DiagonalOp op)
Update the GPU memory with the current values in op.real and op.imag.
Definition: QuEST.c:1531
void destroyDiagonalOp(DiagonalOp op, QuESTEnv env)
Destroys a DiagonalOp created with createDiagonalOp(), freeing its memory.
Definition: QuEST.c:1524
#define NUM_QUBITS
The default number of qubits in the registers created for unit testing (both statevectors and density...
Definition: utilities.hpp:36
void writeToFileSynch(char *fn, const string &contents)
Writes contents to the file with filename fn, which is created and/or overwritten.
Definition: utilities.cpp:1362
int chunkId
The position of the chunk of the operator held by this process in the full operator.
Definition: QuEST.h:306
Information about the environment the program is running in.
Definition: QuEST.h:362
PauliHamil createPauliHamilFromFile(char *fn)
Creates a PauliHamil instance, a real-weighted sum of products of Pauli operators,...
Definition: QuEST.c:1420
void setDiagonalOpElems(DiagonalOp op, long long int startInd, qreal *real, qreal *imag, long long int numElems)
Modifies a subset (starting at index startInd, and ending at index startInd + numElems) of the elemen...
Definition: QuEST.c:1543
Represents a general 2^N by 2^N matrix of complex numbers.
Definition: QuEST.h:186
#define qreal
QMatrix toQMatrix(ComplexMatrix2 src)
Returns a copy of the given 2-by-2 matrix.
Definition: utilities.cpp:1044
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
int numQubitsInStateVec
Number of qubits in the state-vector - this is double the number represented for mixed states.
Definition: QuEST.h:329
void setUniqueFilename(char *outFn, char *prefix)
Modifies outFn to be a filename of format prefix_NUM.txt where NUM is a new unique integer so far.
Definition: utilities.cpp:1358
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:1229
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:310
QVector toQVector(Qureg qureg)
Returns an equal-size copy of the given state-vector qureg.
Definition: utilities.cpp:1113
long long int numAmpsPerChunk
Number of probability amplitudes held in stateVec by this process In the non-MPI version,...
Definition: QuEST.h:332
void applyDiagonalOp(Qureg qureg, DiagonalOp op)
Apply a diagonal operator, which is possibly non-unitary and non-Hermitian, to the entire qureg.
Definition: QuEST.c:1127
qreal * termCoeffs
The real coefficient of each Pauli product. This is an array of length PauliHamil....
Definition: QuEST.h:283
#define qcomp
enum pauliOpType * pauliCodes
The Pauli operators acting on each qubit, flattened over every operator.
Definition: QuEST.h:281
void copyStateToGPU(Qureg qureg)
In GPU mode, this copies the state-vector (or density matrix) from RAM (qureg.stateVec) to VRAM / GPU...
Definition: QuEST_cpu.c:42
#define toComplex(scalar)
int numRanks
Definition: QuEST.h:365
int numQubits
The number of qubits this operator can act on (informing its size)
Definition: QuEST.h:300
int numSumTerms
The number of terms in the weighted sum, or the number of Pauli products.
Definition: QuEST.h:285
Represents a diagonal complex operator on the full Hilbert state of a Qureg.
Definition: QuEST.h:297
@ PAULI_Y
Definition: QuEST.h:96
A Pauli Hamiltonian, expressed as a real-weighted sum of pauli products, and which can hence represen...
Definition: QuEST.h:277
void destroyQureg(Qureg qureg, QuESTEnv env)
Deallocate a Qureg, freeing its memory.
Definition: QuEST.c:77
void setRandomDiagPauliHamil(PauliHamil hamil)
Populates hamil with random coefficients and a random amount number of PAULI_I and PAULI_Z operators.
Definition: utilities.cpp:1241
qreal ** real
Definition: QuEST.h:189
Represents a system of qubits.
Definition: QuEST.h:322
void initDiagonalOp(DiagonalOp op, qreal *real, qreal *imag)
Overwrites the entire DiagonalOp op with the given real and imag complex elements.
Definition: QuEST.c:1537
ComplexArray stateVec
Computational state amplitudes - a subset thereof in the MPI version.
Definition: QuEST.h:341
long long int numElemsPerChunk
The number of the 2^numQubits amplitudes stored on each distributed node.
Definition: QuEST.h:302
void initDiagonalOpFromPauliHamil(DiagonalOp op, PauliHamil hamil)
Populates the diagonal operator op to be equivalent to the given Pauli Hamiltonian hamil,...
Definition: QuEST.c:1550
int isDensityMatrix
Whether this instance is a density-state representation.
Definition: QuEST.h:325
Qureg createCloneQureg(Qureg qureg, QuESTEnv env)
Create a new Qureg which is an exact clone of the passed qureg, which can be either a state-vector or...
Definition: QuEST.c:64
std::vector< std::vector< qcomp > > QMatrix
A complex square matrix.
Definition: utilities.hpp:49
int numQubits
The number of qubits informing the Hilbert dimension of the Hamiltonian.
Definition: QuEST.h:287
int numQubitsRepresented
The number of qubits represented in either the state-vector or density matrix.
Definition: QuEST.h:327
long long int numAmpsTotal
Total number of amplitudes, which are possibly distributed among machines.
Definition: QuEST.h:334
qreal * real
The real values of the 2^numQubits complex elements.
Definition: QuEST.h:308
qreal real
Definition: QuEST.h:105
Qureg createQureg(int numQubits, QuESTEnv env)
Creates a state-vector 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:1414
qreal imag
Definition: QuEST.h:106
void deleteFilesWithPrefixSynch(char *prefix)
Deletes all files with filename starting with prefix.
Definition: utilities.cpp:1375
QMatrix getZeroMatrix(size_t dim)
Returns a dim-by-dim square complex matrix, initialised to all zeroes.
Definition: utilities.cpp:153
Represents one complex number.
Definition: QuEST.h:103
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:398
Qureg createDensityQureg(int numQubits, QuESTEnv env)
Creates a density matrix Qureg object representing a set of qubits which can enter noisy and mixed st...
Definition: QuEST.c:50
PauliHamil createPauliHamil(int numQubits, int numSumTerms)
Dynamically allocates a Hamiltonian expressed as a real-weighted sum of products of Pauli operators.
Definition: QuEST.c:1398
DiagonalOp createDiagonalOp(int numQubits, QuESTEnv env)
Creates a DiagonalOp representing a diagonal operator on the full Hilbert space of a Qureg.
Definition: QuEST.c:1518