QuEST.h
1// Distributed under MIT licence. See https://github.com/QuEST-Kit/QuEST/blob/master/LICENCE.txt for details
2
32# ifndef QUEST_H
33# define QUEST_H
34
35# include "QuEST_precision.h"
36
37
38
39// ensure custatevecHandle_t is defined, even if no GPU
40# ifdef USE_CUQUANTUM
41# include <custatevec.h>
42typedef struct CuQuantumConfig {
43 cudaMemPool_t cuMemPool;
44 cudaStream_t cuStream;
45 custatevecHandle_t cuQuantumHandle;
46 custatevecDeviceMemHandler_t cuMemHandler;
47} CuQuantumConfig;
48# else
49# define CuQuantumConfig void*
50# endif
51
52
53
54// prevent C++ name mangling
55#ifdef __cplusplus
56extern "C" {
57#endif
58
59
60/*
61 * private structures
62 */
63
64// hide these from doxygen
66
72enum phaseGateType {SIGMA_Z=0, S_GATE=1, T_GATE=2};
73
79typedef struct {
80
81 char* buffer; // generated QASM string
82 int bufferSize; // maximum number of chars before overflow
83 int bufferFill; // number of chars currently in buffer
84 int isLogging; // whether gates are being added to buffer
85
86} QASMLogger;
87
94typedef struct ComplexArray
95{
96 qreal *real;
97 qreal *imag;
98} ComplexArray;
99
101
102
103
104/*
105 * public structures
106 */
107
114
120typedef struct Complex
121{
124} Complex;
125
154typedef struct ComplexMatrix2
155{
156 qreal real[2][2];
157 qreal imag[2][2];
159
192typedef struct ComplexMatrix4
193{
194 qreal real[4][4];
195 qreal imag[4][4];
197
203typedef struct ComplexMatrixN
204{
209
215typedef struct Vector
216{
218} Vector;
219
254
289
296typedef struct PauliHamil
297{
307} PauliHamil;
308
316typedef struct DiagonalOp
317{
321 long long int numElemsPerChunk;
331 ComplexArray deviceOperator;
332} DiagonalOp;
333
340typedef struct SubDiagonalOp
341{
345 long long int numElems;
350
352
360typedef struct Qureg
361{
370 long long int numAmpsPerChunk;
372 long long int numAmpsTotal;
377
379 ComplexArray stateVec;
381 ComplexArray pairStateVec;
382
384 ComplexArray deviceStateVec;
387
391 CuQuantumConfig* cuConfig;
392
394 QASMLogger* qasmLog;
395
396} Qureg;
397
405typedef struct QuESTEnv
406{
407 int rank;
409 unsigned long int* seeds;
411
412 // a copy of the QuESTEnv's config, used only in cuQuantum deployment
413 CuQuantumConfig* cuConfig;
414
415} QuESTEnv;
416
417
418
419/*
420 * public functions
421 */
422
579Qureg createQureg(int numQubits, QuESTEnv env);
580
673Qureg createDensityQureg(int numQubits, QuESTEnv env);
674
695
716void destroyQureg(Qureg qureg, QuESTEnv env);
717
776
794
795#ifndef __cplusplus
796#ifndef _WIN32
819void initComplexMatrixN(ComplexMatrixN m, qreal real[][1<<m.numQubits], qreal imag[][1<<m.numQubits]);
820#endif
821#endif
822
858PauliHamil createPauliHamil(int numQubits, int numSumTerms);
859
866void destroyPauliHamil(PauliHamil hamil);
867
915
953void initPauliHamil(PauliHamil hamil, qreal* coeffs, enum pauliOpType* codes);
954
1033DiagonalOp createDiagonalOp(int numQubits, QuESTEnv env);
1034
1048
1067void syncDiagonalOp(DiagonalOp op);
1068
1095void initDiagonalOp(DiagonalOp op, qreal* real, qreal* imag);
1096
1150
1194
1241void setDiagonalOpElems(DiagonalOp op, long long int startInd, qreal* real, qreal* imag, long long int numElems);
1242
1282void applyDiagonalOp(Qureg qureg, DiagonalOp op);
1283
1315
1360SubDiagonalOp createSubDiagonalOp(int numQubits);
1361
1374
1444void diagonalUnitary(Qureg qureg, int* targets, int numTargets, SubDiagonalOp op);
1445
1473void applyGateSubDiagonalOp(Qureg qureg, int* targets, int numTargets, SubDiagonalOp op);
1474
1513void applySubDiagonalOp(Qureg qureg, int* targets, int numTargets, SubDiagonalOp op);
1514
1538void reportState(Qureg qureg);
1539
1547void reportStateToScreen(Qureg qureg, QuESTEnv env, int reportRank);
1548
1555void reportQuregParams(Qureg qureg);
1556
1579void reportPauliHamil(PauliHamil hamil);
1580
1591int getNumQubits(Qureg qureg);
1592
1609long long int getNumAmps(Qureg qureg);
1610
1619void initBlankState(Qureg qureg);
1620
1633void initZeroState(Qureg qureg);
1634
1652void initPlusState(Qureg qureg);
1653
1689void initClassicalState(Qureg qureg, long long int stateInd);
1690
1709void initPureState(Qureg qureg, Qureg pure);
1710
1721void initDebugState(Qureg qureg);
1722
1748void initStateFromAmps(Qureg qureg, qreal* reals, qreal* imags);
1749
1797void setAmps(Qureg qureg, long long int startInd, qreal* reals, qreal* imags, long long int numAmps);
1798
1829void setDensityAmps(Qureg qureg, long long int startRow, long long int startCol, qreal* reals, qreal* imags, long long int numAmps);
1830
1854void setQuregToPauliHamil(Qureg qureg, PauliHamil hamil);
1855
1876void cloneQureg(Qureg targetQureg, Qureg copyQureg);
1877
1916void phaseShift(Qureg qureg, int targetQubit, qreal angle);
1917
1965void controlledPhaseShift(Qureg qureg, int idQubit1, int idQubit2, qreal angle);
1966
2012void multiControlledPhaseShift(Qureg qureg, int *controlQubits, int numControlQubits, qreal angle);
2013
2056void controlledPhaseFlip (Qureg qureg, int idQubit1, int idQubit2);
2057
2105void multiControlledPhaseFlip(Qureg qureg, int *controlQubits, int numControlQubits);
2106
2142void sGate(Qureg qureg, int targetQubit);
2143
2179void tGate(Qureg qureg, int targetQubit);
2180
2197
2209void destroyQuESTEnv(QuESTEnv env);
2210
2220void syncQuESTEnv(QuESTEnv env);
2221
2230int syncQuESTSuccess(int successCode);
2231
2238void reportQuESTEnv(QuESTEnv env);
2239
2257void getEnvironmentString(QuESTEnv env, char str[200]);
2258
2286void copyStateToGPU(Qureg qureg);
2287
2315void copyStateFromGPU(Qureg qureg);
2316
2351void copySubstateToGPU(Qureg qureg, long long int startInd, long long int numAmps);
2352
2383void copySubstateFromGPU(Qureg qureg, long long int startInd, long long int numAmps);
2384
2404Complex getAmp(Qureg qureg, long long int index);
2405
2425qreal getRealAmp(Qureg qureg, long long int index);
2426
2446qreal getImagAmp(Qureg qureg, long long int index);
2447
2467qreal getProbAmp(Qureg qureg, long long int index);
2468
2489Complex getDensityAmp(Qureg qureg, long long int row, long long int col);
2490
2517
2562void compactUnitary(Qureg qureg, int targetQubit, Complex alpha, Complex beta);
2563
2607void unitary(Qureg qureg, int targetQubit, ComplexMatrix2 u);
2608
2649void rotateX(Qureg qureg, int rotQubit, qreal angle);
2650
2692void rotateY(Qureg qureg, int rotQubit, qreal angle);
2693
2736void rotateZ(Qureg qureg, int rotQubit, qreal angle);
2737
2764void rotateAroundAxis(Qureg qureg, int rotQubit, qreal angle, Vector axis);
2765
2808void controlledRotateX(Qureg qureg, int controlQubit, int targetQubit, qreal angle);
2809
2851void controlledRotateY(Qureg qureg, int controlQubit, int targetQubit, qreal angle);
2852
2895void controlledRotateZ(Qureg qureg, int controlQubit, int targetQubit, qreal angle);
2896
2939void controlledRotateAroundAxis(Qureg qureg, int controlQubit, int targetQubit, qreal angle, Vector axis);
2940
2994void controlledCompactUnitary(Qureg qureg, int controlQubit, int targetQubit, Complex alpha, Complex beta);
2995
3049void controlledUnitary(Qureg qureg, int controlQubit, int targetQubit, ComplexMatrix2 u);
3050
3117void multiControlledUnitary(Qureg qureg, int* controlQubits, int numControlQubits, int targetQubit, ComplexMatrix2 u);
3118
3158void pauliX(Qureg qureg, int targetQubit);
3159
3197void pauliY(Qureg qureg, int targetQubit);
3198
3239void pauliZ(Qureg qureg, int targetQubit);
3240
3275void hadamard(Qureg qureg, int targetQubit);
3276
3323void controlledNot(Qureg qureg, int controlQubit, int targetQubit);
3324
3403void multiControlledMultiQubitNot(Qureg qureg, int* ctrls, int numCtrls, int* targs, int numTargs);
3404
3464void multiQubitNot(Qureg qureg, int* targs, int numTargs);
3465
3510void controlledPauliY(Qureg qureg, int controlQubit, int targetQubit);
3511
3544qreal calcProbOfOutcome(Qureg qureg, int measureQubit, int outcome);
3545
3633void calcProbOfAllOutcomes(qreal* outcomeProbs, Qureg qureg, int* qubits, int numQubits);
3634
3668qreal collapseToOutcome(Qureg qureg, int measureQubit, int outcome);
3669
3693int measure(Qureg qureg, int measureQubit);
3694
3719int measureWithStats(Qureg qureg, int measureQubit, qreal *outcomeProb);
3720
3747
3800
3823void seedQuESTDefault(QuESTEnv *env);
3824
3852void seedQuEST(QuESTEnv *env, unsigned long int *seedArray, int numSeeds);
3853
3896void getQuESTSeeds(QuESTEnv env, unsigned long int** seeds, int* numSeeds);
3897
3906void startRecordingQASM(Qureg qureg);
3907
3917void stopRecordingQASM(Qureg qureg);
3918
3925void clearRecordedQASM(Qureg qureg);
3926
3934void printRecordedQASM(Qureg qureg);
3935
3945void writeRecordedQASMToFile(Qureg qureg, char* filename);
3946
3976void mixDephasing(Qureg qureg, int targetQubit, qreal prob);
3977
4008void mixTwoQubitDephasing(Qureg qureg, int qubit1, int qubit2, qreal prob);
4009
4051void mixDepolarising(Qureg qureg, int targetQubit, qreal prob);
4052
4089void mixDamping(Qureg qureg, int targetQubit, qreal prob);
4090
4156void mixTwoQubitDepolarising(Qureg qureg, int qubit1, int qubit2, qreal prob);
4157
4197void mixPauli(Qureg qureg, int targetQubit, qreal probX, qreal probY, qreal probZ);
4198
4219void mixDensityMatrix(Qureg combineQureg, qreal prob, Qureg otherQureg);
4220
4247qreal calcPurity(Qureg qureg);
4248
4283qreal calcFidelity(Qureg qureg, Qureg pureState);
4284
4331void swapGate(Qureg qureg, int qubit1, int qubit2);
4332
4333
4383void sqrtSwapGate(Qureg qureg, int qb1, int qb2);
4384
4448void multiStateControlledUnitary(Qureg qureg, int* controlQubits, int* controlState, int numControlQubits, int targetQubit, ComplexMatrix2 u);
4449
4483void multiRotateZ(Qureg qureg, int* qubits, int numQubits, qreal angle);
4484
4538void multiRotatePauli(Qureg qureg, int* targetQubits, enum pauliOpType* targetPaulis, int numTargets, qreal angle);
4539
4616void multiControlledMultiRotateZ(Qureg qureg, int* controlQubits, int numControls, int* targetQubits, int numTargets, qreal angle);
4617
4726void multiControlledMultiRotatePauli(Qureg qureg, int* controlQubits, int numControls, int* targetQubits, enum pauliOpType* targetPaulis, int numTargets, qreal angle);
4727
4777qreal calcExpecPauliProd(Qureg qureg, int* targetQubits, enum pauliOpType* pauliCodes, int numTargets, Qureg workspace);
4778
4832qreal calcExpecPauliSum(Qureg qureg, enum pauliOpType* allPauliCodes, qreal* termCoeffs, int numSumTerms, Qureg workspace);
4833
4873qreal calcExpecPauliHamil(Qureg qureg, PauliHamil hamil, Qureg workspace);
4874
4945void twoQubitUnitary(Qureg qureg, int targetQubit1, int targetQubit2, ComplexMatrix4 u);
4946
5016void controlledTwoQubitUnitary(Qureg qureg, int controlQubit, int targetQubit1, int targetQubit2, ComplexMatrix4 u);
5017
5099void multiControlledTwoQubitUnitary(Qureg qureg, int* controlQubits, int numControlQubits, int targetQubit1, int targetQubit2, ComplexMatrix4 u);
5100
5193void multiQubitUnitary(Qureg qureg, int* targs, int numTargs, ComplexMatrixN u);
5194
5271void controlledMultiQubitUnitary(Qureg qureg, int ctrl, int* targs, int numTargs, ComplexMatrixN u);
5272
5366void multiControlledMultiQubitUnitary(Qureg qureg, int* ctrls, int numCtrls, int* targs, int numTargs, ComplexMatrixN u);
5367
5412void mixKrausMap(Qureg qureg, int target, ComplexMatrix2 *ops, int numOps);
5413
5453void mixTwoQubitKrausMap(Qureg qureg, int target1, int target2, ComplexMatrix4 *ops, int numOps);
5454
5505void mixMultiQubitKrausMap(Qureg qureg, int* targets, int numTargets, ComplexMatrixN* ops, int numOps);
5506
5540void mixNonTPKrausMap(Qureg qureg, int target, ComplexMatrix2 *ops, int numOps);
5541
5579void mixNonTPTwoQubitKrausMap(Qureg qureg, int target1, int target2, ComplexMatrix4 *ops, int numOps);
5580
5630void mixNonTPMultiQubitKrausMap(Qureg qureg, int* targets, int numTargets, ComplexMatrixN* ops, int numOps);
5631
5664
5688void setWeightedQureg(Complex fac1, Qureg qureg1, Complex fac2, Qureg qureg2, Complex facOut, Qureg out);
5689
5747void applyPauliSum(Qureg inQureg, enum pauliOpType* allPauliCodes, qreal* termCoeffs, int numSumTerms, Qureg outQureg);
5748
5791void applyPauliHamil(Qureg inQureg, PauliHamil hamil, Qureg outQureg);
5792
5871void applyTrotterCircuit(Qureg qureg, PauliHamil hamil, qreal time, int order, int reps);
5872
5892void applyMatrix2(Qureg qureg, int targetQubit, ComplexMatrix2 u);
5893
5944void applyMatrix4(Qureg qureg, int targetQubit1, int targetQubit2, ComplexMatrix4 u);
5945
6013void applyMatrixN(Qureg qureg, int* targs, int numTargs, ComplexMatrixN u);
6014
6043void applyGateMatrixN(Qureg qureg, int* targs, int numTargs, ComplexMatrixN u);
6044
6094void applyMultiControlledGateMatrixN(Qureg qureg, int* ctrls, int numCtrls, int* targs, int numTargs, ComplexMatrixN m);
6095
6147void applyMultiControlledMatrixN(Qureg qureg, int* ctrls, int numCtrls, int* targs, int numTargs, ComplexMatrixN u);
6148
6188void invalidQuESTInputError(const char* errMsg, const char* errFunc);
6189
6190#ifndef __cplusplus
6191#ifndef _WIN32
6192 // hide this function from doxygen
6194
6232ComplexMatrixN bindArraysToStackComplexMatrixN(
6233 int numQubits, qreal re[][1<<numQubits], qreal im[][1<<numQubits],
6234 qreal** reStorage, qreal** imStorage);
6235#endif
6236#endif
6238
6239// hide this function from doxygen
6241#define UNPACK_ARR(...) __VA_ARGS__
6243
6244#ifndef __cplusplus
6292#define getStaticComplexMatrixN(numQubits, re, im) \
6293 bindArraysToStackComplexMatrixN( \
6294 numQubits, \
6295 (qreal[1<<numQubits][1<<numQubits]) UNPACK_ARR re, \
6296 (qreal[1<<numQubits][1<<numQubits]) UNPACK_ARR im, \
6297 (double*[1<<numQubits]) {NULL}, (double*[1<<numQubits]) {NULL} \
6298 )
6299#endif
6300
6407void applyPhaseFunc(Qureg qureg, int* qubits, int numQubits, enum bitEncoding encoding, qreal* coeffs, qreal* exponents, int numTerms);
6408
6518void applyPhaseFuncOverrides(Qureg qureg, int* qubits, int numQubits, enum bitEncoding encoding, qreal* coeffs, qreal* exponents, int numTerms, long long int* overrideInds, qreal* overridePhases, int numOverrides);
6519
6679void applyMultiVarPhaseFunc(Qureg qureg, int* qubits, int* numQubitsPerReg, int numRegs, enum bitEncoding encoding, qreal* coeffs, qreal* exponents, int* numTermsPerReg);
6680
6761void applyMultiVarPhaseFuncOverrides(Qureg qureg, int* qubits, int* numQubitsPerReg, int numRegs, enum bitEncoding encoding, qreal* coeffs, qreal* exponents, int* numTermsPerReg, long long int* overrideInds, qreal* overridePhases, int numOverrides);
6762
6901void applyNamedPhaseFunc(Qureg qureg, int* qubits, int* numQubitsPerReg, int numRegs, enum bitEncoding encoding, enum phaseFunc functionNameCode);
6902
6974void applyNamedPhaseFuncOverrides(Qureg qureg, int* qubits, int* numQubitsPerReg, int numRegs, enum bitEncoding encoding, enum phaseFunc functionNameCode, long long int* overrideInds, qreal* overridePhases, int numOverrides);
6975
7104void applyParamNamedPhaseFunc(Qureg qureg, int* qubits, int* numQubitsPerReg, int numRegs, enum bitEncoding encoding, enum phaseFunc functionNameCode, qreal* params, int numParams);
7105
7179void 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);
7180
7277void applyFullQFT(Qureg qureg);
7278
7397void applyQFT(Qureg qureg, int* qubits, int numQubits);
7398
7421void applyProjector(Qureg qureg, int qubit, int outcome);
7422
7423// end prevention of C++ name mangling
7424#ifdef __cplusplus
7425}
7426#endif
7427
7428#endif // QUEST_H
7429
qreal getRealAmp(Qureg qureg, long long int index)
Get the real component of the complex probability amplitude at an index in the state vector.
Definition: QuEST.c:946
qreal calcExpecPauliProd(Qureg qureg, int *targetQubits, enum pauliOpType *pauliCodes, int numTargets, Qureg workspace)
Computes the expected value of a product of Pauli operators.
Definition: QuEST.c:1297
qreal calcProbOfOutcome(Qureg qureg, int measureQubit, int outcome)
Gives the probability of a specified qubit being measured in the given outcome (0 or 1).
Definition: QuEST.c:1262
Complex getAmp(Qureg qureg, long long int index)
Get the complex amplitude at a given index in the state vector.
Definition: QuEST.c:967
Complex calcInnerProduct(Qureg bra, Qureg ket)
Computes the inner product of two equal-size state vectors, given by.
Definition: QuEST.c:1246
qreal getProbAmp(Qureg qureg, long long int index)
Get the probability of a state-vector at an index in the full state vector.
Definition: QuEST.c:960
qreal calcExpecPauliSum(Qureg qureg, enum pauliOpType *allPauliCodes, qreal *termCoeffs, int numSumTerms, Qureg workspace)
Computes the expected value of a sum of products of Pauli operators.
Definition: QuEST.c:1306
qreal calcExpecPauliHamil(Qureg qureg, PauliHamil hamil, Qureg workspace)
Computes the expected value of qureg under Hermitian operator hamil.
Definition: QuEST.c:1315
qreal calcDensityInnerProduct(Qureg rho1, Qureg rho2)
Computes the Hilbert-Schmidt scalar product (which is equivalent to the Frobenius inner product of ma...
Definition: QuEST.c:1254
qreal calcPurity(Qureg qureg)
Calculates the purity of a density matrix, by the trace of the density matrix squared.
Definition: QuEST.c:1281
int getNumQubits(Qureg qureg)
Returns the number of qubits represented by qureg.
Definition: QuEST.c:936
Complex calcExpecDiagonalOp(Qureg qureg, DiagonalOp op)
Computes the expected value of the diagonal operator op for state qureg.
Definition: QuEST.c:1324
qreal calcHilbertSchmidtDistance(Qureg a, Qureg b)
Computes the Hilbert Schmidt distance between two density matrices a and b, defined as the Frobenius ...
Definition: QuEST.c:1333
long long int getNumAmps(Qureg qureg)
Returns the number of complex amplitudes in a state-vector qureg.
Definition: QuEST.c:940
qreal calcFidelity(Qureg qureg, Qureg pureState)
Calculates the fidelity of qureg (a state-vector or density matrix) against a reference pure state (n...
Definition: QuEST.c:1287
qreal calcTotalProb(Qureg qureg)
A debugging function which calculates the probability of the qubits in qureg being in any state,...
Definition: QuEST.c:1239
void calcProbOfAllOutcomes(qreal *outcomeProbs, Qureg qureg, int *qubits, int numQubits)
Populates outcomeProbs with the probabilities of every outcome of the sub-register contained in qubit...
Definition: QuEST.c:1272
qreal getImagAmp(Qureg qureg, long long int index)
Get the imaginary component of the complex probability amplitude at an index in the state vector.
Definition: QuEST.c:953
Complex getDensityAmp(Qureg qureg, long long int row, long long int col)
Get an amplitude from a density matrix at a given row and column.
Definition: QuEST.c:977
void copyStateToGPU(Qureg qureg)
In GPU mode, this copies the state-vector (or density matrix) from RAM (qureg.stateVec) to VRAM / GPU...
Definition: QuEST_cpu.c:43
void invalidQuESTInputError(const char *errMsg, const char *errFunc)
An internal function called when invalid arguments are passed to a QuEST API call,...
Definition: QuEST_validation.c:231
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.
Definition: QuEST_cpu_distributed.c:1400
void syncQuESTEnv(QuESTEnv env)
Guarantees that all code up to the given point has been executed on all nodes (if running in distribu...
Definition: QuEST_cpu_distributed.c:166
void copyStateFromGPU(Qureg qureg)
In GPU mode, this copies the state-vector (or density matrix) from GPU memory (qureg....
Definition: QuEST_cpu.c:46
void reportQuregParams(Qureg qureg)
Report metainformation about a set of qubits: number of qubits, number of probability amplitudes.
Definition: QuEST_common.c:233
void reportStateToScreen(Qureg qureg, QuESTEnv env, int reportRank)
Print the current state vector of probability amplitudes for a set of qubits to standard out.
Definition: QuEST.c:1722
void copySubstateToGPU(Qureg qureg, long long int startInd, long long int numAmps)
In GPU mode, this copies a substate of the state-vector (or density matrix) from RAM (qureg....
Definition: QuEST.c:1755
void initDebugState(Qureg qureg)
Initialises qureg to be in the un-normalised, non-physical state with with -th complex amplitude give...
Definition: QuEST.c:1718
void getQuESTSeeds(QuESTEnv env, unsigned long int **seeds, int *numSeeds)
Obtain the seeds presently used in random number generation.
Definition: QuEST.c:1750
void reportPauliHamil(PauliHamil hamil)
Print the PauliHamil to screen.
Definition: QuEST.c:1726
void copySubstateFromGPU(Qureg qureg, long long int startInd, long long int numAmps)
In GPU mode, this copies a substate of the state-vector (or density matrix) from GPU VRAM (qureg....
Definition: QuEST.c:1761
int syncQuESTSuccess(int successCode)
Performs a logical AND on all successCodes held by all processes.
Definition: QuEST_cpu_distributed.c:170
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 getEnvironmentString(QuESTEnv env, char str[200])
Sets str to a string containing information about the runtime environment, including whether simulati...
Definition: QuEST_cpu_distributed.c:200
void seedQuESTDefault(QuESTEnv *env)
Seeds the random number generator with the (master node) current time and process ID.
Definition: QuEST.c:1742
void reportQuESTEnv(QuESTEnv env)
Report information about the QuEST environment.
Definition: QuEST_cpu_distributed.c:185
void mixTwoQubitDephasing(Qureg qureg, int qubit1, int qubit2, qreal prob)
Mixes a density matrix qureg to induce two-qubit dephasing noise.
Definition: QuEST.c:1356
void mixTwoQubitKrausMap(Qureg qureg, int target1, int target2, ComplexMatrix4 *ops, int numOps)
Apply a general two-qubit Kraus map to a density matrix, as specified by at most sixteen Kraus operat...
Definition: QuEST.c:1420
void mixNonTPMultiQubitKrausMap(Qureg qureg, int *targets, int numTargets, ComplexMatrixN *ops, int numOps)
Apply a general N-qubit non-trace-preserving Kraus map to a density matrix, as specified by at most (...
Definition: QuEST.c:1460
void mixNonTPKrausMap(Qureg qureg, int target, ComplexMatrix2 *ops, int numOps)
Apply a general non-trace-preserving single-qubit Kraus map to a density matrix, as specified by at m...
Definition: QuEST.c:1440
void mixMultiQubitKrausMap(Qureg qureg, int *targets, int numTargets, ComplexMatrixN *ops, int numOps)
Apply a general N-qubit Kraus map to a density matrix, as specified by at most (2N)^2 Kraus operators...
Definition: QuEST.c:1430
void mixDensityMatrix(Qureg combineQureg, qreal prob, Qureg otherQureg)
Modifies combineQureg to become (1-prob)combineProb + prob otherQureg.
Definition: QuEST.c:1040
void mixPauli(Qureg qureg, int targetQubit, qreal probX, qreal probY, qreal probZ)
Mixes a density matrix qureg to induce general single-qubit Pauli noise.
Definition: QuEST.c:1399
void mixKrausMap(Qureg qureg, int target, ComplexMatrix2 *ops, int numOps)
Apply a general single-qubit Kraus map to a density matrix, as specified by at most four Kraus operat...
Definition: QuEST.c:1410
void mixDephasing(Qureg qureg, int targetQubit, qreal prob)
Mixes a density matrix qureg to induce single-qubit dephasing noise.
Definition: QuEST.c:1346
void mixTwoQubitDepolarising(Qureg qureg, int qubit1, int qubit2, qreal prob)
Mixes a density matrix qureg to induce two-qubit homogeneous depolarising noise.
Definition: QuEST.c:1387
void mixNonTPTwoQubitKrausMap(Qureg qureg, int target1, int target2, ComplexMatrix4 *ops, int numOps)
Apply a general non-trace-preserving two-qubit Kraus map to a density matrix, as specified by at most...
Definition: QuEST.c:1450
void mixDepolarising(Qureg qureg, int targetQubit, qreal prob)
Mixes a density matrix qureg to induce single-qubit homogeneous depolarising noise.
Definition: QuEST.c:1368
void mixDamping(Qureg qureg, int targetQubit, qreal prob)
Mixes a density matrix qureg to induce single-qubit amplitude damping (decay to 0 state).
Definition: QuEST.c:1379
void initClassicalState(Qureg qureg, long long int stateInd)
Initialise qureg into the classical state (also known as a "computational basis state") with index st...
Definition: QuEST.c:134
void initPlusState(Qureg qureg)
Initialise qureg into the plus state.
Definition: QuEST.c:125
void initZeroState(Qureg qureg)
Initialise qureg into the zero state.
Definition: QuEST.c:113
void cloneQureg(Qureg targetQureg, Qureg copyQureg)
Overwrite the amplitudes of targetQureg with those from copyQureg.
Definition: QuEST.c:164
void initStateFromAmps(Qureg qureg, qreal *reals, qreal *imags)
Initialise qureg by specifying all amplitudes.
Definition: QuEST.c:157
void setWeightedQureg(Complex fac1, Qureg qureg1, Complex fac2, Qureg qureg2, Complex facOut, Qureg out)
Modifies qureg out to the result of (facOut out + fac1 qureg1 + fac2 qureg2), imposing no constraints...
Definition: QuEST.c:1068
void setDensityAmps(Qureg qureg, long long int startRow, long long int startCol, qreal *reals, qreal *imags, long long int numAmps)
Overwrites a contiguous subset of the amplitudes in density-matrix qureg, with those passed in reals ...
Definition: QuEST.c:1058
void initPureState(Qureg qureg, Qureg pure)
Initialise qureg into to a given pure state of an equivalent Hilbert dimension.
Definition: QuEST.c:145
void setAmps(Qureg qureg, long long int startInd, qreal *reals, qreal *imags, long long int numAmps)
Overwrites a contiguous subset of the amplitudes in state-vector qureg, with those passed in reals an...
Definition: QuEST.c:1049
void setQuregToPauliHamil(Qureg qureg, PauliHamil hamil)
Overwrites the density-matrix qureg with the Z-basis matrix representation of the given real-weighted...
Definition: QuEST.c:171
void initBlankState(Qureg qureg)
Initialises a qureg to have all-zero-amplitudes.
Definition: QuEST.c:119
int measure(Qureg qureg, int measureQubit)
Measures a single qubit, collapsing it randomly to 0 or 1.
Definition: QuEST.c:1026
qreal collapseToOutcome(Qureg qureg, int measureQubit, int outcome)
Updates qureg to be consistent with measuring measureQubit in the given outcome (0 or 1),...
Definition: QuEST.c:994
int measureWithStats(Qureg qureg, int measureQubit, qreal *outcomeProb)
Measures a single qubit, collapsing it randomly to 0 or 1, and additionally gives the probability of ...
Definition: QuEST.c:1013
void applyMultiVarPhaseFunc(Qureg qureg, int *qubits, int *numQubitsPerReg, int numRegs, enum bitEncoding encoding, qreal *coeffs, qreal *exponents, int *numTermsPerReg)
Induces a phase change upon each amplitude of qureg, determined by a multi-variable exponential polyn...
Definition: QuEST.c:769
void applyParamNamedPhaseFunc(Qureg qureg, int *qubits, int *numQubitsPerReg, int numRegs, enum bitEncoding encoding, enum phaseFunc functionNameCode, qreal *params, int numParams)
Induces a phase change upon each amplitude of qureg, determined by a named, paramaterized (and potent...
Definition: QuEST.c:839
void applyProjector(Qureg qureg, int qubit, int outcome)
Force the target qubit of qureg into the given classical outcome, via a non-renormalising projection.
Definition: QuEST.c:896
void applySubDiagonalOp(Qureg qureg, int *targets, int numTargets, SubDiagonalOp op)
Left-apply a many-qubit a diagonal matrix upon a specific set of qubits of a quantum register.
Definition: QuEST.c:1206
void applyNamedPhaseFuncOverrides(Qureg qureg, int *qubits, int *numQubitsPerReg, int numRegs, enum bitEncoding encoding, enum phaseFunc functionNameCode, long long int *overrideInds, qreal *overridePhases, int numOverrides)
Induces a phase change upon each amplitude of qureg, determined by a named (and potentially multi-var...
Definition: QuEST.c:821
void applyMatrixN(Qureg qureg, int *targs, int numTargs, ComplexMatrixN u)
Apply a general N-by-N matrix, which may be non-unitary, on any number of target qubits.
Definition: QuEST.c:1134
void applyTrotterCircuit(Qureg qureg, PauliHamil hamil, qreal time, int order, int reps)
Applies a trotterisation of unitary evolution to qureg.
Definition: QuEST.c:1101
void applyPhaseFunc(Qureg qureg, int *qubits, int numQubits, enum bitEncoding encoding, qreal *coeffs, qreal *exponents, int numTerms)
Induces a phase change upon each amplitude of qureg, determined by the passed exponential polynomial ...
Definition: QuEST.c:734
void applyPauliHamil(Qureg inQureg, PauliHamil hamil, Qureg outQureg)
Modifies outQureg to be the result of applying PauliHamil (a Hermitian but not necessarily unitary op...
Definition: QuEST.c:1090
void applyGateSubDiagonalOp(Qureg qureg, int *targets, int numTargets, SubDiagonalOp op)
Apply a many-qubit unitary specified as a diagonal matrix upon a specific set of qubits of a quantum ...
Definition: QuEST.c:1215
void applyPauliSum(Qureg inQureg, enum pauliOpType *allPauliCodes, qreal *termCoeffs, int numSumTerms, Qureg outQureg)
Modifies outQureg to be the result of applying the weighted sum of Pauli products (a Hermitian but no...
Definition: QuEST.c:1079
void applyMatrix4(Qureg qureg, int targetQubit1, int targetQubit2, ComplexMatrix4 u)
Apply a general 4-by-4 matrix, which may be non-unitary.
Definition: QuEST.c:1124
void applyMultiVarPhaseFuncOverrides(Qureg qureg, int *qubits, int *numQubitsPerReg, int numRegs, enum bitEncoding encoding, qreal *coeffs, qreal *exponents, int *numTermsPerReg, long long int *overrideInds, qreal *overridePhases, int numOverrides)
Induces a phase change upon each amplitude of qureg, determined by a multi-variable exponential polyn...
Definition: QuEST.c:786
void applyQFT(Qureg qureg, int *qubits, int numQubits)
Applies the quantum Fourier transform (QFT) to a specific subset of qubits of the register qureg.
Definition: QuEST.c:874
void applyGateMatrixN(Qureg qureg, int *targs, int numTargs, ComplexMatrixN u)
Apply a gate specified by a general N-by-N matrix, which may be non-unitary, on any number of target ...
Definition: QuEST.c:1145
void 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)
Induces a phase change upon each amplitude of qureg, determined by a named, parameterised (and potent...
Definition: QuEST.c:856
void applyMultiControlledMatrixN(Qureg qureg, int *ctrls, int numCtrls, int *targs, int numTargs, ComplexMatrixN u)
Apply a general N-by-N matrix, which may be non-unitary, with additional controlled qubits.
Definition: QuEST.c:1182
void applyNamedPhaseFunc(Qureg qureg, int *qubits, int *numQubitsPerReg, int numRegs, enum bitEncoding encoding, enum phaseFunc functionNameCode)
Induces a phase change upon each amplitude of qureg, determined by a named (and potentially multi-var...
Definition: QuEST.c:804
void applyDiagonalOp(Qureg qureg, DiagonalOp op)
Apply a diagonal operator, which is possibly non-unitary and non-Hermitian, to the entire qureg.
Definition: QuEST.c:1195
void applyMatrix2(Qureg qureg, int targetQubit, ComplexMatrix2 u)
Apply a general 2-by-2 matrix, which may be non-unitary.
Definition: QuEST.c:1115
void applyFullQFT(Qureg qureg)
Applies the quantum Fourier transform (QFT) to the entirety of qureg.
Definition: QuEST.c:884
void applyPhaseFuncOverrides(Qureg qureg, int *qubits, int numQubits, enum bitEncoding encoding, qreal *coeffs, qreal *exponents, int numTerms, long long int *overrideInds, qreal *overridePhases, int numOverrides)
Induces a phase change upon each amplitude of qureg, determined by the passed exponential polynomial ...
Definition: QuEST.c:751
void writeRecordedQASMToFile(Qureg qureg, char *filename)
Writes recorded QASM to a file, throwing an error if inaccessible.
Definition: QuEST.c:103
void startRecordingQASM(Qureg qureg)
Enable QASM recording.
Definition: QuEST.c:87
void clearRecordedQASM(Qureg qureg)
Clear all QASM so far recorded.
Definition: QuEST.c:95
void stopRecordingQASM(Qureg qureg)
Disable QASM recording.
Definition: QuEST.c:91
void printRecordedQASM(Qureg qureg)
Print recorded QASM to stdout.
Definition: QuEST.c:99
void initDiagonalOp(DiagonalOp op, qreal *real, qreal *imag)
Overwrites the entire DiagonalOp op with the given real and imag complex elements.
Definition: QuEST.c:1663
pauliOpType
Codes for specifying Pauli operators.
Definition: QuEST.h:113
void destroyDiagonalOp(DiagonalOp op, QuESTEnv env)
Destroys a DiagonalOp created with createDiagonalOp(), freeing its memory.
Definition: QuEST.c:1650
void initDiagonalOpFromPauliHamil(DiagonalOp op, PauliHamil hamil)
Populates the diagonal operator op to be equivalent to the given Pauli Hamiltonian hamil,...
Definition: QuEST.c:1676
SubDiagonalOp createSubDiagonalOp(int numQubits)
Creates a SubDiagonalOp representing a diagonal operator which can act upon a subset of the qubits in...
Definition: QuEST.c:1695
Qureg createQureg(int numQubits, QuESTEnv env)
Creates a state-vector Qureg object representing a set of qubits which will remain in a pure state.
Definition: QuEST.c:36
PauliHamil createPauliHamil(int numQubits, int numSumTerms)
Dynamically allocates a Hamiltonian expressed as a real-weighted sum of products of Pauli operators.
Definition: QuEST.c:1524
void initComplexMatrixN(ComplexMatrixN m, qreal real[][1<< m.numQubits], qreal imag[][1<< m.numQubits])
Initialises a ComplexMatrixN instance to have the passed real and imag values.
Definition: QuEST.c:1512
void destroySubDiagonalOp(SubDiagonalOp op)
Destroy a SubDiagonalOp instance created with createSubDiagonalOp().
Definition: QuEST.c:1709
#define qreal
A precision-agnostic floating point number, as determined by QuEST_PREC.
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:1474
QuESTEnv createQuESTEnv(void)
Create the QuEST execution environment.
Definition: QuEST_cpu_distributed.c:131
void destroyQureg(Qureg qureg, QuESTEnv env)
Deallocate a Qureg, freeing its memory.
Definition: QuEST.c:77
Qureg createDensityQureg(int numQubits, QuESTEnv env)
Creates a density matrix Qureg object representing a set of qubits which can enter noisy and mixed st...
Definition: QuEST.c:50
bitEncoding
Flags for specifying how the bits in sub-register computational basis states are mapped to indices in...
Definition: QuEST.h:288
void destroyComplexMatrixN(ComplexMatrixN matr)
Destroy a ComplexMatrixN instance created with createComplexMatrixN()
Definition: QuEST.c:1495
PauliHamil createPauliHamilFromFile(char *fn)
Creates a PauliHamil instance, a real-weighted sum of products of Pauli operators,...
Definition: QuEST.c:1546
phaseFunc
Flags for specifying named phase functions.
Definition: QuEST.h:249
DiagonalOp createDiagonalOpFromPauliHamilFile(char *fn, QuESTEnv env)
Creates and initialiases a diagonal operator from the Z Pauli Hamiltonian encoded in file with filena...
Definition: QuEST.c:1684
void syncDiagonalOp(DiagonalOp op)
Update the GPU memory with the current values in op.real and op.imag.
Definition: QuEST.c:1657
Qureg createCloneQureg(Qureg qureg, QuESTEnv env)
Create a new Qureg which is an exact clone of the passed qureg, which can be either a state-vector or...
Definition: QuEST.c:64
void destroyPauliHamil(PauliHamil hamil)
Destroy a PauliHamil instance, created with either createPauliHamil() or createPauliHamilFromFile().
Definition: QuEST.c:1540
void initPauliHamil(PauliHamil hamil, qreal *coeffs, enum pauliOpType *codes)
Initialise PauliHamil instance hamil with the given term coefficients and Pauli codes (one for every ...
Definition: QuEST.c:1630
DiagonalOp createDiagonalOp(int numQubits, QuESTEnv env)
Creates a DiagonalOp representing a diagonal operator on the full Hilbert space of a Qureg.
Definition: QuEST.c:1644
void destroyQuESTEnv(QuESTEnv env)
Destroy the QuEST environment.
Definition: QuEST_cpu_distributed.c:176
void setDiagonalOpElems(DiagonalOp op, long long int startInd, qreal *real, qreal *imag, long long int numElems)
Modifies a subset (starting at index startInd, and ending at index startInd + numElems) of the elemen...
Definition: QuEST.c:1669
@ PAULI_Z
Definition: QuEST.h:113
@ PAULI_Y
Definition: QuEST.h:113
@ PAULI_I
Definition: QuEST.h:113
@ PAULI_X
Definition: QuEST.h:113
@ UNSIGNED
Definition: QuEST.h:288
@ TWOS_COMPLEMENT
Definition: QuEST.h:288
@ SCALED_INVERSE_PRODUCT
Definition: QuEST.h:251
@ DISTANCE
Definition: QuEST.h:252
@ SCALED_PRODUCT
Definition: QuEST.h:251
@ SCALED_INVERSE_SHIFTED_DISTANCE
Definition: QuEST.h:252
@ INVERSE_DISTANCE
Definition: QuEST.h:252
@ SCALED_NORM
Definition: QuEST.h:250
@ SCALED_INVERSE_SHIFTED_WEIGHTED_DISTANCE
Definition: QuEST.h:252
@ SCALED_INVERSE_SHIFTED_NORM
Definition: QuEST.h:250
@ INVERSE_PRODUCT
Definition: QuEST.h:251
@ PRODUCT
Definition: QuEST.h:251
@ SCALED_DISTANCE
Definition: QuEST.h:252
@ INVERSE_NORM
Definition: QuEST.h:250
@ NORM
Definition: QuEST.h:250
@ SCALED_INVERSE_DISTANCE
Definition: QuEST.h:252
@ SCALED_INVERSE_NORM
Definition: QuEST.h:250
void multiControlledPhaseShift(Qureg qureg, int *controlQubits, int numControlQubits, qreal angle)
Introduce a phase factor on state of the passed qubits.
Definition: QuEST.c:518
void controlledRotateAroundAxis(Qureg qureg, int controlQubit, int targetQubit, qreal angle, Vector axis)
Applies a controlled rotation by a given angle around a given vector on the Bloch-sphere.
Definition: QuEST.c:622
void controlledRotateX(Qureg qureg, int controlQubit, int targetQubit, qreal angle)
Applies a controlled rotation by a given angle around the X-axis of the Bloch-sphere.
Definition: QuEST.c:228
void multiControlledMultiRotateZ(Qureg qureg, int *controlQubits, int numControls, int *targetQubits, int numTargets, qreal angle)
Apply a multi-controlled multi-target Z rotation, also known as a controlled phase gadget.
Definition: QuEST.c:676
void multiQubitUnitary(Qureg qureg, int *targs, int numTargs, ComplexMatrixN u)
Apply a general multi-qubit unitary (including a global phase factor) with any number of target qubit...
Definition: QuEST.c:304
void phaseShift(Qureg qureg, int targetQubit, qreal angle)
Shift the phase between and of a single qubit by a given angle.
Definition: QuEST.c:495
void diagonalUnitary(Qureg qureg, int *targets, int numTargets, SubDiagonalOp op)
Apply a many-qubit unitary specified as a diagonal matrix upon a specific set of qubits of a quantum ...
Definition: QuEST.c:910
void controlledPhaseFlip(Qureg qureg, int idQubit1, int idQubit2)
Apply the (two-qubit) controlled phase flip gate, also known as the controlled pauliZ gate.
Definition: QuEST.c:583
void multiRotateZ(Qureg qureg, int *qubits, int numQubits, qreal angle)
Apply a multi-qubit Z rotation, also known as a phase gadget, on a selected number of qubits.
Definition: QuEST.c:660
void rotateY(Qureg qureg, int rotQubit, qreal angle)
Rotate a single qubit by a given angle around the Y-axis of the Bloch-sphere.
Definition: QuEST.c:206
void multiRotatePauli(Qureg qureg, int *targetQubits, enum pauliOpType *targetPaulis, int numTargets, qreal angle)
Apply a multi-qubit multi-Pauli rotation, also known as a Pauli gadget, on a selected number of qubit...
Definition: QuEST.c:693
void multiControlledMultiQubitUnitary(Qureg qureg, int *ctrls, int numCtrls, int *targs, int numTargs, ComplexMatrixN u)
Apply a general multi-controlled multi-qubit unitary (including a global phase factor).
Definition: QuEST.c:338
void controlledTwoQubitUnitary(Qureg qureg, int controlQubit, int targetQubit1, int targetQubit2, ComplexMatrix4 u)
Apply a general controlled two-qubit unitary (including a global phase factor).
Definition: QuEST.c:277
void pauliX(Qureg qureg, int targetQubit)
Apply the single-qubit Pauli-X (also known as the X, sigma-X, NOT or bit-flip) gate.
Definition: QuEST.c:440
void multiStateControlledUnitary(Qureg qureg, int *controlQubits, int *controlState, int numControlQubits, int targetQubit, ComplexMatrix2 u)
Apply a general single-qubit unitary with multiple control qubits, conditioned upon a specific bit se...
Definition: QuEST.c:396
void applyMultiControlledGateMatrixN(Qureg qureg, int *ctrls, int numCtrls, int *targs, int numTargs, ComplexMatrixN m)
Apply a general multi-controlled multi-qubit gate specified as an (possibly non-unitary) arbitrary co...
Definition: QuEST.c:1163
void rotateAroundAxis(Qureg qureg, int rotQubit, qreal angle, Vector axis)
Rotate a single qubit by a given angle around a given Vector on the Bloch-sphere.
Definition: QuEST.c:609
void pauliZ(Qureg qureg, int targetQubit)
Apply the single-qubit Pauli-Z (also known as the Z, sigma-Z or phase-flip) gate.
Definition: QuEST.c:462
void controlledPauliY(Qureg qureg, int controlQubit, int targetQubit)
Apply the controlled pauliY (single control, single target) gate, also known as the c-Y and c-sigma-Y...
Definition: QuEST.c:571
void hadamard(Qureg qureg, int targetQubit)
Apply the single-qubit Hadamard gate.
Definition: QuEST.c:184
void sGate(Qureg qureg, int targetQubit)
Apply the single-qubit S gate.
Definition: QuEST.c:473
void multiControlledTwoQubitUnitary(Qureg qureg, int *controlQubits, int numControlQubits, int targetQubit1, int targetQubit2, ComplexMatrix4 u)
Apply a general multi-controlled two-qubit unitary (including a global phase factor).
Definition: QuEST.c:290
void controlledPhaseShift(Qureg qureg, int idQubit1, int idQubit2, qreal angle)
Introduce a phase factor on state of qubits idQubit1 and idQubit2.
Definition: QuEST.c:506
void sqrtSwapGate(Qureg qureg, int qb1, int qb2)
Performs a sqrt SWAP gate between qubit1 and qubit2.
Definition: QuEST.c:647
void unitary(Qureg qureg, int targetQubit, ComplexMatrix2 u)
Apply a general single-qubit unitary (including a global phase factor).
Definition: QuEST.c:356
void multiQubitNot(Qureg qureg, int *targs, int numTargs)
Apply a NOT (or Pauli X) gate with multiple target qubits, which has the same effect as (but is much ...
Definition: QuEST.c:544
void controlledRotateZ(Qureg qureg, int controlQubit, int targetQubit, qreal angle)
Applies a controlled rotation by a given angle around the Z-axis of the Bloch-sphere.
Definition: QuEST.c:252
void swapGate(Qureg qureg, int qubit1, int qubit2)
Performs a SWAP gate between qubit1 and qubit2.
Definition: QuEST.c:635
void rotateX(Qureg qureg, int rotQubit, qreal angle)
Rotate a single qubit by a given angle around the X-axis of the Bloch-sphere.
Definition: QuEST.c:195
void pauliY(Qureg qureg, int targetQubit)
Apply the single-qubit Pauli-Y (also known as the Y or sigma-Y) gate.
Definition: QuEST.c:451
void controlledUnitary(Qureg qureg, int controlQubit, int targetQubit, ComplexMatrix2 u)
Apply a general controlled unitary (single control, single target), which can include a global phase ...
Definition: QuEST.c:368
void controlledNot(Qureg qureg, int controlQubit, int targetQubit)
Apply the controlled not (single control, single target) gate, also known as the c-X,...
Definition: QuEST.c:532
void multiControlledUnitary(Qureg qureg, int *controlQubits, int numControlQubits, int targetQubit, ComplexMatrix2 u)
Apply a general multiple-control single-target unitary, which can include a global phase factor.
Definition: QuEST.c:381
void rotateZ(Qureg qureg, int rotQubit, qreal angle)
Rotate a single qubit by a given angle around the Z-axis of the Bloch-sphere (also known as a phase s...
Definition: QuEST.c:217
void multiControlledMultiRotatePauli(Qureg qureg, int *controlQubits, int numControls, int *targetQubits, enum pauliOpType *targetPaulis, int numTargets, qreal angle)
Apply a multi-controlled multi-target multi-Pauli rotation, also known as a controlled Pauli gadget.
Definition: QuEST.c:713
void twoQubitUnitary(Qureg qureg, int targetQubit1, int targetQubit2, ComplexMatrix4 u)
Apply a general two-qubit unitary (including a global phase factor).
Definition: QuEST.c:264
void tGate(Qureg qureg, int targetQubit)
Apply the single-qubit T gate.
Definition: QuEST.c:484
void multiControlledPhaseFlip(Qureg qureg, int *controlQubits, int numControlQubits)
Apply the multiple-qubit controlled phase flip gate, also known as the multiple-qubit controlled paul...
Definition: QuEST.c:595
void controlledMultiQubitUnitary(Qureg qureg, int ctrl, int *targs, int numTargs, ComplexMatrixN u)
Apply a general controlled multi-qubit unitary (including a global phase factor).
Definition: QuEST.c:321
void compactUnitary(Qureg qureg, int targetQubit, Complex alpha, Complex beta)
Apply a single-qubit unitary parameterised by two given complex scalars.
Definition: QuEST.c:412
void controlledCompactUnitary(Qureg qureg, int controlQubit, int targetQubit, Complex alpha, Complex beta)
Apply a controlled unitary (single control, single target) parameterised by two given complex scalars...
Definition: QuEST.c:425
void controlledRotateY(Qureg qureg, int controlQubit, int targetQubit, qreal angle)
Applies a controlled rotation by a given angle around the Y-axis of the Bloch-sphere.
Definition: QuEST.c:240
void multiControlledMultiQubitNot(Qureg qureg, int *ctrls, int numCtrls, int *targs, int numTargs)
Apply a NOT (or Pauli X) gate with multiple control and target qubits.
Definition: QuEST.c:557
Represents a 2x2 matrix of complex numbers.
Definition: QuEST.h:155
qreal imag[2][2]
Definition: QuEST.h:157
qreal real[2][2]
Definition: QuEST.h:156
Represents a 4x4 matrix of complex numbers.
Definition: QuEST.h:193
qreal imag[4][4]
Definition: QuEST.h:195
qreal real[4][4]
Definition: QuEST.h:194
Represents a general 2^N by 2^N matrix of complex numbers.
Definition: QuEST.h:204
qreal ** real
Definition: QuEST.h:206
int numQubits
Definition: QuEST.h:205
qreal ** imag
Definition: QuEST.h:207
Represents one complex number.
Definition: QuEST.h:121
qreal imag
Definition: QuEST.h:123
qreal real
Definition: QuEST.h:122
Represents a diagonal complex operator on the full Hilbert state of a Qureg.
Definition: QuEST.h:317
qreal * real
The real values of the 2^numQubits complex elements.
Definition: QuEST.h:327
long long int numElemsPerChunk
The number of the 2^numQubits amplitudes stored on each distributed node.
Definition: QuEST.h:321
int numQubits
The number of qubits this operator can act on (informing its size)
Definition: QuEST.h:319
qreal * imag
The imaginary values of the 2^numQubits complex elements.
Definition: QuEST.h:329
int chunkId
The position of the chunk of the operator held by this process in the full operator.
Definition: QuEST.h:325
int numChunks
The number of nodes between which the elements of this operator are split.
Definition: QuEST.h:323
ComplexArray deviceOperator
A copy of the elements stored persistently on the GPU.
Definition: QuEST.h:331
A Pauli Hamiltonian, expressed as a real-weighted sum of pauli products, and which can hence represen...
Definition: QuEST.h:297
enum pauliOpType * pauliCodes
The Pauli operators acting on each qubit, flattened over every operator.
Definition: QuEST.h:300
int numSumTerms
The number of terms in the weighted sum, or the number of Pauli products.
Definition: QuEST.h:304
qreal * termCoeffs
The real coefficient of each Pauli product. This is an array of length PauliHamil....
Definition: QuEST.h:302
int numQubits
The number of qubits informing the Hilbert dimension of the Hamiltonian.
Definition: QuEST.h:306
Information about the environment the program is running in.
Definition: QuEST.h:406
CuQuantumConfig * cuConfig
Definition: QuEST.h:413
int numSeeds
Definition: QuEST.h:410
unsigned long int * seeds
Definition: QuEST.h:409
int rank
Definition: QuEST.h:407
int numRanks
Definition: QuEST.h:408
Represents a system of qubits.
Definition: QuEST.h:361
int numQubitsInStateVec
Number of qubits in the state-vector - this is double the number represented for mixed states.
Definition: QuEST.h:367
cuAmp * deviceCuStateVec
Definition: QuEST.h:390
qreal * firstLevelReduction
Storage for reduction of probabilities on GPU.
Definition: QuEST.h:386
ComplexArray stateVec
Computational state amplitudes - a subset thereof in the MPI version.
Definition: QuEST.h:379
QASMLogger * qasmLog
Storage for generated QASM output.
Definition: QuEST.h:394
CuQuantumConfig * cuConfig
Definition: QuEST.h:391
qreal * secondLevelReduction
Definition: QuEST.h:386
int numChunks
Number of chunks the state vector is broken up into – the number of MPI processes used.
Definition: QuEST.h:376
long long int numAmpsPerChunk
Number of probability amplitudes held in stateVec by this process In the non-MPI version,...
Definition: QuEST.h:370
cuAmp * cuStateVec
Storage for wavefunction amplitues and config (copy of QuESTEnv's handle) in cuQuantum deployment.
Definition: QuEST.h:389
ComplexArray pairStateVec
Temporary storage for a chunk of the state vector received from another process in the MPI version.
Definition: QuEST.h:381
ComplexArray deviceStateVec
Storage for wavefunction amplitudes in the GPU version.
Definition: QuEST.h:384
int chunkId
The position of the chunk of the state vector held by this process in the full state vector.
Definition: QuEST.h:374
int isDensityMatrix
Whether this instance is a density-state representation.
Definition: QuEST.h:363
int numQubitsRepresented
The number of qubits represented in either the state-vector or density matrix.
Definition: QuEST.h:365
long long int numAmpsTotal
Total number of amplitudes, which are possibly distributed among machines.
Definition: QuEST.h:372
Represents a diagonal complex operator of a smaller dimension than the full Hilbert state of a Qureg.
Definition: QuEST.h:341
long long int numElems
The number of diagonal elements, i.e. 2^numQubits.
Definition: QuEST.h:345
qreal * real
The real values of the 2^numQubits complex elements.
Definition: QuEST.h:347
qreal * imag
The imaginary values of the 2^numQubits complex elements.
Definition: QuEST.h:349
int numQubits
The number of target qubits which this SubDiagonalOp can operate upon.
Definition: QuEST.h:343
Represents a 3-vector of real numbers.
Definition: QuEST.h:216
qreal x
Definition: QuEST.h:217
qreal z
Definition: QuEST.h:217
qreal y
Definition: QuEST.h:217