Unit test utilities

Functions used in the unit testing. More...

Typedefs

typedef std::vector< std::vector< qcomp > > QMatrix
 A complex square matrix. More...
 
typedef std::vector< qcompQVector
 A complex vector, which can be zero-initialised with QVector(numAmps). More...
 

Functions

void applyReferenceMatrix (QMatrix &state, int *ctrls, int numCtrls, int *targs, int numTargs, QMatrix op)
 Modifies the density matrix state to be the result of left-multiplying the multi-target operator matrix op, with the specified control and target qubits (in ctrls and targs respectively). More...
 
void applyReferenceMatrix (QMatrix &state, int *targs, int numTargs, QMatrix op)
 Modifies the density matrix state to be the result of left-multiplying the multi-target operator matrix op, with the target qubits (assuming no control qubits). More...
 
void applyReferenceMatrix (QVector &state, int *ctrls, int numCtrls, int *targs, int numTargs, QMatrix op)
 Modifies the state-vector state to be the result of left-multiplying the multi-target operator matrix op, with the specified control and target qubits (in ctrls and targs respectively). More...
 
void applyReferenceMatrix (QVector &state, int *targs, int numTargs, QMatrix op)
 Modifies the state-vector state to be the result of left-multiplying the multi-target operator matrix op, with the specified target qubits (assuming no control qubits). More...
 
void applyReferenceOp (QMatrix &state, int *ctrls, int numCtrls, int *targs, int numTargs, QMatrix op)
 Modifies the density matrix state to be the result of applying the multi-target operator matrix op, with the specified control and target qubits (in ctrls and targs respectively). More...
 
void applyReferenceOp (QMatrix &state, int *ctrls, int numCtrls, int targ1, int targ2, QMatrix op)
 Modifies the density matrix state to be the result of applying the two-target operator matrix op, with the specified control qubits (in ctrls). More...
 
void applyReferenceOp (QMatrix &state, int *ctrls, int numCtrls, int target, QMatrix op)
 Modifies the density matrix state to be the result of applying the single-target operator matrix op, with the specified control qubits (in ctrls). More...
 
void applyReferenceOp (QMatrix &state, int *targs, int numTargs, QMatrix op)
 Modifies the density matrix state to be the result of applying the multi-target operator matrix op, with no control qubits. More...
 
void applyReferenceOp (QMatrix &state, int ctrl, int *targs, int numTargs, QMatrix op)
 Modifies the density matrix state to be the result of applying the multi-target operator matrix op, with a single control qubit ctrl. More...
 
void applyReferenceOp (QMatrix &state, int ctrl, int targ, QMatrix op)
 Modifies the density matrix state to be the result of applying the single-control single-target operator matrix op. More...
 
void applyReferenceOp (QMatrix &state, int ctrl, int targ1, int targ2, QMatrix op)
 Modifies the density matrix state to be the result of applying the two-target operator matrix op, with a single control qubit ctrl. More...
 
void applyReferenceOp (QMatrix &state, int targ, QMatrix op)
 Modifies the density matrix state to be the result of applying the single-target operator matrix op, with no control qubit. More...
 
void applyReferenceOp (QVector &state, int *ctrls, int numCtrls, int *targs, int numTargs, QMatrix op)
 Modifies the state-vector state to be the result of applying the multi-target operator matrix op, with the specified control and target qubits (in ctrls and targs respectively). More...
 
void applyReferenceOp (QVector &state, int *ctrls, int numCtrls, int targ1, int targ2, QMatrix op)
 Modifies the state-vector state to be the result of applying the two-target operator matrix op, with the specified control qubits (in ctrls). More...
 
void applyReferenceOp (QVector &state, int *ctrls, int numCtrls, int target, QMatrix op)
 Modifies the state-vector state to be the result of applying the single-target operator matrix op, with the specified control qubits (in ctrls). More...
 
void applyReferenceOp (QVector &state, int *targs, int numTargs, QMatrix op)
 Modifies the state-vector state to be the result of applying the multi-target operator matrix op, with no contorl qubits. More...
 
void applyReferenceOp (QVector &state, int ctrl, int *targs, int numTargs, QMatrix op)
 Modifies the state-vector state to be the result of applying the multi-target operator matrix op, with a single control qubit (ctrl) This updates state under. More...
 
void applyReferenceOp (QVector &state, int ctrl, int targ, QMatrix op)
 Modifies the state-vector state to be the result of applying the single-target operator matrix op, with a single control qubit (ctrl). More...
 
void applyReferenceOp (QVector &state, int ctrl, int targ1, int targ2, QMatrix op)
 Modifies the state-vector state to be the result of applying the two-target operator matrix op, with a single control qubit (ctrl). More...
 
void applyReferenceOp (QVector &state, int targ, QMatrix op)
 Modifies the state-vector state to be the result of applying the single-target operator matrix op, with no contorl qubits. More...
 
bool areEqual (QMatrix a, QMatrix b)
 Returns true if the absolute value of the difference between every amplitude in matrices a and b is less than REAL_EPS. More...
 
bool areEqual (Qureg qureg, QMatrix matr)
 Performs a hardware-agnostic comparison of density-matrix qureg to matr, checking whether the difference between the real and imaginary components of every amplitude is smaller than the QuEST_PREC-specific REAL_EPS (defined in QuEST_precision) precision. More...
 
bool areEqual (Qureg qureg, QMatrix matr, qreal precision)
 Performs a hardware-agnostic comparison of density-matrix qureg to matr, checking whether the difference between the real and imaginary components of every amplitude is smaller than precision. More...
 
bool areEqual (Qureg qureg, QVector vec)
 Performs a hardware-agnostic comparison of state-vector qureg to vec, checking whether the difference between the real and imaginary components of every amplitude is smaller than the QuEST_PREC-specific REAL_EPS (defined in QuEST_precision) precision. More...
 
bool areEqual (Qureg qureg, QVector vec, qreal precision)
 Performs a hardware-agnostic comparison of state-vector qureg to vec, checking whether the difference between the real and imaginary components of every amplitude is smaller than precision. More...
 
bool areEqual (Qureg qureg1, Qureg qureg2)
 Performs a hardware-agnostic comparison of the given quregs, checking whether the difference between the real and imaginary components of every amplitude is smaller than the QuEST_PREC-specific REAL_EPS (defined in QuEST_precision) precision. More...
 
bool areEqual (Qureg qureg1, Qureg qureg2, qreal precision)
 Performs a hardware-agnostic comparison of the given quregs, checking whether the difference between the real and imaginary components of every amplitude is smaller than precision. More...
 
bool areEqual (QVector a, QVector b)
 Returns true if the absolute value of the difference between every amplitude in vectors a and b is less than REAL_EPS. More...
 
bool areEqual (QVector vec, qreal *reals)
 Returns true if the absolute value of the difference between every element in vec (which must be strictly real) and those implied by reals, is less than REAL_EPS. More...
 
bool areEqual (QVector vec, qreal *reals, qreal *imags)
 Returns true if the absolute value of the difference between every element in vec and those implied by reals and imags, is less than REAL_EPS. More...
 
void assertQuregAndRefInDebugState (Qureg qureg, QMatrix ref)
 Asserts the given density qureg and reference agree, and are properly initialised in the debug state. More...
 
void assertQuregAndRefInDebugState (Qureg qureg, QVector ref)
 Asserts the given statevector qureg and reference agree, and are properly initialised in the debug state. More...
 
CatchGen< int * > bitsets (int numBits)
 Returns a Catch2 generator of every numBits-length bit-set, in increasing lexographic order, where left-most (zero index) bit is treated as LEAST significant (opposite typical convention). More...
 
unsigned int calcLog2 (long unsigned int res)
 Returns log2 of numbers which must be gauranteed to be 2^n. More...
 
void deleteFilesWithPrefixSynch (char *prefix)
 Deletes all files with filename starting with prefix. More...
 
QMatrix getConjugateTranspose (QMatrix a)
 Returns the conjugate transpose of the complex square matrix a. More...
 
QVector getDFT (QVector in)
 Returns the discrete fourier transform of vector in. More...
 
QVector getDFT (QVector in, int *targs, int numTargs)
 Returns the discrete fourier transform of a sub-partition of the vector in. More...
 
QMatrix getExponentialOfDiagonalMatrix (QMatrix a)
 Returns the matrix exponential of a diagonal, square, complex matrix. More...
 
QMatrix getExponentialOfPauliMatrix (qreal angle, QMatrix a)
 Returns the matrix exponential of a kronecker product of pauli matrices (or of any involutory matrices), with exponent factor (-i angle / 2). More...
 
QMatrix getFullOperatorMatrix (int *ctrls, int numCtrls, int *targs, int numTargs, QMatrix op, int numQubits)
 Takes a 2^numTargs-by-2^numTargs matrix op and a returns a 2^numQubits-by-2^numQubits matrix where op is controlled on the given ctrls qubits. More...
 
QMatrix getIdentityMatrix (size_t dim)
 Returns a dim-by-dim identity matrix. More...
 
QMatrix getKetBra (QVector ket, QVector bra)
 Returns the matrix |ket><bra|, with ith-jth element ket(i) conj(bra(j)), since |ket><bra| = sum_i a_i|i> sum_j b_j* <j| = sum_{ij} a_i b_j* |i><j|. More...
 
QMatrix getKroneckerProduct (QMatrix a, QMatrix b)
 Returns the kronecker product of a and b, where a and b are square but possibly differently-sized complex matrices. More...
 
QVector getKroneckerProduct (QVector b, QVector a)
 Returns b (otimes) a. More...
 
QVector getMatrixDiagonal (QMatrix matr)
 Returns the diagonal vector of the given matrix. More...
 
QMatrix getMixedDensityMatrix (std::vector< qreal > probs, std::vector< QVector > states)
 Returns a mixed density matrix formed from mixing the given pure states, which are assumed normalised, but not necessarily orthogonal. More...
 
QVector getNormalised (QVector vec)
 Returns an L2-normalised copy of vec, using Kahan summation for improved accuracy. More...
 
QMatrix getPureDensityMatrix (QVector state)
 Returns a density matrix initialised into the given pure state. More...
 
qcomp getRandomComplex ()
 Returns a random complex number within the square closing (-1-i) and (1+i), from a distribution uniformly randomising the individual real and imaginary components in their domains. More...
 
QMatrix getRandomDensityMatrix (int numQb)
 Returns a random numQb-by-numQb density matrix, from an undisclosed distribution, in a very mixed state. More...
 
int getRandomInt (int min, int max)
 Returns a random integer between min (inclusive) and max (exclusive), from the uniform distribution. More...
 
std::vector< QMatrixgetRandomKrausMap (int numQb, int numOps)
 Returns a random Kraus map of #numOps 2^numQb-by-2^numQb operators, from an undisclosed distribution. More...
 
std::vector< QVectorgetRandomOrthonormalVectors (int numQb, int numStates)
 Returns a list of random orthonormal complex vectors, from an undisclosed distribution. More...
 
std::vector< qrealgetRandomProbabilities (int numProbs)
 Returns a list of random real scalars, each in [0, 1], which sum to unity. More...
 
QMatrix getRandomPureDensityMatrix (int numQb)
 Returns a random numQb-by-numQb density matrix, from an undisclosed distribution, which is pure (corresponds to a random state-vector) More...
 
QMatrix getRandomQMatrix (int dim)
 Returns a dim-by-dim complex matrix, where the real and imaginary value of each element are independently random, under the standard normal distribution (mean 0, standard deviation 1). More...
 
QVector getRandomQVector (int dim)
 Returns a dim-length vector with random complex amplitudes in the square joining {-1-i, 1+i}, of an undisclosed distribution. More...
 
qreal getRandomReal (qreal min, qreal max)
 Returns a random real between min (inclusive) and max (exclusive), from the uniform distribution. More...
 
QVector getRandomStateVector (int numQb)
 Returns a random numQb-length L2-normalised state-vector from an undisclosed distribution. More...
 
QMatrix getRandomUnitary (int numQb)
 Returns a uniformly random (under Haar) 2^numQb-by-2^numQb unitary matrix. More...
 
QMatrix getSwapMatrix (int qb1, int qb2, int numQb)
 Returns the 2^numQb-by-2^numQb unitary matrix which swaps qubits qb1 and qb2; the SWAP gate of not-necessarily-adjacent qubits. More...
 
long long int getTwosComplement (long long int decimal, int numBits)
 Returns the two's complement signed encoding of the unsigned number decimal, which must be a number between 0 and 2^numBits (exclusive). More...
 
long long int getUnsigned (long long int twosComp, int numBits)
 Return the unsigned value of a number, made of #numBits bits, which under two's complement, encodes the signed number twosComp. More...
 
long long int getValueOfTargets (long long int ind, int *targs, int numTargs)
 Returns the integer value of the targeted sub-register for the given full state index ind. More...
 
QMatrix getZeroMatrix (size_t dim)
 Returns a dim-by-dim square complex matrix, initialised to all zeroes. More...
 
CatchGen< pauliOpType * > pauliseqs (int numPaulis)
 Returns a Catch2 generator of every numPaulis-length set of Pauli-matrix types (or base-4 integers). More...
 
CatchGen< int * > sequences (int base, int numDigits)
 Returns a Catch2 generator of every numDigits-length sequence in the given base, in increasing lexographic order, where left-most (zero index) bit is treated as LEAST significant (opposite typical convention). More...
 
void setDiagMatrixOverrides (QMatrix &matr, int *numQubitsPerReg, int numRegs, enum bitEncoding encoding, long long int *overrideInds, qreal *overridePhases, int numOverrides)
 Modifies the given diagonal matrix such that the diagonal elements which correspond to the coordinates in overrideInds are replaced with exp(i phase), as prescribed by overridePhases. More...
 
void setRandomDiagPauliHamil (PauliHamil hamil)
 Populates hamil with random coefficients and a random amount number of PAULI_I and PAULI_Z operators. More...
 
void setRandomPauliSum (PauliHamil hamil)
 Populates hamil with random coefficients and pauli codes. More...
 
void setRandomPauliSum (qreal *coeffs, pauliOpType *codes, int numQubits, int numTerms)
 Populates the coeffs array with random qreals in (-5, 5), and populates codes with random Pauli codes. More...
 
void setRandomTargets (int *targs, int numTargs, int numQb)
 Populates targs with a random selection of numTargs elements from [0,numQb-1]. More...
 
void setSubMatrix (QMatrix &dest, QMatrix sub, size_t r, size_t c)
 Modifies dest by overwriting its submatrix (from top-left corner (r, c) to bottom-right corner (r + dest.size(), c + dest.size()) with the complete elements of sub. More...
 
void setUniqueFilename (char *outFn, char *prefix)
 Modifies outFn to be a filename of format prefix_NUM.txt where NUM is a new unique integer so far. More...
 
CatchGen< int * > sublists (CatchGen< int > &&gen, int numSamps, const int *exclude, int numExclude)
 Returns a Catch2 generator of every length-sublen sublist of the elements generated by gen, which exclude all elements in exclude, in increasing lexographic order. More...
 
CatchGen< int * > sublists (CatchGen< int > &&gen, int numSamps, int excluded)
 Returns a Catch2 generator of every length-sublen sublist of the elements generated by gen which exclude element excluded, in increasing lexographic order. More...
 
CatchGen< int * > sublists (CatchGen< int > &&gen, int sublen)
 Returns a Catch2 generator of every length-sublen sublist of the elements generated by gen, in increasing lexographic order. More...
 
CatchGen< int * > sublists (int *list, int len, int sublen)
 Returns a Catch2 generator of every length-sublen sublist of length-len list, in increasing lexographic order. More...
 
ComplexMatrix2 toComplexMatrix2 (QMatrix qm)
 Returns a ComplexMatrix2 copy of QMatix qm. More...
 
ComplexMatrix4 toComplexMatrix4 (QMatrix qm)
 Returns a ComplexMatrix4 copy of QMatix qm. More...
 
void toComplexMatrixN (QMatrix qm, ComplexMatrixN cm)
 Initialises cm with the values of qm. More...
 
QMatrix toDiagonalQMatrix (QVector vec)
 Returns a diagonal complex matrix formed by the given vector. More...
 
QMatrix toQMatrix (Complex alpha, Complex beta)
 Returns the matrix (where a=alpha, b=beta) {{a, -conj(b)}, {b, conj(a)}} using the qcomp complex type. More...
 
QMatrix toQMatrix (ComplexMatrix2 src)
 Returns a copy of the given 2-by-2 matrix. More...
 
QMatrix toQMatrix (ComplexMatrix4 src)
 Returns a copy of the given 4-by-4 matrix. More...
 
QMatrix toQMatrix (ComplexMatrixN src)
 Returns a copy of the given 2^N-by-2^N matrix. More...
 
QMatrix toQMatrix (DiagonalOp op)
 Returns a 2^N-by-2^N complex diagonal matrix form of the DiagonalOp. More...
 
QMatrix toQMatrix (PauliHamil hamil)
 Returns a 2^N-by-2^N Hermitian matrix form of the PauliHamil. More...
 
QMatrix toQMatrix (qreal *coeffs, pauliOpType *paulis, int numQubits, int numTerms)
 Returns a 2^N-by-2^N Hermitian matrix form of the specified weighted sum of Pauli products. More...
 
QMatrix toQMatrix (Qureg qureg)
 Returns an equal-size copy of the given density matrix qureg. More...
 
QMatrix toQMatrix (SubDiagonalOp op)
 Returns a 2^n-by-2^n complex diagonal matrix form of the SubDiagonalOp, where n = op.numQubits. More...
 
void toQureg (Qureg qureg, QMatrix mat)
 Initialises the density matrix qureg to have the same amplitudes as mat. More...
 
void toQureg (Qureg qureg, QVector vec)
 Initialises the state-vector qureg to have the same amplitudes as vec. More...
 
QVector toQVector (DiagonalOp op)
 Returns a vector with the same of the full diagonal operator, populated with op's elements. More...
 
QVector toQVector (Qureg qureg)
 Returns an equal-size copy of the given state-vector qureg. More...
 
void writeToFileSynch (char *fn, const string &contents)
 Writes contents to the file with filename fn, which is created and/or overwritten. More...
 

Detailed Description

Functions used in the unit testing.

These are mostly unoptimised, analytic implementations of the complex linear algebra that QuEST ultimately effects on quantum states. These are not part of the QuEST API, and require C++14.

Author
Tyson Jones

Typedef Documentation

◆ QMatrix

typedef std::vector<std::vector<qcomp> > QMatrix

A complex square matrix.

Should be initialised with getZeroMatrix(). These have all the natural linear-algebra operator overloads, including left-multiplication onto a vector.

This data-structure is not partitioned between nodes in distributed mode. That is, every node has a complete copy, allowing for safe comparisons.

Author
Tyson Jones

◆ QVector

typedef std::vector<qcomp> QVector

A complex vector, which can be zero-initialised with QVector(numAmps).

These have all the natural linear-algebra operator overloads.

This data-structure is not partitioned between nodes in distributed mode. That is, every node has a complete copy, allowing for safe comparisons.

Author
Tyson Jones

Function Documentation

◆ applyReferenceMatrix() [1/4]

void applyReferenceMatrix ( QMatrix state,
int *  ctrls,
int  numCtrls,
int *  targs,
int  numTargs,
QMatrix  op 
)

Modifies the density matrix state to be the result of left-multiplying the multi-target operator matrix op, with the specified control and target qubits (in ctrls and targs respectively).

Here, op is treated like a simple matrix and is hence left-multiplied onto the state once.

Author
Tyson Jones

◆ applyReferenceMatrix() [2/4]

void applyReferenceMatrix ( QMatrix state,
int *  targs,
int  numTargs,
QMatrix  op 
)

Modifies the density matrix state to be the result of left-multiplying the multi-target operator matrix op, with the target qubits (assuming no control qubits).

Here, op is treated like a simple matrix and is hence left-multiplied onto the state once.

Author
Tyson Jones

◆ applyReferenceMatrix() [3/4]

void applyReferenceMatrix ( QVector state,
int *  ctrls,
int  numCtrls,
int *  targs,
int  numTargs,
QMatrix  op 
)

Modifies the state-vector state to be the result of left-multiplying the multi-target operator matrix op, with the specified control and target qubits (in ctrls and targs respectively).

This is an alias of applyReferenceOp(), since operators are always left-multiplied as matrices onto state-vectors.

Author
Tyson Jones

Referenced by applyReferenceMatrix(), and TEST_CASE().

◆ applyReferenceMatrix() [4/4]

void applyReferenceMatrix ( QVector state,
int *  targs,
int  numTargs,
QMatrix  op 
)

Modifies the state-vector state to be the result of left-multiplying the multi-target operator matrix op, with the specified target qubits (assuming no control qubits).

T

Author
Tyson Jones

◆ applyReferenceOp() [1/16]

void applyReferenceOp ( QMatrix state,
int *  ctrls,
int  numCtrls,
int *  targs,
int  numTargs,
QMatrix  op 
)

Modifies the density matrix state to be the result of applying the multi-target operator matrix op, with the specified control and target qubits (in ctrls and targs respectively).

This updates state under

\[ \text{state} \to \text{op} \, \text{state} \, \text{op}^\dagger \]

even if op is not unitary (which is useful for applying Kraus operators).

op must be a 2^numTargs-by-2^numTargs matrix. Furthermore, every element of targs must not appear in ctrls (and vice-versa), though this is not explicitly checked. Elements of targs and ctrls should be unique.

This function works by computing getFullOperatorMatrix() from the given arguments, and left-multipling it to state, then right-multiplying its conjugate transpose onto the result.

Author
Tyson Jones

◆ applyReferenceOp() [2/16]

void applyReferenceOp ( QMatrix state,
int *  ctrls,
int  numCtrls,
int  targ1,
int  targ2,
QMatrix  op 
)

Modifies the density matrix state to be the result of applying the two-target operator matrix op, with the specified control qubits (in ctrls).

This updates state under

\[ \text{state} \to \text{op} \, \text{state} \, \text{op}^\dagger \]

even if op is not unitary (which is useful for applying Kraus operators).

op must be a 4-by-4 matrix. Both targ1 and targ2 must not appear in ctrls, though this is not explicitly checked. Elements of ctrls, and targ1 and targ2, should be unique.

This function works by computing getFullOperatorMatrix() from the given arguments, and left-multipling it to state, then right-multiplying its conjugate transpose onto the result.

Author
Tyson Jones

◆ applyReferenceOp() [3/16]

void applyReferenceOp ( QMatrix state,
int *  ctrls,
int  numCtrls,
int  target,
QMatrix  op 
)

Modifies the density matrix state to be the result of applying the single-target operator matrix op, with the specified control qubits (in ctrls).

This updates state under

\[ \text{state} \to \text{op} \, \text{state} \, \text{op}^\dagger \]

even if op is not unitary (which is useful for applying Kraus operators).

op must be a 2-by-2 matrix. target must not appear in ctrls, though this is not explicitly checked.

This function works by computing getFullOperatorMatrix() from the given arguments, and left-multipling it to state, then right-multiplying its conjugate transpose onto the result.

Author
Tyson Jones

◆ applyReferenceOp() [4/16]

void applyReferenceOp ( QMatrix state,
int *  targs,
int  numTargs,
QMatrix  op 
)

Modifies the density matrix state to be the result of applying the multi-target operator matrix op, with no control qubits.

This updates state under

\[ \text{state} \to \text{op} \, \text{state} \, \text{op}^\dagger \]

even if op is not unitary (which is useful for applying Kraus operators).

op must be a 2^numTargs-by-2^numTargs matrix. Every element in targs should be unique, though this is not explicitly checked.

This function works by computing getFullOperatorMatrix() from the given arguments, and left-multipling it to state, then right-multiplying its conjugate transpose onto the result.

Author
Tyson Jones

◆ applyReferenceOp() [5/16]

void applyReferenceOp ( QMatrix state,
int  ctrl,
int *  targs,
int  numTargs,
QMatrix  op 
)

Modifies the density matrix state to be the result of applying the multi-target operator matrix op, with a single control qubit ctrl.

This updates state under

\[ \text{state} \to \text{op} \, \text{state} \, \text{op}^\dagger \]

even if op is not unitary (which is useful for applying Kraus operators).

op must be a 2^numTargs-by-2^numTargs matrix, and ctrl must not appear in targs (though this is not explicitly checked).

This function works by computing getFullOperatorMatrix() from the given arguments, and left-multipling it to state, then right-multiplying its conjugate transpose onto the result.

Author
Tyson Jones

◆ applyReferenceOp() [6/16]

void applyReferenceOp ( QMatrix state,
int  ctrl,
int  targ,
QMatrix  op 
)

Modifies the density matrix state to be the result of applying the single-control single-target operator matrix op.

This updates state under

\[ \text{state} \to \text{op} \, \text{state} \, \text{op}^\dagger \]

even if op is not unitary (which is useful for applying Kraus operators).

op must be a 2-by-2 matrix, and ctrl and targ should be different.

This function works by computing getFullOperatorMatrix() from the given arguments, and left-multipling it to state, then right-multiplying its conjugate transpose onto the result.

Author
Tyson Jones

◆ applyReferenceOp() [7/16]

void applyReferenceOp ( QMatrix state,
int  ctrl,
int  targ1,
int  targ2,
QMatrix  op 
)

Modifies the density matrix state to be the result of applying the two-target operator matrix op, with a single control qubit ctrl.

This updates state under

\[ \text{state} \to \text{op} \, \text{state} \, \text{op}^\dagger \]

even if op is not unitary (which is useful for applying Kraus operators).

op must be a 4-by-4 matrix, and ctrl, targ1 and targ2 must be unique.

This function works by computing getFullOperatorMatrix() from the given arguments, and left-multipling it to state, then right-multiplying its conjugate transpose onto the result.

Author
Tyson Jones

◆ applyReferenceOp() [8/16]

void applyReferenceOp ( QMatrix state,
int  targ,
QMatrix  op 
)

Modifies the density matrix state to be the result of applying the single-target operator matrix op, with no control qubit.

This updates state under

\[ \text{state} \to \text{op} \, \text{state} \, \text{op}^\dagger \]

even if op is not unitary (which is useful for applying Kraus operators).

op must be a 2-by-2 matrix.

This function works by computing getFullOperatorMatrix() from the given arguments, and left-multipling it to state, then right-multiplying its conjugate transpose onto the result.

Author
Tyson Jones

◆ applyReferenceOp() [9/16]

void applyReferenceOp ( QVector state,
int *  ctrls,
int  numCtrls,
int *  targs,
int  numTargs,
QMatrix  op 
)

Modifies the state-vector state to be the result of applying the multi-target operator matrix op, with the specified control and target qubits (in ctrls and targs respectively).

This updates state under

\[ \text{state} \to \text{op} \, \text{state} \]

even if op is not unitary.

op must be a 2^numTargs-by-2^numTargs matrix. Furthermore, every element of targs must not appear in ctrls (and vice-versa), though this is not explicitly checked. Elements of targs and ctrls should be unique.

This function works by computing getFullOperatorMatrix() from the given arguments, and left-multiplying it onto state.

Author
Tyson Jones

Referenced by applyReferenceMatrix(), applyReferenceOp(), and TEST_CASE().

◆ applyReferenceOp() [10/16]

void applyReferenceOp ( QVector state,
int *  ctrls,
int  numCtrls,
int  targ1,
int  targ2,
QMatrix  op 
)

Modifies the state-vector state to be the result of applying the two-target operator matrix op, with the specified control qubits (in ctrls).

This updates state under

\[ \text{state} \to \text{op} \, \text{state} \]

even if op is not unitary.

op must be a 4-by-4 matrix. Furthermore, ctrls, targ1 and targ2 should all be unique.

This function works by computing getFullOperatorMatrix() from the given arguments, and left-multiplying it onto state.

Author
Tyson Jones

◆ applyReferenceOp() [11/16]

void applyReferenceOp ( QVector state,
int *  ctrls,
int  numCtrls,
int  target,
QMatrix  op 
)

Modifies the state-vector state to be the result of applying the single-target operator matrix op, with the specified control qubits (in ctrls).

This updates state under

\[ \text{state} \to \text{op} \, \text{state} \]

even if op is not unitary.

op must be a 2-by-2 matrix. Furthermore, elements in ctrls and target should all be unique.

This function works by computing getFullOperatorMatrix() from the given arguments, and left-multiplying it onto state.

Author
Tyson Jones

◆ applyReferenceOp() [12/16]

void applyReferenceOp ( QVector state,
int *  targs,
int  numTargs,
QMatrix  op 
)

Modifies the state-vector state to be the result of applying the multi-target operator matrix op, with no contorl qubits.

This updates state under

\[ \text{state} \to \text{op} \, \text{state} \]

even if op is not unitary.

op must be a 2^numTargs-by-2^numTargs matrix. Furthermore, elements in targs should be unique.

This function works by computing getFullOperatorMatrix() from the given arguments, and left-multiplying it onto state.

Author
Tyson Jones

◆ applyReferenceOp() [13/16]

void applyReferenceOp ( QVector state,
int  ctrl,
int *  targs,
int  numTargs,
QMatrix  op 
)

Modifies the state-vector state to be the result of applying the multi-target operator matrix op, with a single control qubit (ctrl) This updates state under.

\[ \text{state} \to \text{op} \, \text{state} \]

even if op is not unitary.

op must be a 2^numTargs-by-2^numTargs matrix. Furthermore, elements in targs and ctrl should be unique.

This function works by computing getFullOperatorMatrix() from the given arguments, and left-multiplying it onto state.

Author
Tyson Jones

◆ applyReferenceOp() [14/16]

void applyReferenceOp ( QVector state,
int  ctrl,
int  targ,
QMatrix  op 
)

Modifies the state-vector state to be the result of applying the single-target operator matrix op, with a single control qubit (ctrl).

This updates state under

\[ \text{state} \to \text{op} \, \text{state} \]

even if op is not unitary.

op must be a 2-by-2 matrix. Furthermore, ctrl and targ must be different.

This function works by computing getFullOperatorMatrix() from the given arguments, and left-multiplying it onto state.

Author
Tyson Jones

◆ applyReferenceOp() [15/16]

void applyReferenceOp ( QVector state,
int  ctrl,
int  targ1,
int  targ2,
QMatrix  op 
)

Modifies the state-vector state to be the result of applying the two-target operator matrix op, with a single control qubit (ctrl).

This updates state under

\[ \text{state} \to \text{op} \, \text{state} \]

even if op is not unitary.

op must be a 4-by-4 matrix. Furthermore, ctrl, targ1 and targ2 should all be unique.

This function works by computing getFullOperatorMatrix() from the given arguments, and left-multiplying it onto state.

Author
Tyson Jones

◆ applyReferenceOp() [16/16]

void applyReferenceOp ( QVector state,
int  targ,
QMatrix  op 
)

Modifies the state-vector state to be the result of applying the single-target operator matrix op, with no contorl qubits.

This updates state under

\[ \text{state} \to \text{op} \, \text{state} \]

even if op is not unitary.

op must be a 2-by-2 matrix.

This function works by computing getFullOperatorMatrix() from the given arguments, and left-multiplying it onto state.

Author
Tyson Jones

◆ areEqual() [1/10]

bool areEqual ( QMatrix  a,
QMatrix  b 
)

Returns true if the absolute value of the difference between every amplitude in matrices a and b is less than REAL_EPS.

Author
Tyson Jones

◆ areEqual() [2/10]

bool areEqual ( Qureg  qureg,
QMatrix  matr 
)

Performs a hardware-agnostic comparison of density-matrix qureg to matr, checking whether the difference between the real and imaginary components of every amplitude is smaller than the QuEST_PREC-specific REAL_EPS (defined in QuEST_precision) precision.

This function demands qureg is a density matrix, and that qureg and matr have equal dimensions.

In GPU mode, this function involves a GPU to CPU memory copy overhead. In distributed mode, it involves a all-to-all single-int broadcast.

Author
Tyson Jones

◆ areEqual() [3/10]

bool areEqual ( Qureg  qureg,
QMatrix  matr,
qreal  precision 
)

Performs a hardware-agnostic comparison of density-matrix qureg to matr, checking whether the difference between the real and imaginary components of every amplitude is smaller than precision.

This function demands qureg is a density matrix, and that qureg and matr have equal dimensions.

In GPU mode, this function involves a GPU to CPU memory copy overhead. In distributed mode, it involves a all-to-all single-int broadcast.

Author
Tyson Jones

◆ areEqual() [4/10]

bool areEqual ( Qureg  qureg,
QVector  vec 
)

Performs a hardware-agnostic comparison of state-vector qureg to vec, checking whether the difference between the real and imaginary components of every amplitude is smaller than the QuEST_PREC-specific REAL_EPS (defined in QuEST_precision) precision.

This function demands qureg is a state-vector, and that qureg and vec have the same number of amplitudes.

In GPU mode, this function involves a GPU to CPU memory copy overhead. In distributed mode, it involves a all-to-all single-int broadcast.

Author
Tyson Jones

◆ areEqual() [5/10]

bool areEqual ( Qureg  qureg,
QVector  vec,
qreal  precision 
)

Performs a hardware-agnostic comparison of state-vector qureg to vec, checking whether the difference between the real and imaginary components of every amplitude is smaller than precision.

This function demands qureg is a state-vector, and that qureg and vec have the same number of amplitudes.

In GPU mode, this function involves a GPU to CPU memory copy overhead. In distributed mode, it involves a all-to-all single-int broadcast.

Author
Tyson Jones

◆ areEqual() [6/10]

bool areEqual ( Qureg  qureg1,
Qureg  qureg2 
)

Performs a hardware-agnostic comparison of the given quregs, checking whether the difference between the real and imaginary components of every amplitude is smaller than the QuEST_PREC-specific REAL_EPS (defined in QuEST_precision) precision.

This function demands that qureg1 and qureg2 are of the same type (i.e. both state-vectors or both density matrices), and of an equal number of qubits.

In GPU mode, this function involves a GPU to CPU memory copy overhead. In distributed mode, it involves a all-to-all single-int broadcast.

Author
Tyson Jones

◆ areEqual() [7/10]

bool areEqual ( Qureg  qureg1,
Qureg  qureg2,
qreal  precision 
)

Performs a hardware-agnostic comparison of the given quregs, checking whether the difference between the real and imaginary components of every amplitude is smaller than precision.

This function demands that qureg1 and qureg2 are of the same type (i.e. both state-vectors or both density matrices), and of an equal number of qubits.

In GPU mode, this function involves a GPU to CPU memory copy overhead. In distributed mode, it involves a all-to-all single-int broadcast.

Author
Tyson Jones

◆ areEqual() [8/10]

bool areEqual ( QVector  a,
QVector  b 
)

Returns true if the absolute value of the difference between every amplitude in vectors a and b is less than REAL_EPS.

Author
Tyson Jones

Referenced by areEqual(), assertQuregAndRefInDebugState(), getRandomKrausMap(), getRandomUnitary(), and TEST_CASE().

◆ areEqual() [9/10]

bool areEqual ( QVector  vec,
qreal reals 
)

Returns true if the absolute value of the difference between every element in vec (which must be strictly real) and those implied by reals, is less than REAL_EPS.

Author
Tyson Jones

◆ areEqual() [10/10]

bool areEqual ( QVector  vec,
qreal reals,
qreal imags 
)

Returns true if the absolute value of the difference between every element in vec and those implied by reals and imags, is less than REAL_EPS.

Author
Tyson Jones

◆ assertQuregAndRefInDebugState() [1/2]

void assertQuregAndRefInDebugState ( Qureg  qureg,
QMatrix  ref 
)

Asserts the given density qureg and reference agree, and are properly initialised in the debug state.

Assertion uses the DEMAND() macro, calling Catch2's FAIL() if unsatisfied, so does not contribute toward unit test statistics. This should be called within every PREPARE_TEST macro, to ensure that the test states themselves are initially correct, and do not accidentally agree by (e.g.) being all-zero.

Author
Tyson Jones

◆ assertQuregAndRefInDebugState() [2/2]

void assertQuregAndRefInDebugState ( Qureg  qureg,
QVector  ref 
)

Asserts the given statevector qureg and reference agree, and are properly initialised in the debug state.

Assertion uses the DEMAND() macro, calling Catch2's FAIL() if unsatisfied, so does not contribute toward unit test statistics. This should be called within every PREPARE_TEST macro, to ensure that the test states themselves are initially correct, and do not accidentally agree by (e.g.) being all-zero.

Author
Tyson Jones

◆ bitsets()

CatchGen< int * > bitsets ( int  numBits)

Returns a Catch2 generator of every numBits-length bit-set, in increasing lexographic order, where left-most (zero index) bit is treated as LEAST significant (opposite typical convention).

Note that the produced bitset must not be modified during generation.

This function can be used like

int* bits = GENERATE( bitsets(3) );

to produce {0,0,0}, {1,0,0}, {0,1,0}, {1,1,0}, {0,0,1}, {1,0,1}, {0,1,1}, {1,1,1}.

Author
Tyson Jones

Referenced by TEST_CASE().

◆ calcLog2()

unsigned int calcLog2 ( long unsigned int  num)

Returns log2 of numbers which must be gauranteed to be 2^n.

Author
Tyson Jones

Returns log2 of numbers which must be gauranteed to be 2^n.

Referenced by applyReferenceMatrix(), applyReferenceOp(), and TEST_CASE().

◆ deleteFilesWithPrefixSynch()

void deleteFilesWithPrefixSynch ( char *  prefix)

Deletes all files with filename starting with prefix.

In distributed mode, the master node deletes while the other nodes wait until complete.

Author
Tyson Jones

Referenced by TEST_CASE().

◆ getConjugateTranspose()

QMatrix getConjugateTranspose ( QMatrix  a)

Returns the conjugate transpose of the complex square matrix a.

Author
Tyson Jones

Referenced by applyReferenceOp(), getRandomKrausMap(), and getRandomUnitary().

◆ getDFT() [1/2]

QVector getDFT ( QVector  in)

Returns the discrete fourier transform of vector in.

Author
Tyson Jones

Referenced by TEST_CASE().

◆ getDFT() [2/2]

QVector getDFT ( QVector  in,
int *  targs,
int  numTargs 
)

Returns the discrete fourier transform of a sub-partition of the vector in.

Author
Tyson Jones

◆ getExponentialOfDiagonalMatrix()

QMatrix getExponentialOfDiagonalMatrix ( QMatrix  a)

Returns the matrix exponential of a diagonal, square, complex matrix.

This method explicitly checks that the passed matrix a is diagonal.

Author
Tyson Jones

Referenced by TEST_CASE().

◆ getExponentialOfPauliMatrix()

QMatrix getExponentialOfPauliMatrix ( qreal  angle,
QMatrix  a 
)

Returns the matrix exponential of a kronecker product of pauli matrices (or of any involutory matrices), with exponent factor (-i angle / 2).

This method will not explicitly check that the passed matrix a is kronecker product of involutory matrices, but will otherwise return an incorrect exponential.

Author
Tyson Jones

Referenced by TEST_CASE().

◆ getFullOperatorMatrix()

QMatrix getFullOperatorMatrix ( int *  ctrls,
int  numCtrls,
int *  targs,
int  numTargs,
QMatrix  op,
int  numQubits 
)

Takes a 2^numTargs-by-2^numTargs matrix op and a returns a 2^numQubits-by-2^numQubits matrix where op is controlled on the given ctrls qubits.

The union of {ctrls} and {targs} must be unique (though this is not explicitly checked), and every element must be >= 0 (not checked). The passed {ctrls} and {targs} arrays are unmodified.

This funciton works by first swapping {targs} and {ctrls} (via swap unitaries) to be strictly increasing {0,1,...}, building controlled(op), tensoring it to the full Hilbert space, and then 'unswapping'. The returned matrix has form: swap1 ... swapN . controlled(op) . swapN ... swap1

Author
Tyson Jones

Referenced by applyReferenceMatrix(), applyReferenceOp(), and TEST_CASE().

◆ getIdentityMatrix()

QMatrix getIdentityMatrix ( size_t  dim)

Returns a dim-by-dim identity matrix.

Author
Tyson Jones

Referenced by getExponentialOfPauliMatrix(), getFullOperatorMatrix(), getRandomKrausMap(), getRandomUnitary(), and getSwapMatrix().

◆ getKetBra()

QMatrix getKetBra ( QVector  ket,
QVector  bra 
)

Returns the matrix |ket><bra|, with ith-jth element ket(i) conj(bra(j)), since |ket><bra| = sum_i a_i|i> sum_j b_j* <j| = sum_{ij} a_i b_j* |i><j|.

The dimensions of bra and ket must agree, and the returned square complex matrix has dimensions size(bra) x size(bra).

Author
Tyson Jones

Referenced by getPureDensityMatrix(), getRandomDensityMatrix(), and TEST_CASE().

◆ getKroneckerProduct() [1/2]

QMatrix getKroneckerProduct ( QMatrix  a,
QMatrix  b 
)

Returns the kronecker product of a and b, where a and b are square but possibly differently-sized complex matrices.

Author
Tyson Jones

◆ getKroneckerProduct() [2/2]

QVector getKroneckerProduct ( QVector  b,
QVector  a 
)

Returns b (otimes) a.

If b and a are state-vectors, the resulting kronecker product is the seperable state formed by joining the qubits in the state-vectors, producing |b>|a> (a is least significant)

Author
Tyson Jones

Referenced by getFullOperatorMatrix(), getSwapMatrix(), TEST_CASE(), and toQMatrix().

◆ getMatrixDiagonal()

QVector getMatrixDiagonal ( QMatrix  matr)

Returns the diagonal vector of the given matrix.

Author
Tyson Jones

◆ getMixedDensityMatrix()

QMatrix getMixedDensityMatrix ( std::vector< qreal probs,
std::vector< QVector states 
)

Returns a mixed density matrix formed from mixing the given pure states, which are assumed normalised, but not necessarily orthogonal.

Author
Tyson Jones

Referenced by TEST_CASE().

◆ getNormalised()

QVector getNormalised ( QVector  vec)

Returns an L2-normalised copy of vec, using Kahan summation for improved accuracy.

Author
Tyson Jones

Referenced by getRandomOrthonormalVectors(), and getRandomStateVector().

◆ getPureDensityMatrix()

QMatrix getPureDensityMatrix ( QVector  state)

Returns a density matrix initialised into the given pure state.

Author
Tyson Jones

Referenced by getMixedDensityMatrix(), getRandomPureDensityMatrix(), and TEST_CASE().

◆ getRandomComplex()

qcomp getRandomComplex ( )

Returns a random complex number within the square closing (-1-i) and (1+i), from a distribution uniformly randomising the individual real and imaginary components in their domains.

Author
Tyson Jones

Referenced by getRandomQVector(), and TEST_CASE().

◆ getRandomDensityMatrix()

QMatrix getRandomDensityMatrix ( int  numQb)

Returns a random numQb-by-numQb density matrix, from an undisclosed distribution, in a very mixed state.

This function works by generating 2^numQb random pure states, and mixing them with random probabilities.

Author
Tyson Jones

Referenced by TEST_CASE().

◆ getRandomInt()

int getRandomInt ( int  min,
int  max 
)

Returns a random integer between min (inclusive) and max (exclusive), from the uniform distribution.

Demands that max > min.

Author
Tyson Jones

Referenced by setRandomPauliSum(), and TEST_CASE().

◆ getRandomKrausMap()

std::vector< QMatrix > getRandomKrausMap ( int  numQb,
int  numOps 
)

Returns a random Kraus map of #numOps 2^numQb-by-2^numQb operators, from an undisclosed distribution.

Note this method is very simple and cannot generate all possible Kraus maps. It works by generating numOps random unitary matrices, and randomly re-normalising them, such that the sum of ops[j]^dagger ops[j] = 1

Author
Tyson Jones

Referenced by TEST_CASE().

◆ getRandomOrthonormalVectors()

std::vector< QVector > getRandomOrthonormalVectors ( int  numQb,
int  numStates 
)

Returns a list of random orthonormal complex vectors, from an undisclosed distribution.

Author
Tyson Jones

Referenced by TEST_CASE().

◆ getRandomProbabilities()

std::vector< qreal > getRandomProbabilities ( int  numProbs)

Returns a list of random real scalars, each in [0, 1], which sum to unity.

Author
Tyson Jones

Referenced by getRandomDensityMatrix(), and TEST_CASE().

◆ getRandomPureDensityMatrix()

QMatrix getRandomPureDensityMatrix ( int  numQb)

Returns a random numQb-by-numQb density matrix, from an undisclosed distribution, which is pure (corresponds to a random state-vector)

Author
Tyson Jones

◆ getRandomQMatrix()

QMatrix getRandomQMatrix ( int  dim)

Returns a dim-by-dim complex matrix, where the real and imaginary value of each element are independently random, under the standard normal distribution (mean 0, standard deviation 1).

Author
Tyson Jones

Referenced by getRandomUnitary(), and TEST_CASE().

◆ getRandomQVector()

QVector getRandomQVector ( int  dim)

Returns a dim-length vector with random complex amplitudes in the square joining {-1-i, 1+i}, of an undisclosed distribution.

The resulting vector is NOT L2-normalised.

Author
Tyson Jones

Referenced by getRandomStateVector(), and TEST_CASE().

◆ getRandomReal()

qreal getRandomReal ( qreal  min,
qreal  max 
)

Returns a random real between min (inclusive) and max (exclusive), from the uniform distribution.

Demands that max > min.

Author
Tyson Jones

Referenced by getRandomComplex(), getRandomInt(), getRandomKrausMap(), getRandomProbabilities(), setRandomDiagPauliHamil(), setRandomPauliSum(), and TEST_CASE().

◆ getRandomStateVector()

QVector getRandomStateVector ( int  numQb)

Returns a random numQb-length L2-normalised state-vector from an undisclosed distribution.

This function works by randomly generating each complex amplitude, then L2-normalising.

Author
Tyson Jones

Referenced by getRandomDensityMatrix(), getRandomOrthonormalVectors(), getRandomPureDensityMatrix(), and TEST_CASE().

◆ getRandomUnitary()

QMatrix getRandomUnitary ( int  numQb)

Returns a uniformly random (under Haar) 2^numQb-by-2^numQb unitary matrix.

This function works by first generating a complex matrix where each element is independently random; the real and imaginary component thereof are independent standard normally-distributed (mean 0, standard-dev 1). Then, the matrix is orthonormalised via the Gram Schmidt algorithm. The resulting unitary matrix MAY be uniformly distributed under the Haar measure, but we make no assurance. This routine may return an identity matrix if it was unable to sufficiently precisely produce a unitary of the given size.

Author
Tyson Jones

Referenced by getRandomKrausMap(), and TEST_CASE().

◆ getSwapMatrix()

QMatrix getSwapMatrix ( int  qb1,
int  qb2,
int  numQb 
)

Returns the 2^numQb-by-2^numQb unitary matrix which swaps qubits qb1 and qb2; the SWAP gate of not-necessarily-adjacent qubits.

If qb1 == qb2, returns the identity matrix.

Author
Tyson Jones

Referenced by getFullOperatorMatrix().

◆ getTwosComplement()

long long int getTwosComplement ( long long int  decimal,
int  numBits 
)

Returns the two's complement signed encoding of the unsigned number decimal, which must be a number between 0 and 2^numBits (exclusive).

The returned number lies in [-2^(numBits-1), 2^(numBits-1)-1]

Author
Tyson Jones

Referenced by TEST_CASE().

◆ getUnsigned()

long long int getUnsigned ( long long int  twosComp,
int  numBits 
)

Return the unsigned value of a number, made of #numBits bits, which under two's complement, encodes the signed number twosComp.

The returned number lies in [0, 2^(numBits)-1]

Author
Tyson Jones

Referenced by setDiagMatrixOverrides().

◆ getValueOfTargets()

long long int getValueOfTargets ( long long int  ind,
int *  targs,
int  numTargs 
)

Returns the integer value of the targeted sub-register for the given full state index ind.

Author
Tyson Jones

Referenced by getDFT().

◆ getZeroMatrix()

QMatrix getZeroMatrix ( size_t  dim)

◆ pauliseqs()

CatchGen< pauliOpType * > pauliseqs ( int  numPaulis)

Returns a Catch2 generator of every numPaulis-length set of Pauli-matrix types (or base-4 integers).

Generates in increasing lexographic order, where the left-most (zero index) pauli is treated as LEAST significant. Note that the sequence must not be modified during generation.

This function can be used like

pauliOpType* set = GENERATE( pauliseqs(2) );

to produce {I,I}, {X,I}, {Y,I}, {Z,I}, {I,X}, {X,X}, {Y,X}, {Z,X}, {I,Y}, {X,Y}, {Y,Y}, {Z,Y}, {I,Z}, {X,Z}, {Y,Z}, {Z,Z}/

Author
Tyson Jones

◆ sequences()

CatchGen< int * > sequences ( int  base,
int  numDigits 
)

Returns a Catch2 generator of every numDigits-length sequence in the given base, in increasing lexographic order, where left-most (zero index) bit is treated as LEAST significant (opposite typical convention).

Note that the sequence must not be modified during generation.

This function can be used like

int base = 3;
int numDigits = 2;
int* seq = GENERATE_COPY( sequences(base, numDigits) );

to produce {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {2,1}, {0,2}, {1,2}, {2,2}.

Author
Tyson Jones

◆ setDiagMatrixOverrides()

void setDiagMatrixOverrides ( QMatrix matr,
int *  numQubitsPerReg,
int  numRegs,
enum bitEncoding  encoding,
long long int *  overrideInds,
qreal overridePhases,
int  numOverrides 
)

Modifies the given diagonal matrix such that the diagonal elements which correspond to the coordinates in overrideInds are replaced with exp(i phase), as prescribed by overridePhases.

This function assumes that the given registers are contiguous, are in order of increasing significance, and that the matrix is proportionately sized and structured to act on the space of all registers combined. Overrides can be repeated, and only the first encountered for a given index will be effected (much like applyMultiVarPhaseFuncOverrides()).

Author
Tyson Jones

Referenced by TEST_CASE().

◆ setRandomDiagPauliHamil()

void setRandomDiagPauliHamil ( PauliHamil  hamil)

Populates hamil with random coefficients and a random amount number of PAULI_I and PAULI_Z operators.

Author
Tyson Jones

Referenced by TEST_CASE().

◆ setRandomPauliSum() [1/2]

void setRandomPauliSum ( PauliHamil  hamil)

Populates hamil with random coefficients and pauli codes.

Author
Tyson Jones

◆ setRandomPauliSum() [2/2]

void setRandomPauliSum ( qreal coeffs,
pauliOpType codes,
int  numQubits,
int  numTerms 
)

Populates the coeffs array with random qreals in (-5, 5), and populates codes with random Pauli codes.

Author
Tyson Jones

Referenced by setRandomPauliSum(), and TEST_CASE().

◆ setRandomTargets()

void setRandomTargets ( int *  targs,
int  numTargs,
int  numQb 
)

Populates targs with a random selection of numTargs elements from [0,numQb-1].

List targs does not need to be initialised and its elements are overwritten.

Author
Tyson Jones

Referenced by TEST_CASE().

◆ setSubMatrix()

void setSubMatrix ( QMatrix dest,
QMatrix  sub,
size_t  r,
size_t  c 
)

Modifies dest by overwriting its submatrix (from top-left corner (r, c) to bottom-right corner (r + dest.size(), c + dest.size()) with the complete elements of sub.

This demands that dest.size() >= sub.size() + max(r,c).

Author
Tyson Jones

Referenced by getFullOperatorMatrix(), and getSwapMatrix().

◆ setUniqueFilename()

void setUniqueFilename ( char *  outFn,
char *  prefix 
)

Modifies outFn to be a filename of format prefix_NUM.txt where NUM is a new unique integer so far.

This is useful for getting unique filenames for independent test cases of functions requiring reading/writing to file, to avoid IO locks (especially common in distributed mode).

Author
Tyson Jones

Referenced by TEST_CASE().

◆ sublists() [1/4]

CatchGen< int * > sublists ( CatchGen< int > &&  gen,
int  numSamps,
const int *  exclude,
int  numExclude 
)

Returns a Catch2 generator of every length-sublen sublist of the elements generated by gen, which exclude all elements in exclude, in increasing lexographic order.

This generates every fixed-length combination of gen's elements the nominated exclusions, and every permutation of each.

There is on need for the elements of exclude to actually appear in those of gen. sublen must less than or equal to the number of elements in gen, after the nominated exclusions.

Note that the sublist must not be modified, else further generation may break (QuEST's internal functions will indeed modify but restore the qubit index lists given to them, which is ok). Assumes list contains no duplicates, otherwise the generated sublists may be duplicated.

This function can be used like

int sublen = 2;
int exclude[2] = {3,4};
int* sublist = GENERATE_COPY( sublists(range(1,6), sublen, exclude, 2) );

to generate {1,2}, {1,5}, {2,1}, {2,5}, {5,1}, {5,2}

Author
Tyson Jones

◆ sublists() [2/4]

CatchGen< int * > sublists ( CatchGen< int > &&  gen,
int  numSamps,
int  excluded 
)

Returns a Catch2 generator of every length-sublen sublist of the elements generated by gen which exclude element excluded, in increasing lexographic order.

This generates every fixed-length combination of gen's elements the nominated exclusions, and every permutation of each.

sublen must less than or equal to the number of elements in gen, after the nominated exclusion. There is no need for excluded to actually appear in the elements of gen.

Note that the sublist must not be modified, else further generation may break (QuEST's internal functions will indeed modify but restore the qubit index lists given to them, which is ok). Assumes list contains no duplicates, otherwise the generated sublists may be duplicated.

This function can be used like

int sublen = 2;
int excluded = 1;
int* sublist = GENERATE_COPY( sublists(range(1,4), sublen, excluded) );

to generate {2,3}, {3,2}.

Author
Tyson Jones

◆ sublists() [3/4]

CatchGen< int * > sublists ( CatchGen< int > &&  gen,
int  sublen 
)

Returns a Catch2 generator of every length-sublen sublist of the elements generated by gen, in increasing lexographic order.

This generates every fixed-length combination of gen's elements, and every permutation of each. Note that the produced sublist must not be modified, else further generation may break (QuEST's internal functions will indeed modify but restore the qubit index lists given to them, which is ok). Assumes list contains no duplicates, otherwise the generated sublists may be duplicated.

This function can be used like

int sublen = 2;
int* sublist = GENERATE_COPY( sublists(list, 4, sublen) );

to generate {1,2}, {1,3}, {1,4}, {2,1}, {2,3}, {2,4}, {3,1}, {3,2}, {3, 4}, {4,1}, {4,2}, {4, 3}.

Author
Tyson Jones

◆ sublists() [4/4]

CatchGen< int * > sublists ( int *  list,
int  len,
int  sublen 
)

Returns a Catch2 generator of every length-sublen sublist of length-len list, in increasing lexographic order.

This generates every fixed-length combination of the given list and every permutation of each. & If the sublist length is the full list length, this generator produces every permutation correctly. Note that the sublist must not be modified, else further & generation may break (QuEST's internal functions will indeed modify but restore the qubit index lists given to them, which is ok). Assumes list contains no duplicates, otherwise the generated sublists may be duplicated.

This function can be used like

int list[4] = {1,2,3,4};
int sublen = 2;
int* sublist = GENERATE_COPY( sublists(list, 4, sublen) );

to generate {1,2}, {1,3}, {1,4}, {2,1}, {2,3}, {2,4}, {3,1}, {3,2}, {3, 4}, {4,1}, {4,2}, {4, 3}.

Author
Tyson Jones

Referenced by TEST_CASE().

◆ toComplexMatrix2()

ComplexMatrix2 toComplexMatrix2 ( QMatrix  qm)

Returns a ComplexMatrix2 copy of QMatix qm.

Demands that qm is a 2-by-2 matrix.

Author
Tyson Jones

Referenced by TEST_CASE().

◆ toComplexMatrix4()

ComplexMatrix4 toComplexMatrix4 ( QMatrix  qm)

Returns a ComplexMatrix4 copy of QMatix qm.

Demands that qm is a 4-by-4 matrix.

Author
Tyson Jones

Referenced by TEST_CASE().

◆ toComplexMatrixN()

void toComplexMatrixN ( QMatrix  qm,
ComplexMatrixN  cm 
)

Initialises cm with the values of qm.

Demands that cm is a previously created ComplexMatrixN instance, with the same dimensions as qm.

Author
Tyson Jones

Referenced by TEST_CASE().

◆ toDiagonalQMatrix()

QMatrix toDiagonalQMatrix ( QVector  vec)

Returns a diagonal complex matrix formed by the given vector.

Author
Tyson Jones

◆ toQMatrix() [1/9]

QMatrix toQMatrix ( Complex  alpha,
Complex  beta 
)

Returns the matrix (where a=alpha, b=beta) {{a, -conj(b)}, {b, conj(a)}} using the qcomp complex type.

Author
Tyson Jones

◆ toQMatrix() [2/9]

QMatrix toQMatrix ( ComplexMatrix2  src)

Returns a copy of the given 2-by-2 matrix.

Author
Tyson Jones

Referenced by TEST_CASE(), and toQMatrix().

◆ toQMatrix() [3/9]

QMatrix toQMatrix ( ComplexMatrix4  src)

Returns a copy of the given 4-by-4 matrix.

Author
Tyson Jones

◆ toQMatrix() [4/9]

QMatrix toQMatrix ( ComplexMatrixN  src)

Returns a copy of the given 2^N-by-2^N matrix.

Author
Tyson Jones

◆ toQMatrix() [5/9]

QMatrix toQMatrix ( DiagonalOp  op)

Returns a 2^N-by-2^N complex diagonal matrix form of the DiagonalOp.

Author
Tyson Jones

◆ toQMatrix() [6/9]

QMatrix toQMatrix ( PauliHamil  hamil)

Returns a 2^N-by-2^N Hermitian matrix form of the PauliHamil.

Author
Tyson Jones

◆ toQMatrix() [7/9]

QMatrix toQMatrix ( qreal coeffs,
pauliOpType paulis,
int  numQubits,
int  numTerms 
)

Returns a 2^N-by-2^N Hermitian matrix form of the specified weighted sum of Pauli products.

Author
Tyson Jones

◆ toQMatrix() [8/9]

QMatrix toQMatrix ( Qureg  qureg)

Returns an equal-size copy of the given density matrix qureg.

In GPU mode, this function involves a copy of qureg from GPU memory to RAM. In distributed mode, this involves an all-to-all broadcast of qureg.

Author
Tyson Jones

◆ toQMatrix() [9/9]

QMatrix toQMatrix ( SubDiagonalOp  op)

Returns a 2^n-by-2^n complex diagonal matrix form of the SubDiagonalOp, where n = op.numQubits.

Author
Tyson Jones

◆ toQureg() [1/2]

void toQureg ( Qureg  qureg,
QMatrix  mat 
)

Initialises the density matrix qureg to have the same amplitudes as mat.

Demands qureg is a density matrix of equal dimensions to mat. In GPU mode, this function involves a copy from RAM to GPU memory. This function has no communication cost in distributed mode.

Author
Tyson Jones

◆ toQureg() [2/2]

void toQureg ( Qureg  qureg,
QVector  vec 
)

Initialises the state-vector qureg to have the same amplitudes as vec.

Demands qureg is a state-vector of an equal size to vec. In GPU mode, this function involves a copy from RAM to GPU memory. This function has no communication cost in distributed mode.

Author
Tyson Jones

Referenced by TEST_CASE().

◆ toQVector() [1/2]

QVector toQVector ( DiagonalOp  op)

Returns a vector with the same of the full diagonal operator, populated with op's elements.

In distributed mode, this involves an all-to-all broadcast of op.

Author
Tyson Jones

◆ toQVector() [2/2]

QVector toQVector ( Qureg  qureg)

Returns an equal-size copy of the given state-vector qureg.

In GPU mode, this function involves a copy of qureg from GPU memory to RAM. In distributed mode, this involves an all-to-all broadcast of qureg.

Author
Tyson Jones

Referenced by TEST_CASE(), and toQMatrix().

◆ writeToFileSynch()

void writeToFileSynch ( char *  fn,
const string &  contents 
)

Writes contents to the file with filename fn, which is created and/or overwritten.

In distributed mode, the master node writes while the other nodes wait until complete.

Author
Tyson Jones

Referenced by TEST_CASE().