QuEST_cpu_local.c File Reference
#include "QuEST.h"
#include "QuEST_internal.h"
#include "QuEST_precision.h"
#include "mt19937ar.h"
#include "QuEST_cpu_internal.h"
#include <stdlib.h>
#include <stdio.h>
#include <math.h>
#include <time.h>
#include <sys/types.h>

Go to the source code of this file.

Functions

QuESTEnv createQuESTEnv (void)
 Create the QuEST execution environment. More...
 
void densmatr_applyDiagonalOp (Qureg qureg, DiagonalOp op)
 
Complex densmatr_calcExpecDiagonalOp (Qureg qureg, DiagonalOp op)
 
qreal densmatr_calcFidelity (Qureg qureg, Qureg pureState)
 
qreal densmatr_calcHilbertSchmidtDistance (Qureg a, Qureg b)
 
qreal densmatr_calcInnerProduct (Qureg a, Qureg b)
 
void densmatr_calcProbOfAllOutcomes (qreal *retProbs, Qureg qureg, int *qubits, int numQubits)
 
qreal densmatr_calcProbOfOutcome (Qureg qureg, int measureQubit, int outcome)
 
qreal densmatr_calcPurity (Qureg qureg)
 
qreal densmatr_calcTotalProb (Qureg qureg)
 
void densmatr_initPureState (Qureg qureg, Qureg pureState)
 
void densmatr_mixDamping (Qureg qureg, int targetQubit, qreal damping)
 
void densmatr_mixDepolarising (Qureg qureg, int targetQubit, qreal depolLevel)
 
void densmatr_mixTwoQubitDepolarising (Qureg qureg, int qubit1, int qubit2, qreal depolLevel)
 
void destroyQuESTEnv (QuESTEnv env)
 Destroy the QuEST environment. More...
 
void getEnvironmentString (QuESTEnv env, char str[200])
 Sets str to a string containing information about the runtime environment, including whether simulation is using CUDA (for GPU), OpenMP (for multithreading) and/or MPI (for distribution). More...
 
void reportQuESTEnv (QuESTEnv env)
 Report information about the QuEST environment. More...
 
void seedQuEST (QuESTEnv *env, unsigned long int *seedArray, int numSeeds)
 Seeds the random number generator with a custom array of key(s), overriding the default keys. More...
 
Complex statevec_calcExpecDiagonalOp (Qureg qureg, DiagonalOp op)
 
Complex statevec_calcInnerProduct (Qureg bra, Qureg ket)
 
void statevec_calcProbOfAllOutcomes (qreal *retProbs, Qureg qureg, int *qubits, int numQubits)
 
qreal statevec_calcProbOfOutcome (Qureg qureg, int measureQubit, int outcome)
 
qreal statevec_calcTotalProb (Qureg qureg)
 
void statevec_collapseToKnownProbOutcome (Qureg qureg, int measureQubit, int outcome, qreal stateProb)
 
void statevec_compactUnitary (Qureg qureg, int targetQubit, Complex alpha, Complex beta)
 
void statevec_controlledCompactUnitary (Qureg qureg, int controlQubit, int targetQubit, Complex alpha, Complex beta)
 
void statevec_controlledNot (Qureg qureg, int controlQubit, int targetQubit)
 
void statevec_controlledPauliY (Qureg qureg, int controlQubit, int targetQubit)
 
void statevec_controlledPauliYConj (Qureg qureg, int controlQubit, int targetQubit)
 
void statevec_controlledUnitary (Qureg qureg, int controlQubit, int targetQubit, ComplexMatrix2 u)
 
qreal statevec_getImagAmp (Qureg qureg, long long int index)
 
qreal statevec_getRealAmp (Qureg qureg, long long int index)
 
void statevec_hadamard (Qureg qureg, int targetQubit)
 
void statevec_multiControlledMultiQubitNot (Qureg qureg, int ctrlMask, int targMask)
 
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 already fit in the node, it operates the unitary direct. More...
 
void statevec_multiControlledTwoQubitUnitary (Qureg qureg, long long int ctrlMask, int q1, int q2, ComplexMatrix4 u)
 This calls swapQubitAmps only when it would involve a distributed communication; if the qubit chunks already fit in the node, it operates the unitary direct. More...
 
void statevec_multiControlledUnitary (Qureg qureg, long long int ctrlQubitsMask, long long int ctrlFlipMask, int targetQubit, ComplexMatrix2 u)
 
void statevec_pauliX (Qureg qureg, int targetQubit)
 
void statevec_pauliY (Qureg qureg, int targetQubit)
 
void statevec_pauliYConj (Qureg qureg, int targetQubit)
 
void statevec_swapQubitAmps (Qureg qureg, int qb1, int qb2)
 
void statevec_unitary (Qureg qureg, int targetQubit, ComplexMatrix2 u)
 
void syncQuESTEnv (QuESTEnv env)
 Guarantees that all code up to the given point has been executed on all nodes (if running in distributed mode) More...
 
int syncQuESTSuccess (int successCode)
 Performs a logical AND on all successCodes held by all processes. More...
 

Detailed Description

An implementation of the pure backend in ../QuEST_ops_pure.h for a local (non-MPI, non-GPU) environment. Mostly pure-state wrappers for the local/distributed functions implemented in QuEST_cpu

Author
Ania Brown
Tyson Jones
Balint Koczor

Definition in file QuEST_cpu_local.c.

Function Documentation

◆ densmatr_applyDiagonalOp()

void densmatr_applyDiagonalOp ( Qureg  qureg,
DiagonalOp  op 
)

Definition at line 358 of file QuEST_cpu_local.c.

358  {
359 
360  // we must preload qureg.pairStateVec with the elements of op.
361  // instead of needless cloning, we'll just temporarily swap the pointers
362  qreal* rePtr = qureg.pairStateVec.real;
363  qreal* imPtr = qureg.pairStateVec.imag;
364  qureg.pairStateVec.real = op.real;
365  qureg.pairStateVec.imag = op.imag;
366 
368 
369  qureg.pairStateVec.real = rePtr;
370  qureg.pairStateVec.imag = imPtr;
371 }

References densmatr_applyDiagonalOpLocal(), DiagonalOp::imag, Qureg::pairStateVec, qreal, and DiagonalOp::real.

◆ densmatr_calcExpecDiagonalOp()

Complex densmatr_calcExpecDiagonalOp ( Qureg  qureg,
DiagonalOp  op 
)

Definition at line 378 of file QuEST_cpu_local.c.

378  {
379 
380  return densmatr_calcExpecDiagonalOpLocal(qureg, op);
381 }

References densmatr_calcExpecDiagonalOpLocal().

◆ densmatr_calcFidelity()

qreal densmatr_calcFidelity ( Qureg  qureg,
Qureg  pureState 
)

Definition at line 75 of file QuEST_cpu_local.c.

75  {
76 
77  // save pointers to qureg's pair state
78  qreal* quregPairRePtr = qureg.pairStateVec.real;
79  qreal* quregPairImPtr = qureg.pairStateVec.imag;
80 
81  // populate qureg pair state with pure state (by repointing)
82  qureg.pairStateVec.real = pureState.stateVec.real;
83  qureg.pairStateVec.imag = pureState.stateVec.imag;
84 
85  // calculate fidelity using pairState
86  qreal fid = densmatr_calcFidelityLocal(qureg, pureState);
87 
88  // restore pointers
89  qureg.pairStateVec.real = quregPairRePtr;
90  qureg.pairStateVec.imag = quregPairImPtr;
91 
92  return fid;
93 }

References densmatr_calcFidelityLocal(), Qureg::pairStateVec, qreal, and Qureg::stateVec.

◆ densmatr_calcHilbertSchmidtDistance()

qreal densmatr_calcHilbertSchmidtDistance ( Qureg  a,
Qureg  b 
)

Definition at line 62 of file QuEST_cpu_local.c.

62  {
63 
65  qreal dist = sqrt(distSquared);
66  return dist;
67 }

References densmatr_calcHilbertSchmidtDistanceSquaredLocal(), and qreal.

◆ densmatr_calcInnerProduct()

qreal densmatr_calcInnerProduct ( Qureg  a,
Qureg  b 
)

Definition at line 69 of file QuEST_cpu_local.c.

69  {
70 
72  return scalar;
73 }

References densmatr_calcInnerProductLocal(), and qreal.

◆ densmatr_calcProbOfAllOutcomes()

void densmatr_calcProbOfAllOutcomes ( qreal retProbs,
Qureg  qureg,
int *  qubits,
int  numQubits 
)

Definition at line 317 of file QuEST_cpu_local.c.

317  {
318 
319  densmatr_calcProbOfAllOutcomesLocal(retProbs, qureg, qubits, numQubits);
320 }

References densmatr_calcProbOfAllOutcomesLocal().

◆ densmatr_calcProbOfOutcome()

qreal densmatr_calcProbOfOutcome ( Qureg  qureg,
int  measureQubit,
int  outcome 
)

Definition at line 304 of file QuEST_cpu_local.c.

304  {
305 
306  qreal outcomeProb = densmatr_findProbabilityOfZeroLocal(qureg, measureQubit);
307  if (outcome == 1)
308  outcomeProb = 1.0 - outcomeProb;
309  return outcomeProb;
310 }

References densmatr_findProbabilityOfZeroLocal(), and qreal.

◆ densmatr_calcPurity()

qreal densmatr_calcPurity ( Qureg  qureg)

Definition at line 57 of file QuEST_cpu_local.c.

57  {
58 
59  return densmatr_calcPurityLocal(qureg);
60 }

References densmatr_calcPurityLocal().

◆ densmatr_calcTotalProb()

qreal densmatr_calcTotalProb ( Qureg  qureg)

Definition at line 118 of file QuEST_cpu_local.c.

118  {
119 
120  // computes the trace using Kahan summation
121  qreal pTotal=0;
122  qreal y, t, c;
123  c = 0;
124 
125  long long int numCols = 1LL << qureg.numQubitsRepresented;
126  long long diagIndex;
127 
128  for (int col=0; col< numCols; col++) {
129  diagIndex = col*(numCols + 1);
130  y = qureg.stateVec.real[diagIndex] - c;
131  t = pTotal + y;
132  c = ( t - pTotal ) - y; // brackets are important
133  pTotal = t;
134  }
135 
136  // does not check imaginary component, by design
137 
138  return pTotal;
139 }

References Qureg::numQubitsRepresented, qreal, and Qureg::stateVec.

◆ densmatr_initPureState()

void densmatr_initPureState ( Qureg  qureg,
Qureg  pureState 
)

Definition at line 95 of file QuEST_cpu_local.c.

95  {
96 
97  // save pointers to qureg's pair state
98  qreal* quregPairRePtr = qureg.pairStateVec.real;
99  qreal* quregPairImPtr = qureg.pairStateVec.imag;
100 
101  // populate qureg pair state with pure state (by repointing)
102  qureg.pairStateVec.real = pureState.stateVec.real;
103  qureg.pairStateVec.imag = pureState.stateVec.imag;
104 
105  // populate density matrix via it's pairState
106  densmatr_initPureStateLocal(qureg, pureState);
107 
108  // restore pointers
109  qureg.pairStateVec.real = quregPairRePtr;
110  qureg.pairStateVec.imag = quregPairImPtr;
111 }

References densmatr_initPureStateLocal(), Qureg::pairStateVec, qreal, and Qureg::stateVec.

◆ densmatr_mixDamping()

void densmatr_mixDamping ( Qureg  qureg,
int  targetQubit,
qreal  damping 
)

Definition at line 37 of file QuEST_cpu_local.c.

37  {
38  if (damping == 0)
39  return;
40 
41  densmatr_mixDampingLocal(qureg, targetQubit, damping);
42 }

References densmatr_mixDampingLocal().

◆ densmatr_mixDepolarising()

void densmatr_mixDepolarising ( Qureg  qureg,
int  targetQubit,
qreal  depolLevel 
)

Definition at line 30 of file QuEST_cpu_local.c.

30  {
31  if (depolLevel == 0)
32  return;
33 
34  densmatr_mixDepolarisingLocal(qureg, targetQubit, depolLevel);
35 }

References densmatr_mixDepolarisingLocal().

◆ densmatr_mixTwoQubitDepolarising()

void densmatr_mixTwoQubitDepolarising ( Qureg  qureg,
int  qubit1,
int  qubit2,
qreal  depolLevel 
)

Definition at line 44 of file QuEST_cpu_local.c.

44  {
45  if (depolLevel == 0)
46  return;
47  qreal eta = 2/depolLevel;
48  qreal delta = eta - 1 - sqrt( (eta-1)*(eta-1) - 1 );
49  qreal gamma = 1+delta;
50  // TODO -- test delta too small
51 
52  gamma = 1/(gamma*gamma*gamma);
53  densmatr_mixTwoQubitDephasing(qureg, qubit1, qubit2, depolLevel);
54  densmatr_mixTwoQubitDepolarisingLocal(qureg, qubit1, qubit2, delta, gamma);
55 }

References densmatr_mixTwoQubitDephasing(), densmatr_mixTwoQubitDepolarisingLocal(), and qreal.

◆ statevec_calcExpecDiagonalOp()

Complex statevec_calcExpecDiagonalOp ( Qureg  qureg,
DiagonalOp  op 
)

Definition at line 373 of file QuEST_cpu_local.c.

373  {
374 
375  return statevec_calcExpecDiagonalOpLocal(qureg, op);
376 }

References statevec_calcExpecDiagonalOpLocal().

◆ statevec_calcInnerProduct()

Complex statevec_calcInnerProduct ( Qureg  bra,
Qureg  ket 
)

Definition at line 113 of file QuEST_cpu_local.c.

113  {
114 
115  return statevec_calcInnerProductLocal(bra, ket);
116 }

References statevec_calcInnerProductLocal().

◆ statevec_calcProbOfAllOutcomes()

void statevec_calcProbOfAllOutcomes ( qreal retProbs,
Qureg  qureg,
int *  qubits,
int  numQubits 
)

Definition at line 312 of file QuEST_cpu_local.c.

312  {
313 
314  statevec_calcProbOfAllOutcomesLocal(retProbs, qureg, qubits, numQubits);
315 }

References statevec_calcProbOfAllOutcomesLocal().

◆ statevec_calcProbOfOutcome()

qreal statevec_calcProbOfOutcome ( Qureg  qureg,
int  measureQubit,
int  outcome 
)

Definition at line 296 of file QuEST_cpu_local.c.

297 {
298  qreal outcomeProb = statevec_findProbabilityOfZeroLocal(qureg, measureQubit);
299  if (outcome==1)
300  outcomeProb = 1.0 - outcomeProb;
301  return outcomeProb;
302 }

References qreal, and statevec_findProbabilityOfZeroLocal().

◆ statevec_calcTotalProb()

qreal statevec_calcTotalProb ( Qureg  qureg)

Definition at line 141 of file QuEST_cpu_local.c.

141  {
142  // implemented using Kahan summation for greater accuracy at a slight floating
143  // point operation overhead. For more details see https://en.wikipedia.org/wiki/Kahan_summation_algorithm
144  qreal pTotal=0;
145  qreal y, t, c;
146  long long int index;
147  long long int numAmpsPerRank = qureg.numAmpsPerChunk;
148  c = 0.0;
149  for (index=0; index<numAmpsPerRank; index++){
150  // Perform pTotal+=qureg.stateVec.real[index]*qureg.stateVec.real[index]; by Kahan
151 
152  y = qureg.stateVec.real[index]*qureg.stateVec.real[index] - c;
153  t = pTotal + y;
154  // Don't change the bracketing on the following line
155  c = ( t - pTotal ) - y;
156  pTotal = t;
157 
158  // Perform pTotal+=qureg.stateVec.imag[index]*qureg.stateVec.imag[index]; by Kahan
159 
160  y = qureg.stateVec.imag[index]*qureg.stateVec.imag[index] - c;
161  t = pTotal + y;
162  // Don't change the bracketing on the following line
163  c = ( t - pTotal ) - y;
164  pTotal = t;
165  }
166  return pTotal;
167 }

References Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.

◆ statevec_collapseToKnownProbOutcome()

void statevec_collapseToKnownProbOutcome ( Qureg  qureg,
int  measureQubit,
int  outcome,
qreal  stateProb 
)

Definition at line 322 of file QuEST_cpu_local.c.

323 {
324  statevec_collapseToKnownProbOutcomeLocal(qureg, measureQubit, outcome, stateProb);
325 }

References statevec_collapseToKnownProbOutcomeLocal().

◆ statevec_compactUnitary()

void statevec_compactUnitary ( Qureg  qureg,
int  targetQubit,
Complex  alpha,
Complex  beta 
)

Definition at line 227 of file QuEST_cpu_local.c.

228 {
229  statevec_compactUnitaryLocal(qureg, targetQubit, alpha, beta);
230 }

References statevec_compactUnitaryLocal().

◆ statevec_controlledCompactUnitary()

void statevec_controlledCompactUnitary ( Qureg  qureg,
int  controlQubit,
int  targetQubit,
Complex  alpha,
Complex  beta 
)

Definition at line 237 of file QuEST_cpu_local.c.

238 {
239  statevec_controlledCompactUnitaryLocal(qureg, controlQubit, targetQubit, alpha, beta);
240 }

References statevec_controlledCompactUnitaryLocal().

◆ statevec_controlledNot()

void statevec_controlledNot ( Qureg  qureg,
int  controlQubit,
int  targetQubit 
)

Definition at line 286 of file QuEST_cpu_local.c.

287 {
288  statevec_controlledNotLocal(qureg, controlQubit, targetQubit);
289 }

References statevec_controlledNotLocal().

◆ statevec_controlledPauliY()

void statevec_controlledPauliY ( Qureg  qureg,
int  controlQubit,
int  targetQubit 
)

Definition at line 269 of file QuEST_cpu_local.c.

270 {
271  int conjFac = 1;
272  statevec_controlledPauliYLocal(qureg, controlQubit, targetQubit, conjFac);
273 }

References statevec_controlledPauliYLocal().

◆ statevec_controlledPauliYConj()

void statevec_controlledPauliYConj ( Qureg  qureg,
int  controlQubit,
int  targetQubit 
)

Definition at line 275 of file QuEST_cpu_local.c.

276 {
277  int conjFac = -1;
278  statevec_controlledPauliYLocal(qureg, controlQubit, targetQubit, conjFac);
279 }

References statevec_controlledPauliYLocal().

◆ statevec_controlledUnitary()

void statevec_controlledUnitary ( Qureg  qureg,
int  controlQubit,
int  targetQubit,
ComplexMatrix2  u 
)

Definition at line 242 of file QuEST_cpu_local.c.

243 {
244  statevec_controlledUnitaryLocal(qureg, controlQubit, targetQubit, u);
245 }

References statevec_controlledUnitaryLocal().

◆ statevec_getImagAmp()

qreal statevec_getImagAmp ( Qureg  qureg,
long long int  index 
)

Definition at line 223 of file QuEST_cpu_local.c.

223  {
224  return qureg.stateVec.imag[index];
225 }

References Qureg::stateVec.

◆ statevec_getRealAmp()

qreal statevec_getRealAmp ( Qureg  qureg,
long long int  index 
)

Definition at line 219 of file QuEST_cpu_local.c.

219  {
220  return qureg.stateVec.real[index];
221 }

References Qureg::stateVec.

◆ statevec_hadamard()

void statevec_hadamard ( Qureg  qureg,
int  targetQubit 
)

Definition at line 281 of file QuEST_cpu_local.c.

282 {
283  statevec_hadamardLocal(qureg, targetQubit);
284 }

References statevec_hadamardLocal().

◆ statevec_multiControlledMultiQubitNot()

void statevec_multiControlledMultiQubitNot ( Qureg  qureg,
int  ctrlMask,
int  targMask 
)

Definition at line 291 of file QuEST_cpu_local.c.

291  {
292 
293  statevec_multiControlledMultiQubitNotLocal(qureg, ctrlMask, targMask);
294 }

References statevec_multiControlledMultiQubitNotLocal().

◆ statevec_multiControlledMultiQubitUnitary()

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 already fit in the node, it operates the unitary direct.

It is already gauranteed here that all target qubits can fit on each node (this is validated in the front-end)

Todo:
refactor so that the 'swap back' isn't performed; instead the qubit locations are updated.

Definition at line 348 of file QuEST_cpu_local.c.

349 {
350  statevec_multiControlledMultiQubitUnitaryLocal(qureg, ctrlMask, targs, numTargs, u);
351 }

References statevec_multiControlledMultiQubitUnitaryLocal().

◆ statevec_multiControlledTwoQubitUnitary()

void statevec_multiControlledTwoQubitUnitary ( Qureg  qureg,
long long int  ctrlMask,
int  q1,
int  q2,
ComplexMatrix4  u 
)

This calls swapQubitAmps only when it would involve a distributed communication; if the qubit chunks already fit in the node, it operates the unitary direct.

Note the order of q1 and q2 in the call to twoQubitUnitaryLocal is important.

Todo:

refactor so that the 'swap back' isn't performed; instead the qubit locations are updated.

the double swap (q1,q2 to 0,1) may be possible simultaneously by a bespoke swap routine.

Definition at line 343 of file QuEST_cpu_local.c.

344 {
345  statevec_multiControlledTwoQubitUnitaryLocal(qureg, ctrlMask, q1, q2, u);
346 }

References statevec_multiControlledTwoQubitUnitaryLocal().

◆ statevec_multiControlledUnitary()

void statevec_multiControlledUnitary ( Qureg  qureg,
long long int  ctrlQubitsMask,
long long int  ctrlFlipMask,
int  targetQubit,
ComplexMatrix2  u 
)

Definition at line 247 of file QuEST_cpu_local.c.

248 {
249  statevec_multiControlledUnitaryLocal(qureg, targetQubit, ctrlQubitsMask, ctrlFlipMask, u);
250 }

References statevec_multiControlledUnitaryLocal().

◆ statevec_pauliX()

void statevec_pauliX ( Qureg  qureg,
int  targetQubit 
)

Definition at line 252 of file QuEST_cpu_local.c.

253 {
254  statevec_pauliXLocal(qureg, targetQubit);
255 }

References statevec_pauliXLocal().

◆ statevec_pauliY()

void statevec_pauliY ( Qureg  qureg,
int  targetQubit 
)

Definition at line 257 of file QuEST_cpu_local.c.

258 {
259  int conjFac = 1;
260  statevec_pauliYLocal(qureg, targetQubit, conjFac);
261 }

References statevec_pauliYLocal().

◆ statevec_pauliYConj()

void statevec_pauliYConj ( Qureg  qureg,
int  targetQubit 
)

Definition at line 263 of file QuEST_cpu_local.c.

264 {
265  int conjFac = -1;
266  statevec_pauliYLocal(qureg, targetQubit, conjFac);
267 }

References statevec_pauliYLocal().

◆ statevec_swapQubitAmps()

void statevec_swapQubitAmps ( Qureg  qureg,
int  qb1,
int  qb2 
)

Definition at line 353 of file QuEST_cpu_local.c.

354 {
355  statevec_swapQubitAmpsLocal(qureg, qb1, qb2);
356 }

References statevec_swapQubitAmpsLocal().

◆ statevec_unitary()

void statevec_unitary ( Qureg  qureg,
int  targetQubit,
ComplexMatrix2  u 
)

Definition at line 232 of file QuEST_cpu_local.c.

233 {
234  statevec_unitaryLocal(qureg, targetQubit, u);
235 }

References statevec_unitaryLocal().

void statevec_controlledNotLocal(Qureg qureg, int controlQubit, int targetQubit)
Definition: QuEST_cpu.c:2679
void statevec_pauliYLocal(Qureg qureg, int targetQubit, int conjFac)
Definition: QuEST_cpu.c:2887
qreal densmatr_calcHilbertSchmidtDistanceSquaredLocal(Qureg a, Qureg b)
computes Tr((a-b) conjTrans(a-b)) = sum of abs values of (a-b)
Definition: QuEST_cpu.c:934
void densmatr_calcProbOfAllOutcomesLocal(qreal *outcomeProbs, Qureg qureg, int *qubits, int numQubits)
Definition: QuEST_cpu.c:3616
ComplexArray pairStateVec
Temporary storage for a chunk of the state vector received from another process in the MPI version.
Definition: QuEST.h:343
void statevec_swapQubitAmpsLocal(Qureg qureg, int qb1, int qb2)
It is ensured that all amplitudes needing to be swapped are on this node.
Definition: QuEST_cpu.c:3922
void statevec_multiControlledMultiQubitNotLocal(Qureg qureg, int ctrlMask, int targMask)
Definition: QuEST_cpu.c:2778
void statevec_collapseToKnownProbOutcomeLocal(Qureg qureg, int measureQubit, int outcome, qreal totalProbability)
Update the state vector to be consistent with measuring measureQubit=0 if outcome=0 and measureQubit=...
Definition: QuEST_cpu.c:3767
void statevec_compactUnitaryLocal(Qureg qureg, int targetQubit, Complex alpha, Complex beta)
Definition: QuEST_cpu.c:1754
qreal densmatr_calcPurityLocal(Qureg qureg)
Definition: QuEST_cpu.c:872
Complex statevec_calcExpecDiagonalOpLocal(Qureg qureg, DiagonalOp op)
Definition: QuEST_cpu.c:4124
void statevec_multiControlledMultiQubitUnitaryLocal(Qureg qureg, long long int ctrlMask, int *targs, int numTargs, ComplexMatrixN u)
Definition: QuEST_cpu.c:1912
void statevec_multiControlledTwoQubitUnitaryLocal(Qureg qureg, long long int ctrlMask, int q1, int q2, ComplexMatrix4 u)
Definition: QuEST_cpu.c:1813
#define qreal
void statevec_calcProbOfAllOutcomesLocal(qreal *outcomeProbs, Qureg qureg, int *qubits, int numQubits)
Definition: QuEST_cpu.c:3549
qreal * imag
The imaginary values of the 2^numQubits complex elements.
Definition: QuEST.h:310
long long int numAmpsPerChunk
Number of probability amplitudes held in stateVec by this process In the non-MPI version,...
Definition: QuEST.h:332
void densmatr_mixTwoQubitDepolarisingLocal(Qureg qureg, int qubit1, int qubit2, qreal delta, qreal gamma)
Definition: QuEST_cpu.c:393
void statevec_controlledPauliYLocal(Qureg qureg, int controlQubit, int targetQubit, int conjFac)
Definition: QuEST_cpu.c:2982
void densmatr_mixTwoQubitDephasing(Qureg qureg, int qubit1, int qubit2, qreal dephase)
Definition: QuEST_cpu.c:90
void statevec_multiControlledUnitaryLocal(Qureg qureg, int targetQubit, long long int ctrlQubitsMask, long long int ctrlFlipMask, ComplexMatrix2 u)
Definition: QuEST_cpu.c:2268
Complex densmatr_calcExpecDiagonalOpLocal(Qureg qureg, DiagonalOp op)
Definition: QuEST_cpu.c:4167
void statevec_pauliXLocal(Qureg qureg, int targetQubit)
Definition: QuEST_cpu.c:2593
void densmatr_initPureStateLocal(Qureg targetQureg, Qureg copyQureg)
Definition: QuEST_cpu.c:1195
void densmatr_mixDampingLocal(Qureg qureg, int targetQubit, qreal damping)
Definition: QuEST_cpu.c:180
qreal densmatr_calcInnerProductLocal(Qureg a, Qureg b)
computes Tr(conjTrans(a) b) = sum of (a_ij^* b_ij)
Definition: QuEST_cpu.c:969
ComplexArray stateVec
Computational state amplitudes - a subset thereof in the MPI version.
Definition: QuEST.h:341
void densmatr_mixDepolarisingLocal(Qureg qureg, int targetQubit, qreal depolLevel)
Definition: QuEST_cpu.c:131
void statevec_hadamardLocal(Qureg qureg, int targetQubit)
Definition: QuEST_cpu.c:3078
void densmatr_applyDiagonalOpLocal(Qureg qureg, DiagonalOp op)
Definition: QuEST_cpu.c:4082
void statevec_controlledUnitaryLocal(Qureg qureg, int controlQubit, int targetQubit, ComplexMatrix2 u)
Definition: QuEST_cpu.c:2336
int numQubitsRepresented
The number of qubits represented in either the state-vector or density matrix.
Definition: QuEST.h:327
qreal * real
The real values of the 2^numQubits complex elements.
Definition: QuEST.h:308
void statevec_unitaryLocal(Qureg qureg, int targetQubit, ComplexMatrix2 u)
Definition: QuEST_cpu.c:2026
qreal densmatr_findProbabilityOfZeroLocal(Qureg qureg, int measureQubit)
Definition: QuEST_cpu.c:3402
void statevec_controlledCompactUnitaryLocal(Qureg qureg, int controlQubit, int targetQubit, Complex alpha, Complex beta)
Definition: QuEST_cpu.c:2196
Complex statevec_calcInnerProductLocal(Qureg bra, Qureg ket)
Definition: QuEST_cpu.c:1082
qreal densmatr_calcFidelityLocal(Qureg qureg, Qureg pureState)
computes a few dens-columns-worth of (vec^*T) dens * vec
Definition: QuEST_cpu.c:1001
qreal statevec_findProbabilityOfZeroLocal(Qureg qureg, int measureQubit)
Measure the total probability of a specified qubit being in the zero state across all amplitudes in t...
Definition: QuEST_cpu.c:3457