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)
 
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 reportQuESTEnv (QuESTEnv env)
 Report information about the QuEST environment. More...
 
void seedQuESTDefault (void)
 Seed the Mersenne Twister used for random number generation in the QuEST environment with an example defualt seed. More...
 
Complex statevec_calcExpecDiagonalOp (Qureg qureg, DiagonalOp op)
 
Complex statevec_calcInnerProduct (Qureg bra, Qureg ket)
 
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_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 325 of file QuEST_cpu_local.c.

325  {
326 
327  // we must preload qureg.pairStateVec with the elements of op.
328  // instead of needless cloning, we'll just temporarily swap the pointers
329  qreal* rePtr = qureg.pairStateVec.real;
330  qreal* imPtr = qureg.pairStateVec.imag;
331  qureg.pairStateVec.real = op.real;
332  qureg.pairStateVec.imag = op.imag;
333 
335 
336  qureg.pairStateVec.real = rePtr;
337  qureg.pairStateVec.imag = imPtr;
338 }

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

◆ densmatr_calcExpecDiagonalOp()

Complex densmatr_calcExpecDiagonalOp ( Qureg  qureg,
DiagonalOp  op 
)

Definition at line 345 of file QuEST_cpu_local.c.

345  {
346 
347  return densmatr_calcExpecDiagonalOpLocal(qureg, op);
348 }

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_calcProbOfOutcome()

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

Definition at line 287 of file QuEST_cpu_local.c.

287  {
288 
289  qreal outcomeProb = densmatr_findProbabilityOfZeroLocal(qureg, measureQubit);
290  if (outcome == 1)
291  outcomeProb = 1.0 - outcomeProb;
292  return outcomeProb;
293 }

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 340 of file QuEST_cpu_local.c.

340  {
341 
342  return statevec_calcExpecDiagonalOpLocal(qureg, op);
343 }

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_calcProbOfOutcome()

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

Definition at line 279 of file QuEST_cpu_local.c.

280 {
281  qreal stateProb=0;
282  stateProb = statevec_findProbabilityOfZeroLocal(qureg, measureQubit);
283  if (outcome==1) stateProb = 1.0 - stateProb;
284  return stateProb;
285 }

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 295 of file QuEST_cpu_local.c.

296 {
297  statevec_collapseToKnownProbOutcomeLocal(qureg, measureQubit, outcome, stateProb);
298 }

References statevec_collapseToKnownProbOutcomeLocal().

◆ statevec_compactUnitary()

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

Definition at line 215 of file QuEST_cpu_local.c.

216 {
217  statevec_compactUnitaryLocal(qureg, targetQubit, alpha, beta);
218 }

References statevec_compactUnitaryLocal().

◆ statevec_controlledCompactUnitary()

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

Definition at line 225 of file QuEST_cpu_local.c.

226 {
227  statevec_controlledCompactUnitaryLocal(qureg, controlQubit, targetQubit, alpha, beta);
228 }

References statevec_controlledCompactUnitaryLocal().

◆ statevec_controlledNot()

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

Definition at line 274 of file QuEST_cpu_local.c.

275 {
276  statevec_controlledNotLocal(qureg, controlQubit, targetQubit);
277 }

References statevec_controlledNotLocal().

◆ statevec_controlledPauliY()

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

Definition at line 257 of file QuEST_cpu_local.c.

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

References statevec_controlledPauliYLocal().

◆ statevec_controlledPauliYConj()

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

Definition at line 263 of file QuEST_cpu_local.c.

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

References statevec_controlledPauliYLocal().

◆ statevec_controlledUnitary()

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

Definition at line 230 of file QuEST_cpu_local.c.

231 {
232  statevec_controlledUnitaryLocal(qureg, controlQubit, targetQubit, u);
233 }

References statevec_controlledUnitaryLocal().

◆ statevec_getImagAmp()

qreal statevec_getImagAmp ( Qureg  qureg,
long long int  index 
)

Definition at line 211 of file QuEST_cpu_local.c.

211  {
212  return qureg.stateVec.imag[index];
213 }

References Qureg::stateVec.

◆ statevec_getRealAmp()

qreal statevec_getRealAmp ( Qureg  qureg,
long long int  index 
)

Definition at line 207 of file QuEST_cpu_local.c.

207  {
208  return qureg.stateVec.real[index];
209 }

References Qureg::stateVec.

◆ statevec_hadamard()

void statevec_hadamard ( Qureg  qureg,
int  targetQubit 
)

Definition at line 269 of file QuEST_cpu_local.c.

270 {
271  statevec_hadamardLocal(qureg, targetQubit);
272 }

References statevec_hadamardLocal().

◆ 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 315 of file QuEST_cpu_local.c.

316 {
317  statevec_multiControlledMultiQubitUnitaryLocal(qureg, ctrlMask, targs, numTargs, u);
318 }

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. @TODO: the double swap (q1,q2 to 0,1) may be possible simultaneously by a bespoke swap routine.

Definition at line 310 of file QuEST_cpu_local.c.

311 {
312  statevec_multiControlledTwoQubitUnitaryLocal(qureg, ctrlMask, q1, q2, u);
313 }

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 235 of file QuEST_cpu_local.c.

236 {
237  statevec_multiControlledUnitaryLocal(qureg, targetQubit, ctrlQubitsMask, ctrlFlipMask, u);
238 }

References statevec_multiControlledUnitaryLocal().

◆ statevec_pauliX()

void statevec_pauliX ( Qureg  qureg,
int  targetQubit 
)

Definition at line 240 of file QuEST_cpu_local.c.

241 {
242  statevec_pauliXLocal(qureg, targetQubit);
243 }

References statevec_pauliXLocal().

◆ statevec_pauliY()

void statevec_pauliY ( Qureg  qureg,
int  targetQubit 
)

Definition at line 245 of file QuEST_cpu_local.c.

246 {
247  int conjFac = 1;
248  statevec_pauliYLocal(qureg, targetQubit, conjFac);
249 }

References statevec_pauliYLocal().

◆ statevec_pauliYConj()

void statevec_pauliYConj ( Qureg  qureg,
int  targetQubit 
)

Definition at line 251 of file QuEST_cpu_local.c.

252 {
253  int conjFac = -1;
254  statevec_pauliYLocal(qureg, targetQubit, conjFac);
255 }

References statevec_pauliYLocal().

◆ statevec_swapQubitAmps()

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

Definition at line 320 of file QuEST_cpu_local.c.

321 {
322  statevec_swapQubitAmpsLocal(qureg, qb1, qb2);
323 }

References statevec_swapQubitAmpsLocal().

◆ statevec_unitary()

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

Definition at line 220 of file QuEST_cpu_local.c.

221 {
222  statevec_unitaryLocal(qureg, targetQubit, u);
223 }

References statevec_unitaryLocal().

void statevec_controlledNotLocal(Qureg qureg, int controlQubit, int targetQubit)
Definition: QuEST_cpu.c:2584
void statevec_pauliYLocal(Qureg qureg, int targetQubit, int conjFac)
Definition: QuEST_cpu.c:2682
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:923
ComplexArray pairStateVec
Temporary storage for a chunk of the state vector received from another process in the MPI version.
Definition: QuEST.h:224
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:3536
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:3380
void statevec_compactUnitaryLocal(Qureg qureg, int targetQubit, Complex alpha, Complex beta)
Definition: QuEST_cpu.c:1688
qreal densmatr_calcPurityLocal(Qureg qureg)
Definition: QuEST_cpu.c:861
Complex statevec_calcExpecDiagonalOpLocal(Qureg qureg, DiagonalOp op)
Definition: QuEST_cpu.c:3738
void statevec_multiControlledMultiQubitUnitaryLocal(Qureg qureg, long long int ctrlMask, int *targs, int numTargs, ComplexMatrixN u)
Definition: QuEST_cpu.c:1846
void statevec_multiControlledTwoQubitUnitaryLocal(Qureg qureg, long long int ctrlMask, int q1, int q2, ComplexMatrix4 u)
Definition: QuEST_cpu.c:1747
#define qreal
qreal * imag
The imaginary values of the 2^numQubits complex elements.
Definition: QuEST.h:191
long long int numAmpsPerChunk
Number of probability amplitudes held in stateVec by this process In the non-MPI version,...
Definition: QuEST.h:213
void densmatr_mixTwoQubitDepolarisingLocal(Qureg qureg, int qubit1, int qubit2, qreal delta, qreal gamma)
Definition: QuEST_cpu.c:387
void statevec_controlledPauliYLocal(Qureg qureg, int controlQubit, int targetQubit, int conjFac)
Definition: QuEST_cpu.c:2776
void densmatr_mixTwoQubitDephasing(Qureg qureg, int qubit1, int qubit2, qreal dephase)
Definition: QuEST_cpu.c:84
void statevec_multiControlledUnitaryLocal(Qureg qureg, int targetQubit, long long int ctrlQubitsMask, long long int ctrlFlipMask, ComplexMatrix2 u)
Definition: QuEST_cpu.c:2173
Complex densmatr_calcExpecDiagonalOpLocal(Qureg qureg, DiagonalOp op)
Definition: QuEST_cpu.c:3781
void statevec_pauliXLocal(Qureg qureg, int targetQubit)
Definition: QuEST_cpu.c:2498
void densmatr_initPureStateLocal(Qureg targetQureg, Qureg copyQureg)
Definition: QuEST_cpu.c:1184
void densmatr_mixDampingLocal(Qureg qureg, int targetQubit, qreal damping)
Definition: QuEST_cpu.c:174
qreal densmatr_calcInnerProductLocal(Qureg a, Qureg b)
computes Tr(conjTrans(a) b) = sum of (a_ij^* b_ij)
Definition: QuEST_cpu.c:958
ComplexArray stateVec
Computational state amplitudes - a subset thereof in the MPI version.
Definition: QuEST.h:222
void densmatr_mixDepolarisingLocal(Qureg qureg, int targetQubit, qreal depolLevel)
Definition: QuEST_cpu.c:125
void statevec_hadamardLocal(Qureg qureg, int targetQubit)
Definition: QuEST_cpu.c:2872
void densmatr_applyDiagonalOpLocal(Qureg qureg, DiagonalOp op)
Definition: QuEST_cpu.c:3696
void statevec_controlledUnitaryLocal(Qureg qureg, int controlQubit, int targetQubit, ComplexMatrix2 u)
Definition: QuEST_cpu.c:2241
int numQubitsRepresented
The number of qubits represented in either the state-vector or density matrix.
Definition: QuEST.h:208
qreal * real
The real values of the 2^numQubits complex elements.
Definition: QuEST.h:189
void statevec_unitaryLocal(Qureg qureg, int targetQubit, ComplexMatrix2 u)
Definition: QuEST_cpu.c:1932
qreal densmatr_findProbabilityOfZeroLocal(Qureg qureg, int measureQubit)
Definition: QuEST_cpu.c:3151
void statevec_controlledCompactUnitaryLocal(Qureg qureg, int controlQubit, int targetQubit, Complex alpha, Complex beta)
Definition: QuEST_cpu.c:2101
Complex statevec_calcInnerProductLocal(Qureg bra, Qureg ket)
Definition: QuEST_cpu.c:1071
qreal densmatr_calcFidelityLocal(Qureg qureg, Qureg pureState)
computes a few dens-columns-worth of (vec^*T) dens * vec
Definition: QuEST_cpu.c:990
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:3206