QuEST_common.c
Go to the documentation of this file.
1 // Distributed under MIT licence. See https://github.com/QuEST-Kit/QuEST/blob/master/LICENCE.txt for details
2 
18 # include "QuEST.h"
19 # include "QuEST_internal.h"
20 # include "QuEST_precision.h"
21 # include "QuEST_validation.h"
22 # include "QuEST_qasm.h"
23 # include "mt19937ar.h"
24 
25 #if defined(_WIN32) && ! defined(__MINGW32__)
26  #include <Windows.h>
27  #include <io.h>
28  #include <process.h>
29 #else
30  #include <unistd.h>
31  #include <sys/time.h>
32 #endif
33 
34 # include <sys/types.h>
35 # include <stdio.h>
36 # include <stdlib.h>
37 
38 
39 // expose PI on GPU build
40 #ifndef M_PI
41 #define M_PI 3.1415926535897932384626433832795028841971
42 #endif
43 
44 
45 #ifdef __cplusplus
46 extern "C" {
47 #endif
48 
49 /* builds a bit-string where 1 indicates a qubit is present in this list */
50 long long int getQubitBitMask(int* qubits, int numQubits) {
51 
52  long long int mask=0;
53  for (int i=0; i<numQubits; i++)
54  mask = mask | (1LL << qubits[i]);
55 
56  return mask;
57 }
58 
59 /* builds a bit-string where 1 indicates control qubits conditioned on 0 ('flipped') */
60 long long int getControlFlipMask(int* controlQubits, int* controlState, int numControlQubits) {
61 
62  long long int mask=0;
63  for (int i=0; i<numControlQubits; i++)
64  if (controlState[i] == 0)
65  mask = mask | (1LL << controlQubits[i]);
66 
67  return mask;
68 }
69 
70 void ensureIndsIncrease(int* ind1, int* ind2) {
71 
72  if (*ind1 > *ind2) {
73  int copy = *ind1;
74  *ind1 = *ind2;
75  *ind2 = copy;
76  }
77 }
78 
80 
81  return sqrt(vec.x*vec.x + vec.y*vec.y + vec.z*vec.z);
82 }
83 
85 
86  qreal mag = getVectorMagnitude(vec);
87  Vector unitVec = (Vector) {.x=vec.x/mag, .y=vec.y/mag, .z=vec.z/mag};
88  return unitVec;
89 }
90 
92 
93  Complex conjScalar;
94  conjScalar.real = scalar.real;
95  conjScalar.imag = - scalar.imag;
96  return conjScalar;
97 }
98 
99 #define macro_setConjugateMatrix(dest, src, dim) \
100  for (int i=0; i<dim; i++) \
101  for (int j=0; j<dim; j++) { \
102  dest.real[i][j] = src.real[i][j]; \
103  dest.imag[i][j] = - src.imag[i][j]; /* negative for conjugate */ \
104  }
106  ComplexMatrix2 conj;
107  macro_setConjugateMatrix(conj, src, 2);
108  return conj;
109 }
111  ComplexMatrix4 conj;
112  macro_setConjugateMatrix(conj, src, 4);
113  return conj;
114 }
116  int len = 1 << m.numQubits;
117  macro_setConjugateMatrix(m, m, len);
118 }
119 
120 void getComplexPairFromRotation(qreal angle, Vector axis, Complex* alpha, Complex* beta) {
121 
122  Vector unitAxis = getUnitVector(axis);
123  alpha->real = cos(angle/2.0);
124  alpha->imag = - sin(angle/2.0)*unitAxis.z;
125  beta->real = sin(angle/2.0)*unitAxis.y;
126  beta->imag = - sin(angle/2.0)*unitAxis.x;
127 }
128 
130 void getZYZRotAnglesFromComplexPair(Complex alpha, Complex beta, qreal* rz2, qreal* ry, qreal* rz1) {
131 
132  qreal alphaMag = sqrt(alpha.real*alpha.real + alpha.imag*alpha.imag);
133  *ry = 2.0 * acos(alphaMag);
134 
135  qreal alphaPhase = atan2(alpha.imag, alpha.real);
136  qreal betaPhase = atan2(beta.imag, beta.real);
137  *rz2 = - alphaPhase + betaPhase;
138  *rz1 = - alphaPhase - betaPhase;
139 }
140 
143 
144  qreal r0c0Phase = atan2(u.imag[0][0], u.real[0][0]);
145  qreal r1c1Phase = atan2(u.imag[1][1], u.real[1][1]);
146  *globalPhase = (r0c0Phase + r1c1Phase)/2.0;
147 
148  qreal cosPhase = cos(*globalPhase);
149  qreal sinPhase = sin(*globalPhase);
150  alpha->real = u.real[0][0]*cosPhase + u.imag[0][0]*sinPhase;
151  alpha->imag = u.imag[0][0]*cosPhase - u.real[0][0]*sinPhase;
152  beta->real = u.real[1][0]*cosPhase + u.imag[1][0]*sinPhase;
153  beta->imag = u.imag[1][0]*cosPhase - u.real[1][0]*sinPhase;
154 }
155 
156 void shiftIndices(int* indices, int numIndices, int shift) {
157  for (int j=0; j < numIndices; j++)
158  indices[j] += shift;
159 }
160 
161 void shiftSubregIndices(int* allInds, int* numIndsPerReg, int numRegs, int shift) {
162  int i=0;
163  for (int r=0; r<numRegs; r++)
164  for (int j=0; j<numIndsPerReg[r]; j++)
165  allInds[i++] += shift;
166 }
167 
168 int generateMeasurementOutcome(qreal zeroProb, qreal *outcomeProb) {
169 
170  // randomly choose outcome
171  int outcome;
172  if (zeroProb < REAL_EPS)
173  outcome = 1;
174  else if (1-zeroProb < REAL_EPS)
175  outcome = 0;
176  else
177  outcome = (genrand_real1() > zeroProb);
178 
179  // set probability of outcome
180  *outcomeProb = (outcome==0)? zeroProb : 1-zeroProb;
181 
182  return outcome;
183 }
184 
185 unsigned long int hashString(char *str){
186  unsigned long int hash = 5381;
187  int c;
188 
189  while ((c = *str++))
190  hash = ((hash << 5) + hash) + c; /* hash * 33 + c */
191 
192  return hash;
193 }
194 
195 void getQuESTDefaultSeedKey(unsigned long int *key){
196  // init MT random number generator with two keys -- time and pid
197  // for the MPI version, it is ok that all procs will get the same seed as random numbers will only be
198  // used by the master process
199 #if defined(_WIN32) && ! defined(__MINGW32__)
200 
201  unsigned long int pid = (unsigned long int) _getpid();
202  unsigned long int msecs = (unsigned long int) GetTickCount64();
203 
204  key[0] = msecs; key[1] = pid;
205 #else
206  struct timeval tv;
207  gettimeofday(&tv, NULL);
208 
209  double time_in_mill =
210  (tv.tv_sec) * 1000 + (tv.tv_usec) / 1000 ; // convert tv_sec & tv_usec to millisecond
211 
212  unsigned long int pid = getpid();
213  unsigned long int msecs = (unsigned long int) time_in_mill;
214 
215  key[0] = msecs; key[1] = pid;
216 #endif
217 }
218 
219 void reportState(Qureg qureg){
220  FILE *state;
221  char filename[100];
222  long long int index;
223  sprintf(filename, "state_rank_%d.csv", qureg.chunkId);
224  state = fopen(filename, "w");
225  if (qureg.chunkId==0) fprintf(state, "real, imag\n");
226 
227  for(index=0; index<qureg.numAmpsPerChunk; index++){
228  # if QuEST_PREC==1 || QuEST_PREC==2
229  fprintf(state, "%.12f, %.12f\n", qureg.stateVec.real[index], qureg.stateVec.imag[index]);
230  # elif QuEST_PREC == 4
231  fprintf(state, "%.12Lf, %.12Lf\n", qureg.stateVec.real[index], qureg.stateVec.imag[index]);
232  #endif
233  }
234  fclose(state);
235 }
236 
238  long long int numAmps = 1LL << qureg.numQubitsInStateVec;
239  long long int numAmpsPerRank = numAmps/qureg.numChunks;
240  if (qureg.chunkId==0){
241  printf("QUBITS:\n");
242  printf("Number of qubits is %d.\n", qureg.numQubitsInStateVec);
243  printf("Number of amps is %lld.\n", numAmps);
244  printf("Number of amps per rank is %lld.\n", numAmpsPerRank);
245  }
246 }
247 
248 qreal statevec_getProbAmp(Qureg qureg, long long int index){
249  qreal real = statevec_getRealAmp(qureg, index);
250  qreal imag = statevec_getImagAmp(qureg, index);
251  return real*real + imag*imag;
252 }
253 
254 void statevec_phaseShift(Qureg qureg, int targetQubit, qreal angle) {
255  Complex term;
256  term.real = cos(angle);
257  term.imag = sin(angle);
258  statevec_phaseShiftByTerm(qureg, targetQubit, term);
259 }
260 
261 void statevec_pauliZ(Qureg qureg, int targetQubit) {
262  Complex term;
263  term.real = -1;
264  term.imag = 0;
265  statevec_phaseShiftByTerm(qureg, targetQubit, term);
266 }
267 
268 void statevec_sGate(Qureg qureg, int targetQubit) {
269  Complex term;
270  term.real = 0;
271  term.imag = 1;
272  statevec_phaseShiftByTerm(qureg, targetQubit, term);
273 }
274 
275 void statevec_tGate(Qureg qureg, int targetQubit) {
276  Complex term;
277  term.real = 1/sqrt(2);
278  term.imag = 1/sqrt(2);
279  statevec_phaseShiftByTerm(qureg, targetQubit, term);
280 }
281 
282 void statevec_sGateConj(Qureg qureg, int targetQubit) {
283  Complex term;
284  term.real = 0;
285  term.imag = -1;
286  statevec_phaseShiftByTerm(qureg, targetQubit, term);
287 }
288 
289 void statevec_tGateConj(Qureg qureg, int targetQubit) {
290  Complex term;
291  term.real = 1/sqrt(2);
292  term.imag = -1/sqrt(2);
293  statevec_phaseShiftByTerm(qureg, targetQubit, term);
294 }
295 
296 void statevec_rotateX(Qureg qureg, int rotQubit, qreal angle){
297 
298  Vector unitAxis = {1, 0, 0};
299  statevec_rotateAroundAxis(qureg, rotQubit, angle, unitAxis);
300 }
301 
302 void statevec_rotateY(Qureg qureg, int rotQubit, qreal angle){
303 
304  Vector unitAxis = {0, 1, 0};
305  statevec_rotateAroundAxis(qureg, rotQubit, angle, unitAxis);
306 }
307 
308 void statevec_rotateZ(Qureg qureg, int rotQubit, qreal angle){
309 
310  Vector unitAxis = {0, 0, 1};
311  statevec_rotateAroundAxis(qureg, rotQubit, angle, unitAxis);
312 }
313 
314 void statevec_rotateAroundAxis(Qureg qureg, int rotQubit, qreal angle, Vector axis){
315 
316  Complex alpha, beta;
317  getComplexPairFromRotation(angle, axis, &alpha, &beta);
318  statevec_compactUnitary(qureg, rotQubit, alpha, beta);
319 }
320 
321 void statevec_rotateAroundAxisConj(Qureg qureg, int rotQubit, qreal angle, Vector axis){
322 
323  Complex alpha, beta;
324  getComplexPairFromRotation(angle, axis, &alpha, &beta);
325  alpha.imag *= -1;
326  beta.imag *= -1;
327  statevec_compactUnitary(qureg, rotQubit, alpha, beta);
328 }
329 
330 void statevec_controlledRotateAroundAxis(Qureg qureg, int controlQubit, int targetQubit, qreal angle, Vector axis){
331 
332  Complex alpha, beta;
333  getComplexPairFromRotation(angle, axis, &alpha, &beta);
334  statevec_controlledCompactUnitary(qureg, controlQubit, targetQubit, alpha, beta);
335 }
336 
337 void statevec_controlledRotateAroundAxisConj(Qureg qureg, int controlQubit, int targetQubit, qreal angle, Vector axis){
338 
339  Complex alpha, beta;
340  getComplexPairFromRotation(angle, axis, &alpha, &beta);
341  alpha.imag *= -1;
342  beta.imag *= -1;
343  statevec_controlledCompactUnitary(qureg, controlQubit, targetQubit, alpha, beta);
344 }
345 
346 void statevec_controlledRotateX(Qureg qureg, int controlQubit, int targetQubit, qreal angle){
347 
348  Vector unitAxis = {1, 0, 0};
349  statevec_controlledRotateAroundAxis(qureg, controlQubit, targetQubit, angle, unitAxis);
350 }
351 
352 void statevec_controlledRotateY(Qureg qureg, int controlQubit, int targetQubit, qreal angle){
353 
354  Vector unitAxis = {0, 1, 0};
355  statevec_controlledRotateAroundAxis(qureg, controlQubit, targetQubit, angle, unitAxis);
356 }
357 
358 void statevec_controlledRotateZ(Qureg qureg, int controlQubit, int targetQubit, qreal angle){
359 
360  Vector unitAxis = {0, 0, 1};
361  statevec_controlledRotateAroundAxis(qureg, controlQubit, targetQubit, angle, unitAxis);
362 }
363 
364 int statevec_measureWithStats(Qureg qureg, int measureQubit, qreal *outcomeProb) {
365 
366  qreal zeroProb = statevec_calcProbOfOutcome(qureg, measureQubit, 0);
367  int outcome = generateMeasurementOutcome(zeroProb, outcomeProb);
368  statevec_collapseToKnownProbOutcome(qureg, measureQubit, outcome, *outcomeProb);
369  return outcome;
370 }
371 
372 int densmatr_measureWithStats(Qureg qureg, int measureQubit, qreal *outcomeProb) {
373 
374  qreal zeroProb = densmatr_calcProbOfOutcome(qureg, measureQubit, 0);
375  int outcome = generateMeasurementOutcome(zeroProb, outcomeProb);
376  densmatr_collapseToKnownProbOutcome(qureg, measureQubit, outcome, *outcomeProb);
377  return outcome;
378 }
379 
381 
382  Complex innerProd = statevec_calcInnerProduct(qureg, pureState);
383  qreal innerProdMag = innerProd.real*innerProd.real + innerProd.imag*innerProd.imag;
384  return innerProdMag;
385 }
386 
387 void statevec_sqrtSwapGate(Qureg qureg, int qb1, int qb2) {
388 
389  ComplexMatrix4 u = (ComplexMatrix4) {.real={{0}}, .imag={{0}}};
390  u.real[0][0]=1;
391  u.real[3][3]=1;
392  u.real[1][1] = .5; u.imag[1][1] = .5;
393  u.real[1][2] = .5; u.imag[1][2] =-.5;
394  u.real[2][1] = .5; u.imag[2][1] =-.5;
395  u.real[2][2] = .5; u.imag[2][2] = .5;
396 
397  statevec_twoQubitUnitary(qureg, qb1, qb2, u);
398 }
399 
400 void statevec_sqrtSwapGateConj(Qureg qureg, int qb1, int qb2) {
401 
402  ComplexMatrix4 u = (ComplexMatrix4) {.real={{0}}, .imag={{0}}};
403  u.real[0][0]=1;
404  u.real[3][3]=1;
405  u.real[1][1] = .5; u.imag[1][1] =-.5;
406  u.real[1][2] = .5; u.imag[1][2] = .5;
407  u.real[2][1] = .5; u.imag[2][1] = .5;
408  u.real[2][2] = .5; u.imag[2][2] =-.5;
409 
410  statevec_twoQubitUnitary(qureg, qb1, qb2, u);
411 }
412 
415  Qureg qureg, int* targetQubits, enum pauliOpType* targetPaulis, int numTargets, qreal angle,
416  int applyConj
417 ) {
418  qreal fac = 1/sqrt(2);
419  Complex uRxAlpha = {.real = fac, .imag = 0}; // Rx(pi/2)* rotates Z -> Y
420  Complex uRxBeta = {.real = 0, .imag = (applyConj)? fac : -fac};
421  Complex uRyAlpha = {.real = fac, .imag = 0}; // Ry(-pi/2) rotates Z -> X
422  Complex uRyBeta = {.real = -fac, .imag = 0};
423 
424  // mask may be modified to remove superfluous Identity ops
425  long long int mask = getQubitBitMask(targetQubits, numTargets);
426 
427  // rotate basis so that exp(Z) will effect exp(Y) and exp(X)
428  for (int t=0; t < numTargets; t++) {
429  if (targetPaulis[t] == PAULI_I)
430  mask -= 1LL << targetQubits[t]; // remove target from mask
431  if (targetPaulis[t] == PAULI_X)
432  statevec_compactUnitary(qureg, targetQubits[t], uRyAlpha, uRyBeta);
433  if (targetPaulis[t] == PAULI_Y)
434  statevec_compactUnitary(qureg, targetQubits[t], uRxAlpha, uRxBeta);
435  // (targetPaulis[t] == 3) is Z basis
436  }
437 
438  // does nothing if there are no qubits to 'rotate'
439  if (mask != 0)
440  statevec_multiRotateZ(qureg, mask, (applyConj)? -angle : angle);
441 
442  // undo X and Y basis rotations
443  uRxBeta.imag *= -1;
444  uRyBeta.real *= -1;
445  for (int t=0; t < numTargets; t++) {
446  if (targetPaulis[t] == PAULI_X)
447  statevec_compactUnitary(qureg, targetQubits[t], uRyAlpha, uRyBeta);
448  if (targetPaulis[t] == PAULI_Y)
449  statevec_compactUnitary(qureg, targetQubits[t], uRxAlpha, uRxBeta);
450  }
451 }
452 
454  Qureg qureg, long long int ctrlMask, int* targetQubits, enum pauliOpType* targetPaulis, int numTargets, qreal angle,
455  int applyConj
456 ) {
457  qreal fac = 1/sqrt(2);
458  qreal sgn = (applyConj)? 1 : -1;
459  ComplexMatrix2 uRx = {.real={{fac,0},{0,fac}}, .imag={{0,sgn*fac},{sgn*fac,0}}}; // Rx(pi/2)* rotates Z -> Y
460  ComplexMatrix2 uRy = {.real={{fac,fac},{-fac,fac}}, .imag={{0,0},{0,0}}}; // Ry(-pi/2) rotates Z -> X
461 
462  // this function is controlled on the all-one state, so no ctrl flips
463  long long int ctrlFlipMask = 0;
464 
465  // mask may be modified to remove superfluous Identity ops
466  long long int targMask = getQubitBitMask(targetQubits, numTargets);
467 
468  // rotate basis so that exp(Z) will effect exp(Y) and exp(X)
469  for (int t=0; t < numTargets; t++) {
470  if (targetPaulis[t] == PAULI_I)
471  targMask -= 1LL << targetQubits[t]; // remove target from mask
472  if (targetPaulis[t] == PAULI_X)
473  statevec_multiControlledUnitary(qureg, ctrlMask, ctrlFlipMask, targetQubits[t], uRy);
474  if (targetPaulis[t] == PAULI_Y)
475  statevec_multiControlledUnitary(qureg, ctrlMask, ctrlFlipMask, targetQubits[t], uRx);
476  // (targetPaulis[t] == 3) is Z basis
477  }
478 
479  // does nothing if there are no qubits to 'rotate'
480  if (targMask != 0)
481  statevec_multiControlledMultiRotateZ(qureg, ctrlMask, targMask, (applyConj)? -angle : angle);
482 
483  // undo X and Y basis rotations
484  uRx.imag[0][1] *= -1; uRx.imag[1][0] *= -1;
485  uRy.real[0][1] *= -1; uRy.real[1][0] *= -1;
486  for (int t=0; t < numTargets; t++) {
487  if (targetPaulis[t] == PAULI_X)
488  statevec_multiControlledUnitary(qureg, ctrlMask, ctrlFlipMask, targetQubits[t], uRy);
489  if (targetPaulis[t] == PAULI_Y)
490  statevec_multiControlledUnitary(qureg, ctrlMask, ctrlFlipMask, targetQubits[t], uRx);
491  }
492 }
493 
494 /* produces both pauli|qureg> or pauli * qureg (as a density matrix) */
495 void statevec_applyPauliProd(Qureg workspace, int* targetQubits, enum pauliOpType* pauliCodes, int numTargets) {
496 
497  for (int i=0; i < numTargets; i++) {
498  // (pauliCodes[i] == PAULI_I) applies no operation
499  if (pauliCodes[i] == PAULI_X)
500  statevec_pauliX(workspace, targetQubits[i]);
501  if (pauliCodes[i] == PAULI_Y)
502  statevec_pauliY(workspace, targetQubits[i]);
503  if (pauliCodes[i] == PAULI_Z)
504  statevec_pauliZ(workspace, targetQubits[i]);
505  }
506 }
507 
508 // <pauli> = <qureg|pauli|qureg> = qureg . pauli(qureg)
509 qreal statevec_calcExpecPauliProd(Qureg qureg, int* targetQubits, enum pauliOpType* pauliCodes, int numTargets, Qureg workspace) {
510 
511  statevec_cloneQureg(workspace, qureg);
512  statevec_applyPauliProd(workspace, targetQubits, pauliCodes, numTargets);
513 
514  // compute the expected value
515  qreal value;
516  if (qureg.isDensityMatrix)
517  value = densmatr_calcTotalProb(workspace); // Trace(ops qureg)
518  else
519  value = statevec_calcInnerProduct(workspace, qureg).real; // <qureg|ops|qureg>
520 
521  return value;
522 }
523 
524 qreal statevec_calcExpecPauliSum(Qureg qureg, enum pauliOpType* allCodes, qreal* termCoeffs, int numSumTerms, Qureg workspace) {
525 
526  int numQb = qureg.numQubitsRepresented;
527  int targs[100]; // [numQb];
528  for (int q=0; q < numQb; q++)
529  targs[q] = q;
530 
531  qreal value = 0;
532  for (int t=0; t < numSumTerms; t++)
533  value += termCoeffs[t] * statevec_calcExpecPauliProd(qureg, targs, &allCodes[t*numQb], numQb, workspace);
534 
535  return value;
536 }
537 
538 void statevec_applyPauliSum(Qureg inQureg, enum pauliOpType* allCodes, qreal* termCoeffs, int numSumTerms, Qureg outQureg) {
539 
540  int numQb = inQureg.numQubitsRepresented;
541  int targs[100]; // [numQb];
542  for (int q=0; q < numQb; q++)
543  targs[q] = q;
544 
545  statevec_initBlankState(outQureg);
546 
547  for (int t=0; t < numSumTerms; t++) {
548  Complex coef = (Complex) {.real=termCoeffs[t], .imag=0};
549  Complex iden = (Complex) {.real=1, .imag=0};
550  Complex zero = (Complex) {.real=0, .imag=0};
551 
552  // outQureg += coef paulis(inQureg)
553  statevec_applyPauliProd(inQureg, targs, &allCodes[t*numQb], numQb);
554  statevec_setWeightedQureg(coef, inQureg, iden, outQureg, zero, outQureg);
555 
556  // undero paulis(inQureg), exploiting XX=YY=ZZ=I
557  statevec_applyPauliProd(inQureg, targs, &allCodes[t*numQb], numQb);
558  }
559 }
560 
561 void statevec_twoQubitUnitary(Qureg qureg, int targetQubit1, int targetQubit2, ComplexMatrix4 u) {
562 
563  long long int ctrlMask = 0;
564  statevec_multiControlledTwoQubitUnitary(qureg, ctrlMask, targetQubit1, targetQubit2, u);
565 }
566 
567 void statevec_controlledTwoQubitUnitary(Qureg qureg, int controlQubit, int targetQubit1, int targetQubit2, ComplexMatrix4 u) {
568 
569  long long int ctrlMask = 1LL << controlQubit;
570  statevec_multiControlledTwoQubitUnitary(qureg, ctrlMask, targetQubit1, targetQubit2, u);
571 }
572 
573 void statevec_multiQubitUnitary(Qureg qureg, int* targets, int numTargets, ComplexMatrixN u) {
574 
575  long long int ctrlMask = 0;
576  statevec_multiControlledMultiQubitUnitary(qureg, ctrlMask, targets, numTargets, u);
577 }
578 
579 void statevec_controlledMultiQubitUnitary(Qureg qureg, int ctrl, int* targets, int numTargets, ComplexMatrixN u) {
580 
581  long long int ctrlMask = 1LL << ctrl;
582  statevec_multiControlledMultiQubitUnitary(qureg, ctrlMask, targets, numTargets, u);
583 }
584 
585 #define macro_populateKrausOperator(superOp, ops, numOps, opDim) \
586  /* clear the superop */ \
587  for (int r=0; r < (opDim)*(opDim); r++) \
588  for (int c=0; c < (opDim)*(opDim); c++) { \
589  superOp->real[r][c] = 0; \
590  superOp->imag[r][c] = 0; \
591  } \
592  /* add each op's contribution to the superop */ \
593  for (int n=0; n<(numOps); n++) \
594  /* superop += conjugate(op) (x) op, where (x) is a tensor product */ \
595  for (int i = 0; i < (opDim); i++) \
596  for (int j = 0; j < (opDim); j++) \
597  for (int k = 0; k < (opDim); k++) \
598  for (int l = 0; l < (opDim); l++) { \
599  superOp->real[i*(opDim) + k][j*(opDim) + l] += \
600  ops[n].real[i][j]*ops[n].real[k][l] + \
601  ops[n].imag[i][j]*ops[n].imag[k][l]; \
602  superOp->imag[i*(opDim) + k][j*(opDim) + l] += \
603  ops[n].real[i][j]*ops[n].imag[k][l] - \
604  ops[n].imag[i][j]*ops[n].real[k][l]; \
605  }
606 
608  int opDim = 2;
609  macro_populateKrausOperator(superOp, ops, numOps, opDim);
610 }
612  int opDim = 4;
613  macro_populateKrausOperator(superOp, ops, numOps, opDim);
614 }
616  int opDim = 1 << ops[0].numQubits;
617  macro_populateKrausOperator(superOp, ops, numOps, opDim);
618 }
619 
620 void densmatr_applyKrausSuperoperator(Qureg qureg, int target, ComplexMatrix4 superOp) {
621 
622  long long int ctrlMask = 0;
623  statevec_multiControlledTwoQubitUnitary(qureg, ctrlMask, target, target + qureg.numQubitsRepresented, superOp);
624 }
625 
626 void densmatr_applyTwoQubitKrausSuperoperator(Qureg qureg, int target1, int target2, ComplexMatrixN superOp) {
627 
628  long long int ctrlMask = 0;
629  int numQb = qureg.numQubitsRepresented;
630  int allTargets[4] = {target1, target2, target1+numQb, target2+numQb};
631  statevec_multiControlledMultiQubitUnitary(qureg, ctrlMask, allTargets, 4, superOp);
632 }
633 
634 void densmatr_applyMultiQubitKrausSuperoperator(Qureg qureg, int *targets, int numTargets, ComplexMatrixN superOp) {
635  long long int ctrlMask = 0;
636  int allTargets[200]; // [2*numTargets];
637  for (int t=0; t < numTargets; t++) {
638  allTargets[t] = targets[t];
639  allTargets[t+numTargets] = targets[t] + qureg.numQubitsRepresented;
640  }
641  statevec_multiControlledMultiQubitUnitary(qureg, ctrlMask, allTargets, 2*numTargets, superOp);
642 }
643 
644 void densmatr_mixKrausMap(Qureg qureg, int target, ComplexMatrix2 *ops, int numOps) {
645 
646  ComplexMatrix4 superOp;
647  populateKrausSuperOperator2(&superOp, ops, numOps);
648  densmatr_applyKrausSuperoperator(qureg, target, superOp);
649 }
650 
651 #ifndef _WIN32
653  int numQubits, qreal re[][1<<numQubits], qreal im[][1<<numQubits],
654  qreal** reStorage, qreal** imStorage
655 ) {
656  ComplexMatrixN m;
657  m.numQubits = numQubits;
658  m.real = reStorage;
659  m.imag = imStorage;
660 
661  int len = 1<<numQubits;
662  for (int i=0; i<len; i++) {
663  m.real[i] = re[i];
664  m.imag[i] = im[i];
665  }
666  return m;
667 }
668 
669 #define macro_initialiseStackComplexMatrixN(matrix, numQubits, real, imag) \
670  /* reStorage_ and imStorage_ must not exist in calling scope */ \
671  qreal* reStorage_[1<<(numQubits)]; \
672  qreal* imStorage_[1<<(numQubits)]; \
673  matrix = bindArraysToStackComplexMatrixN((numQubits), real, imag, reStorage_, imStorage_);
674 
675 #define macro_allocStackComplexMatrixN(matrix, numQubits) \
676  /* reArr_, imArr_, reStorage_, and imStorage_ must not exist in calling scope */ \
677  qreal reArr_[1<<(numQubits)][1<<(numQubits)]; \
678  qreal imArr_[1<<(numQubits)][1<<(numQubits)]; \
679  macro_initialiseStackComplexMatrixN(matrix, (numQubits), reArr_, imArr_);
680 #endif
681 
682 void densmatr_mixTwoQubitKrausMap(Qureg qureg, int target1, int target2, ComplexMatrix4 *ops, int numOps) {
683 
684  // if NOT on Windows, allocate ComplexN on stack
685  #ifndef _WIN32
686  ComplexMatrixN superOp;
687  macro_allocStackComplexMatrixN(superOp, 4);
688  populateKrausSuperOperator4(&superOp, ops, numOps);
689  densmatr_applyTwoQubitKrausSuperoperator(qureg, target1, target2, superOp);
690 
691  // but on Windows, we MUST allocated dynamically
692  #else
694  populateKrausSuperOperator4(&superOp, ops, numOps);
695  densmatr_applyTwoQubitKrausSuperoperator(qureg, target1, target2, superOp);
696  destroyComplexMatrixN(superOp);
697 
698  #endif
699 }
700 
701 void densmatr_mixMultiQubitKrausMap(Qureg qureg, int* targets, int numTargets, ComplexMatrixN* ops, int numOps) {
702 
703  ComplexMatrixN superOp;
704 
705  /* superOp will contain 2^(4 numTargets) complex numbers.
706  * At double precision, superOp will cost additional memory:
707  * numTargs=1 -> 0.25 KiB
708  * numTargs=2 -> 4 KiB
709  * numTargs=3 -> 64 KiB
710  * numTargs=4 -> 1 MiB
711  * numTargs=5 -> 16 MiB.
712  * At quad precision (usually 10 B per number, but possibly 16 B due to alignment),
713  * this costs at most double.
714  *
715  * Hence, if superOp is kept in the stack, numTargs >= 4 would exceed Windows' 1 MB
716  * maximum stack-space allocation (numTargs >= 5 exceeding Linux' 8 MB). Therefore,
717  * for numTargets < 4, superOp will be kept in the stack, else in the heap
718  */
719 
720  // if NOT on Windows, allocate ComplexN on stack depending on size
721  #ifndef _WIN32
722  if (numTargets < 4) {
723  // everything must live in 'if' since this macro declares local vars
724  macro_allocStackComplexMatrixN(superOp, 2*numTargets);
725  populateKrausSuperOperatorN(&superOp, ops, numOps);
726  densmatr_applyMultiQubitKrausSuperoperator(qureg, targets, numTargets, superOp);
727  }
728  else {
729  superOp = createComplexMatrixN(2*numTargets);
730  populateKrausSuperOperatorN(&superOp, ops, numOps);
731  densmatr_applyMultiQubitKrausSuperoperator(qureg, targets, numTargets, superOp);
732  destroyComplexMatrixN(superOp);
733  }
734  // on Windows, we must always create in heap
735  #else
736  superOp = createComplexMatrixN(2*numTargets);
737  populateKrausSuperOperatorN(&superOp, ops, numOps);
738  densmatr_applyMultiQubitKrausSuperoperator(qureg, targets, numTargets, superOp);
739  destroyComplexMatrixN(superOp);
740  #endif
741 }
742 
743 void densmatr_mixPauli(Qureg qureg, int qubit, qreal probX, qreal probY, qreal probZ) {
744 
745  // convert pauli probabilities into Kraus map
746  const int numOps = 4;
747  ComplexMatrix2 ops[4]; // [numOps];
748  for (int n=0; n < numOps; n++)
749  ops[n] = (ComplexMatrix2) {.real={{0}}, .imag={{0}}};
750 
751  qreal facs[4] = { // literal numOps=4 for valid initialisation
752  sqrt(1-(probX + probY + probZ)),
753  sqrt(probX),
754  sqrt(probY),
755  sqrt(probZ)
756  };
757  ops[0].real[0][0] = facs[0]; ops[0].real[1][1] = facs[0];
758  ops[1].real[0][1] = facs[1]; ops[1].real[1][0] = facs[1];
759  ops[2].imag[0][1] = -facs[2]; ops[2].imag[1][0] = facs[2];
760  ops[3].real[0][0] = facs[3]; ops[3].real[1][1] = -facs[3];
761 
762  densmatr_mixKrausMap(qureg, qubit, ops, numOps);
763 }
764 
765 void applyExponentiatedPauliHamil(Qureg qureg, PauliHamil hamil, qreal fac, int reverse) {
766 
767  /* applies a first-order one-repetition approximation of exp(-i fac H)
768  * to qureg. Letting H = sum_j c_j h_j, it does this via
769  * exp(-i fac H) ~ prod_j exp(-i fac c_j h_j), where each inner exp
770  * is performed with multiRotatePauli (with pre-factor 2).
771  */
772 
773  // prepare targets for multiRotatePauli
774  // (all qubits; actual targets are determined by Pauli codes)
775  int vecTargs[100]; // [hamil.numQubits];
776  int densTargs[100]; // [hamil.numQubits];
777  for (int q=0; q<hamil.numQubits; q++) {
778  vecTargs[q] = q;
779  densTargs[q] = q + hamil.numQubits;
780  }
781 
782  for (int i=0; i<hamil.numSumTerms; i++) {
783 
784  int t=i;
785  if (reverse)
786  t=hamil.numSumTerms-1-i;
787 
788  qreal angle = 2*fac*hamil.termCoeffs[t];
789 
791  qureg, vecTargs, &(hamil.pauliCodes[t*hamil.numQubits]),
792  hamil.numQubits, angle, 0);
793 
794  if (qureg.isDensityMatrix)
796  qureg, densTargs, &(hamil.pauliCodes[t*hamil.numQubits]),
797  hamil.numQubits, angle, 1);
798 
799  // record qasm
800  char buff[1024];
801  int b=0;
802  for (int q=0; q<hamil.numQubits; q++) {
803  enum pauliOpType op = hamil.pauliCodes[q + t*hamil.numQubits];
804 
805  char p = 'I';
806  if (op == PAULI_X) p = 'X';
807  if (op == PAULI_Y) p = 'Y';
808  if (op == PAULI_Z) p = 'Z';
809  buff[b++] = p;
810  buff[b++] = ' ';
811  }
812  buff[b] = '\0';
813 
814  qasm_recordComment(qureg,
815  "Here, a multiRotatePauli with angle %g and paulis %s was applied.",
816  angle, buff);
817  }
818 }
819 
820 void applySymmetrizedTrotterCircuit(Qureg qureg, PauliHamil hamil, qreal time, int order) {
821 
822  if (order == 1) {
823  applyExponentiatedPauliHamil(qureg, hamil, time, 0);
824  }
825  else if (order == 2) {
826  applyExponentiatedPauliHamil(qureg, hamil, time/2., 0);
827  applyExponentiatedPauliHamil(qureg, hamil, time/2., 1);
828  }
829  else {
830  qreal p = 1. / (4 - pow(4, 1./(order-1)));
831  int lower = order-2;
832  applySymmetrizedTrotterCircuit(qureg, hamil, p*time, lower);
833  applySymmetrizedTrotterCircuit(qureg, hamil, p*time, lower);
834  applySymmetrizedTrotterCircuit(qureg, hamil, (1-4*p)*time, lower);
835  applySymmetrizedTrotterCircuit(qureg, hamil, p*time, lower);
836  applySymmetrizedTrotterCircuit(qureg, hamil, p*time, lower);
837  }
838 }
839 
840 void agnostic_applyTrotterCircuit(Qureg qureg, PauliHamil hamil, qreal time, int order, int reps) {
841 
842  if (time == 0)
843  return;
844 
845  for (int r=0; r<reps; r++)
846  applySymmetrizedTrotterCircuit(qureg, hamil, time/reps, order);
847 }
848 
849 void agnostic_applyQFT(Qureg qureg, int* qubits, int numQubits) {
850 
851  int densShift = qureg.numQubitsRepresented;
852 
853  // start with top/left-most qubit, work down
854  for (int q=numQubits-1; q >= 0; q--) {
855 
856  // H
857  statevec_hadamard(qureg, qubits[q]);
858  if (qureg.isDensityMatrix)
859  statevec_hadamard(qureg, qubits[q] + densShift);
860  qasm_recordGate(qureg, GATE_HADAMARD, qubits[q]);
861 
862  if (q == 0)
863  break;
864 
865  // succession of C-phases, control on qubits[q], targeting each of
866  // qubits[q-1], qubits[q-1], ... qubits[0]. This can be expressed by
867  // a single invocation of applyNamedPhaseFunc product
868 
869  int numRegs = 2;
870  int numQubitsPerReg[2] = {q, 1};
871  int regs[100]; // [q+1];
872  for (int i=0; i<q+1; i++)
873  regs[i] = qubits[i]; // qubits[q] is in own register
874 
875  int numParams = 1;
876  qreal params[1] = { M_PI / (1 << q) };
877 
878  int conj = 0;
880  qureg, regs, numQubitsPerReg, numRegs,
881  UNSIGNED, SCALED_PRODUCT, params, numParams,
882  NULL, NULL, 0,
883  conj);
884  if (qureg.isDensityMatrix) {
885  conj = 1;
886  shiftSubregIndices(regs, numQubitsPerReg, numRegs, densShift);
888  qureg, regs, numQubitsPerReg, numRegs,
889  UNSIGNED, SCALED_PRODUCT, params, numParams,
890  NULL, NULL, 0,
891  conj);
892  shiftSubregIndices(regs, numQubitsPerReg, numRegs, - densShift);
893  }
895  qureg, regs, numQubitsPerReg, numRegs,
896  UNSIGNED, SCALED_PRODUCT, params, numParams,
897  NULL, NULL, 0);
898  }
899 
900  // final swaps
901  for (int i=0; i<(numQubits/2); i++) {
902 
903  int qb1 = qubits[i];
904  int qb2 = qubits[numQubits-i-1];
905 
906  statevec_swapQubitAmps(qureg, qb1, qb2);
907  if (qureg.isDensityMatrix)
908  statevec_swapQubitAmps(qureg, qb1 + densShift, qb2 + densShift);
909  qasm_recordControlledGate(qureg, GATE_SWAP, qb1, qb2);
910  }
911 }
912 
913 #ifdef __cplusplus
914 }
915 #endif
void statevec_phaseShift(Qureg qureg, int targetQubit, qreal angle)
Definition: QuEST_common.c:254
Represents a 3-vector of real numbers.
Definition: QuEST.h:198
void statevec_sGate(Qureg qureg, int targetQubit)
Definition: QuEST_common.c:268
void statevec_pauliZ(Qureg qureg, int targetQubit)
Definition: QuEST_common.c:261
pauliOpType
Codes for specifying Pauli operators.
Definition: QuEST.h:96
void statevec_rotateX(Qureg qureg, int rotQubit, qreal angle)
Definition: QuEST_common.c:296
void reportQuregParams(Qureg qureg)
Report metainformation about a set of qubits: number of qubits, number of probability amplitudes.
Definition: QuEST_common.c:237
void agnostic_applyQFT(Qureg qureg, int *qubits, int numQubits)
Definition: QuEST_common.c:849
#define macro_setConjugateMatrix(dest, src, dim)
Definition: QuEST_common.c:99
qreal real[4][4]
Definition: QuEST.h:177
void densmatr_mixKrausMap(Qureg qureg, int target, ComplexMatrix2 *ops, int numOps)
Definition: QuEST_common.c:644
@ PAULI_Z
Definition: QuEST.h:96
void applyExponentiatedPauliHamil(Qureg qureg, PauliHamil hamil, qreal fac, int reverse)
Definition: QuEST_common.c:765
void populateKrausSuperOperator4(ComplexMatrixN *superOp, ComplexMatrix4 *ops, int numOps)
Definition: QuEST_common.c:611
void populateKrausSuperOperatorN(ComplexMatrixN *superOp, ComplexMatrixN *ops, int numOps)
Definition: QuEST_common.c:615
void statevec_multiControlledMultiRotatePauli(Qureg qureg, long long int ctrlMask, int *targetQubits, enum pauliOpType *targetPaulis, int numTargets, qreal angle, int applyConj)
Definition: QuEST_common.c:453
void destroyComplexMatrixN(ComplexMatrixN m)
Destroy a ComplexMatrixN instance created with createComplexMatrixN()
Definition: QuEST.c:1369
void shiftIndices(int *indices, int numIndices, int shift)
Definition: QuEST_common.c:156
qreal statevec_calcExpecPauliProd(Qureg qureg, int *targetQubits, enum pauliOpType *pauliCodes, int numTargets, Qureg workspace)
Definition: QuEST_common.c:509
void statevec_controlledTwoQubitUnitary(Qureg qureg, int controlQubit, int targetQubit1, int targetQubit2, ComplexMatrix4 u)
Definition: QuEST_common.c:567
@ PAULI_I
Definition: QuEST.h:96
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
qreal statevec_calcExpecPauliSum(Qureg qureg, enum pauliOpType *allCodes, qreal *termCoeffs, int numSumTerms, Qureg workspace)
Definition: QuEST_common.c:524
void statevec_twoQubitUnitary(Qureg qureg, int targetQubit1, int targetQubit2, ComplexMatrix4 u)
Definition: QuEST_common.c:561
void statevec_tGateConj(Qureg qureg, int targetQubit)
Definition: QuEST_common.c:289
qreal z
Definition: QuEST.h:200
ComplexMatrix4 getConjugateMatrix4(ComplexMatrix4 src)
Definition: QuEST_common.c:110
int numChunks
Number of chunks the state vector is broken up into – the number of MPI processes used.
Definition: QuEST.h:338
void statevec_applyPauliProd(Qureg workspace, int *targetQubits, enum pauliOpType *pauliCodes, int numTargets)
Definition: QuEST_common.c:495
void statevec_multiControlledMultiQubitUnitary(Qureg qureg, long long int ctrlMask, int *targs, int numTargs, ComplexMatrixN u)
This calls swapQubitAmps only when it would involve a distributed communication; if the qubit chunks ...
void statevec_multiRotateZ(Qureg qureg, long long int mask, qreal angle)
Definition: QuEST_cpu.c:3316
void getQuESTDefaultSeedKey(unsigned long int *key)
Definition: QuEST_common.c:195
qreal statevec_calcProbOfOutcome(Qureg qureg, int measureQubit, int outcome)
void densmatr_mixMultiQubitKrausMap(Qureg qureg, int *targets, int numTargets, ComplexMatrixN *ops, int numOps)
Definition: QuEST_common.c:701
Complex getConjugateScalar(Complex scalar)
Definition: QuEST_common.c:91
@ GATE_HADAMARD
Definition: QuEST_qasm.h:26
ComplexMatrix2 getConjugateMatrix2(ComplexMatrix2 src)
Definition: QuEST_common.c:105
Vector getUnitVector(Vector vec)
Definition: QuEST_common.c:84
Represents a 4x4 matrix of complex numbers.
Definition: QuEST.h:175
@ UNSIGNED
Definition: QuEST.h:269
void getComplexPairFromRotation(qreal angle, Vector axis, Complex *alpha, Complex *beta)
Definition: QuEST_common.c:120
void statevec_pauliY(Qureg qureg, int targetQubit)
ComplexMatrixN bindArraysToStackComplexMatrixN(int numQubits, qreal re[][1<< numQubits], qreal im[][1<< numQubits], qreal **reStorage, qreal **imStorage)
Definition: QuEST_common.c:652
Represents a general 2^N by 2^N matrix of complex numbers.
Definition: QuEST.h:186
void statevec_controlledMultiQubitUnitary(Qureg qureg, int ctrl, int *targets, int numTargets, ComplexMatrixN u)
Definition: QuEST_common.c:579
qreal statevec_calcFidelity(Qureg qureg, Qureg pureState)
Definition: QuEST_common.c:380
#define qreal
#define macro_allocStackComplexMatrixN(matrix, numQubits)
Definition: QuEST_common.c:675
void statevec_sqrtSwapGate(Qureg qureg, int qb1, int qb2)
Definition: QuEST_common.c:387
@ PAULI_X
Definition: QuEST.h:96
qreal statevec_getProbAmp(Qureg qureg, long long int index)
Definition: QuEST_common.c:248
void reportState(Qureg qureg)
Print the current state vector of probability amplitudes for a set of qubits to file.
Definition: QuEST_common.c:219
void statevec_swapQubitAmps(Qureg qureg, int qb1, int qb2)
int numQubitsInStateVec
Number of qubits in the state-vector - this is double the number represented for mixed states.
Definition: QuEST.h:329
qreal densmatr_calcTotalProb(Qureg qureg)
void densmatr_collapseToKnownProbOutcome(Qureg qureg, int measureQubit, int outcome, qreal outcomeProb)
Renorms (/prob) every | * outcome * >< * outcome * | state, setting all others to zero.
Definition: QuEST_cpu.c:791
int chunkId
The position of the chunk of the state vector held by this process in the full state vector.
Definition: QuEST.h:336
void qasm_recordControlledGate(Qureg qureg, TargetGate gate, int controlQubit, int targetQubit)
Definition: QuEST_qasm.c:239
qreal y
Definition: QuEST.h:200
unsigned long int hashString(char *str)
Definition: QuEST_common.c:185
qreal imag[2][2]
Definition: QuEST.h:140
void statevec_applyParamNamedPhaseFuncOverrides(Qureg qureg, int *qubits, int *numQubitsPerReg, int numRegs, enum bitEncoding encoding, enum phaseFunc functionNameCode, qreal *params, int numParams, long long int *overrideInds, qreal *overridePhases, int numOverrides, int conj)
Definition: QuEST_cpu.c:4446
void statevec_multiControlledMultiRotateZ(Qureg qureg, long long int ctrlMask, long long int targMask, qreal angle)
Definition: QuEST_cpu.c:3358
qreal x
Definition: QuEST.h:200
int generateMeasurementOutcome(qreal zeroProb, qreal *outcomeProb)
Definition: QuEST_common.c:168
long long int numAmpsPerChunk
Number of probability amplitudes held in stateVec by this process In the non-MPI version,...
Definition: QuEST.h:332
qreal * termCoeffs
The real coefficient of each Pauli product. This is an array of length PauliHamil....
Definition: QuEST.h:283
void agnostic_applyTrotterCircuit(Qureg qureg, PauliHamil hamil, qreal time, int order, int reps)
Definition: QuEST_common.c:840
void densmatr_applyKrausSuperoperator(Qureg qureg, int target, ComplexMatrix4 superOp)
Definition: QuEST_common.c:620
enum pauliOpType * pauliCodes
The Pauli operators acting on each qubit, flattened over every operator.
Definition: QuEST.h:281
void statevec_tGate(Qureg qureg, int targetQubit)
Definition: QuEST_common.c:275
void statevec_initBlankState(Qureg qureg)
Definition: QuEST_cpu.c:1464
void populateKrausSuperOperator2(ComplexMatrix4 *superOp, ComplexMatrix2 *ops, int numOps)
Definition: QuEST_common.c:607
@ SCALED_PRODUCT
Definition: QuEST.h:233
void statevec_cloneQureg(Qureg targetQureg, Qureg copyQureg)
works for both statevectors and density matrices
Definition: QuEST_cpu.c:1572
void statevec_rotateAroundAxis(Qureg qureg, int rotQubit, qreal angle, Vector axis)
Definition: QuEST_common.c:314
qreal imag[4][4]
Definition: QuEST.h:178
void setConjugateMatrixN(ComplexMatrixN m)
Definition: QuEST_common.c:115
void densmatr_mixTwoQubitKrausMap(Qureg qureg, int target1, int target2, ComplexMatrix4 *ops, int numOps)
Definition: QuEST_common.c:682
void statevec_multiControlledTwoQubitUnitary(Qureg qureg, long long int ctrlMask, int targetQubit1, int targetQubit2, ComplexMatrix4 u)
This calls swapQubitAmps only when it would involve a distributed communication; if the qubit chunks ...
int numSumTerms
The number of terms in the weighted sum, or the number of Pauli products.
Definition: QuEST.h:285
long long int getQubitBitMask(int *qubits, int numQubits)
Definition: QuEST_common.c:50
void statevec_compactUnitary(Qureg qureg, int targetQubit, Complex alpha, Complex beta)
@ 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
double genrand_real1(void)
Definition: mt19937ar.c:150
void statevec_multiRotatePauli(Qureg qureg, int *targetQubits, enum pauliOpType *targetPaulis, int numTargets, qreal angle, int applyConj)
applyConj=1 will apply conjugate operation, else applyConj=0
Definition: QuEST_common.c:414
void qasm_recordComment(Qureg qureg, char *comment,...)
Definition: QuEST_qasm.c:121
void statevec_multiQubitUnitary(Qureg qureg, int *targets, int numTargets, ComplexMatrixN u)
Definition: QuEST_common.c:573
qreal ** real
Definition: QuEST.h:189
void statevec_controlledRotateAroundAxis(Qureg qureg, int controlQubit, int targetQubit, qreal angle, Vector axis)
Definition: QuEST_common.c:330
void statevec_applyPauliSum(Qureg inQureg, enum pauliOpType *allCodes, qreal *termCoeffs, int numSumTerms, Qureg outQureg)
Definition: QuEST_common.c:538
Complex statevec_calcInnerProduct(Qureg bra, Qureg ket)
Terrible code which unnecessarily individually computes and sums the real and imaginary components of...
void statevec_rotateAroundAxisConj(Qureg qureg, int rotQubit, qreal angle, Vector axis)
Definition: QuEST_common.c:321
void statevec_sGateConj(Qureg qureg, int targetQubit)
Definition: QuEST_common.c:282
void statevec_setWeightedQureg(Complex fac1, Qureg qureg1, Complex fac2, Qureg qureg2, Complex facOut, Qureg out)
Definition: QuEST_cpu.c:4005
void statevec_collapseToKnownProbOutcome(Qureg qureg, int measureQubit, int outcome, qreal outcomeProb)
void shiftSubregIndices(int *allInds, int *numIndsPerReg, int numRegs, int shift)
Definition: QuEST_common.c:161
long long int getControlFlipMask(int *controlQubits, int *controlState, int numControlQubits)
Definition: QuEST_common.c:60
qreal statevec_getImagAmp(Qureg qureg, long long int index)
Represents a system of qubits.
Definition: QuEST.h:322
int statevec_measureWithStats(Qureg qureg, int measureQubit, qreal *outcomeProb)
Definition: QuEST_common.c:364
qreal ** imag
Definition: QuEST.h:190
void getZYZRotAnglesFromComplexPair(Complex alpha, Complex beta, qreal *rz2, qreal *ry, qreal *rz1)
maps U(alpha, beta) to Rz(rz2) Ry(ry) Rz(rz1)
Definition: QuEST_common.c:130
void statevec_controlledRotateX(Qureg qureg, int controlQubit, int targetQubit, qreal angle)
Definition: QuEST_common.c:346
ComplexArray stateVec
Computational state amplitudes - a subset thereof in the MPI version.
Definition: QuEST.h:341
qreal real[2][2]
Definition: QuEST.h:139
int isDensityMatrix
Whether this instance is a density-state representation.
Definition: QuEST.h:325
void statevec_hadamard(Qureg qureg, int targetQubit)
void statevec_rotateY(Qureg qureg, int rotQubit, qreal angle)
Definition: QuEST_common.c:302
int numQubits
Definition: QuEST.h:188
void statevec_multiControlledUnitary(Qureg qureg, long long int ctrlQubitsMask, long long int ctrlFlipMask, int targetQubit, ComplexMatrix2 u)
qreal getVectorMagnitude(Vector vec)
Definition: QuEST_common.c:79
void densmatr_applyTwoQubitKrausSuperoperator(Qureg qureg, int target1, int target2, ComplexMatrixN superOp)
Definition: QuEST_common.c:626
void statevec_phaseShiftByTerm(Qureg qureg, int targetQubit, Complex term)
Definition: QuEST_cpu.c:3185
int numQubits
The number of qubits informing the Hilbert dimension of the Hamiltonian.
Definition: QuEST.h:287
void statevec_sqrtSwapGateConj(Qureg qureg, int qb1, int qb2)
Definition: QuEST_common.c:400
void statevec_controlledRotateZ(Qureg qureg, int controlQubit, int targetQubit, qreal angle)
Definition: QuEST_common.c:358
int numQubitsRepresented
The number of qubits represented in either the state-vector or density matrix.
Definition: QuEST.h:327
void statevec_pauliX(Qureg qureg, int targetQubit)
@ GATE_SWAP
Definition: QuEST_qasm.h:33
qreal real
Definition: QuEST.h:105
void densmatr_applyMultiQubitKrausSuperoperator(Qureg qureg, int *targets, int numTargets, ComplexMatrixN superOp)
Definition: QuEST_common.c:634
void densmatr_mixPauli(Qureg qureg, int qubit, qreal probX, qreal probY, qreal probZ)
Definition: QuEST_common.c:743
qreal imag
Definition: QuEST.h:106
void ensureIndsIncrease(int *ind1, int *ind2)
Definition: QuEST_common.c:70
Represents one complex number.
Definition: QuEST.h:103
void statevec_controlledRotateY(Qureg qureg, int controlQubit, int targetQubit, qreal angle)
Definition: QuEST_common.c:352
void statevec_controlledRotateAroundAxisConj(Qureg qureg, int controlQubit, int targetQubit, qreal angle, Vector axis)
Definition: QuEST_common.c:337
#define macro_populateKrausOperator(superOp, ops, numOps, opDim)
Definition: QuEST_common.c:585
void qasm_recordGate(Qureg qureg, TargetGate gate, int targetQubit)
Definition: QuEST_qasm.c:179
qreal densmatr_calcProbOfOutcome(Qureg qureg, int measureQubit, int outcome)
void statevec_rotateZ(Qureg qureg, int rotQubit, qreal angle)
Definition: QuEST_common.c:308
#define M_PI
Definition: QuEST_common.c:41
void applySymmetrizedTrotterCircuit(Qureg qureg, PauliHamil hamil, qreal time, int order)
Definition: QuEST_common.c:820
void qasm_recordNamedPhaseFunc(Qureg qureg, int *qubits, int *numQubitsPerReg, int numRegs, enum bitEncoding encoding, enum phaseFunc funcName, qreal *params, int numParams, long long int *overrideInds, qreal *overridePhases, int numOverrides)
Definition: QuEST_qasm.c:726
void statevec_controlledCompactUnitary(Qureg qureg, int controlQubit, int targetQubit, Complex alpha, Complex beta)
qreal statevec_getRealAmp(Qureg qureg, long long int index)
int densmatr_measureWithStats(Qureg qureg, int measureQubit, qreal *outcomeProb)
Definition: QuEST_common.c:372
Represents a 2x2 matrix of complex numbers.
Definition: QuEST.h:137
void getComplexPairAndPhaseFromUnitary(ComplexMatrix2 u, Complex *alpha, Complex *beta, qreal *globalPhase)
maps U(r0c0, r0c1, r1c0, r1c1) to exp(i globalPhase) U(alpha, beta)
Definition: QuEST_common.c:142