Unit test utilities

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. 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 (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 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, 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...
 
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...
 
QMatrix getConjugateTranspose (QMatrix a)
 Returns the conjugate transpose of the complex square matrix a. 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 getNormalised (QVector vec)
 Returns an L2-normalised copy of vec, using Kahan summation for improved accuracy. 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...
 
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...
 
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 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 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...
 
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 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...
 
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...
 

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

Definition at line 49 of file utilities.hpp.

◆ 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

Definition at line 60 of file utilities.hpp.

Function Documentation

◆ applyReferenceMatrix() [1/2]

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

Definition at line 692 of file utilities.cpp.

694  {
695  // for density matrices, op is left-multiplied only
696  int numQubits = calcLog2(state.size());
697  QMatrix leftOp = getFullOperatorMatrix(ctrls, numCtrls, targs, numTargs, op, numQubits);
698  state = leftOp * state;
699 }

References calcLog2(), and getFullOperatorMatrix().

◆ applyReferenceMatrix() [2/2]

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

Definition at line 686 of file utilities.cpp.

688  {
689  // for state-vectors, the op is always just left-multiplied
690  applyReferenceOp(state, ctrls, numCtrls, targs, numTargs, op);
691 }

References applyReferenceOp().

Referenced by TEST_CASE().

◆ 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

Definition at line 629 of file utilities.cpp.

631  {
632  int numQubits = calcLog2(state.size());
633  QMatrix leftOp = getFullOperatorMatrix(ctrls, numCtrls, targs, numTargs, op, numQubits);
634  QMatrix rightOp = getConjugateTranspose(leftOp);
635  state = leftOp * state * rightOp;
636 }

References calcLog2(), getConjugateTranspose(), and getFullOperatorMatrix().

◆ 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

Definition at line 637 of file utilities.cpp.

639  {
640  int targs[2] = {targ1, targ2};
641  applyReferenceOp(state, ctrls, numCtrls, targs, 2, op);
642 }

References applyReferenceOp().

◆ 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

Definition at line 643 of file utilities.cpp.

645  {
646  int targs[1] = {target};
647  applyReferenceOp(state, ctrls, numCtrls, targs, 1, op);
648 }

References applyReferenceOp().

◆ 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

Definition at line 649 of file utilities.cpp.

651  {
652  applyReferenceOp(state, NULL, 0, targs, numTargs, op);
653 }

References applyReferenceOp().

◆ 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

Definition at line 661 of file utilities.cpp.

663  {
664  int ctrls[1] = {ctrl};
665  applyReferenceOp(state, ctrls, 1, targs, numTargs, op);
666 }

References applyReferenceOp().

◆ 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

Definition at line 654 of file utilities.cpp.

656  {
657  int ctrls[1] = {ctrl};
658  int targs[1] = {targ};
659  applyReferenceOp(state, ctrls, 1, targs, 1, op);
660 }

References applyReferenceOp().

◆ 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

Definition at line 667 of file utilities.cpp.

669  {
670  int ctrls[1] = {ctrl};
671  int targs[2] = {targ1, targ2};
672  applyReferenceOp(state, ctrls, 1, targs, 2, op);
673 }

References applyReferenceOp().

◆ 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

Definition at line 674 of file utilities.cpp.

676  {
677  int targs[1] = {targ};
678  applyReferenceOp(state, NULL, 0, targs, 1, op);
679 }

References applyReferenceOp().

◆ 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

Definition at line 573 of file utilities.cpp.

575  {
576  int numQubits = calcLog2(state.size());
577  QMatrix fullOp = getFullOperatorMatrix(ctrls, numCtrls, targs, numTargs, op, numQubits);
578  state = fullOp * state;
579 }

References calcLog2(), and getFullOperatorMatrix().

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

Definition at line 580 of file utilities.cpp.

582  {
583  int targs[2] = {targ1, targ2};
584  applyReferenceOp(state, ctrls, numCtrls, targs, 2, op);
585 }

References applyReferenceOp().

◆ 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

Definition at line 586 of file utilities.cpp.

588  {
589  int targs[1] = {target};
590  applyReferenceOp(state, ctrls, numCtrls, targs, 1, op);
591 }

References applyReferenceOp().

◆ 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

Definition at line 592 of file utilities.cpp.

594  {
595  applyReferenceOp(state, NULL, 0, targs, numTargs, op);
596 }

References applyReferenceOp().

◆ 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

Definition at line 604 of file utilities.cpp.

606  {
607  int ctrls[1] = {ctrl};
608  applyReferenceOp(state, ctrls, 1, targs, numTargs, op);
609 }

References applyReferenceOp().

◆ 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

Definition at line 597 of file utilities.cpp.

599  {
600  int ctrls[1] = {ctrl};
601  int targs[1] = {targ};
602  applyReferenceOp(state, ctrls, 1, targs, 1, op);
603 }

References applyReferenceOp().

◆ 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

Definition at line 610 of file utilities.cpp.

612  {
613  int ctrls[1] = {ctrl};
614  int targs[2] = {targ1, targ2};
615  applyReferenceOp(state, ctrls, 1, targs, 2, op);
616 }

References applyReferenceOp().

◆ 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

Definition at line 617 of file utilities.cpp.

619  {
620  int targs[1] = {targ};
621  applyReferenceOp(state, NULL, 0, targs, 1, op);
622 }

References applyReferenceOp().

◆ areEqual() [1/9]

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

Definition at line 396 of file utilities.cpp.

396  {
397  DEMAND( a.size() == b.size() );
398 
399  for (size_t i=0; i<a.size(); i++)
400  for (size_t j=0; j<b.size(); j++)
401  if (abs(a[i][j] - b[i][j]) > REAL_EPS)
402  return false;
403  return true;
404 }

References DEMAND.

◆ areEqual() [2/9]

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

Definition at line 819 of file utilities.cpp.

819  {
820  return areEqual(qureg, matr, REAL_EPS);
821 }

References areEqual().

◆ areEqual() [3/9]

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

Definition at line 767 of file utilities.cpp.

767  {
768  DEMAND( qureg.isDensityMatrix );
769  DEMAND( (int) (matr.size()*matr.size()) == qureg.numAmpsTotal );
770 
771  // ensure local qureg.stateVec is up to date
772  copyStateFromGPU(qureg);
774 
775  // the starting index in vec of this node's qureg partition.
776  long long int startInd = qureg.chunkId * qureg.numAmpsPerChunk;
777  long long int globalInd, row, col, i;
778  int ampsAgree;
779 
780  // compare each of this node's amplitude to the corresponding matr sub-matrix
781  for (i=0; i<qureg.numAmpsPerChunk; i++) {
782  globalInd = startInd + i;
783  row = globalInd % matr.size();
784  col = globalInd / matr.size();
785  qreal realDif = absReal(qureg.stateVec.real[i] - real(matr[row][col]));
786  qreal imagDif = absReal(qureg.stateVec.imag[i] - imag(matr[row][col]));
787  ampsAgree = (realDif < precision && imagDif < precision);
788 
789  // DEBUG
790  if (!ampsAgree) {
791  printf("[msg from utilities.cpp] node %d has a disagreement at (global) index %lld of (%g) + i(%g)\n",
792  qureg.chunkId, globalInd, realDif, imagDif
793  );
794  }
795 
796  // break loop as soon as amplitudes disagree
797  if (!ampsAgree)
798  break;
799 
800  /* TODO:
801  * of the nodes which disagree, the lowest-rank should send its
802  * disagreeing (i, row, col, stateVec[i]) to rank 0 which should
803  * report it immediately (before the impending DEMAND failure)
804  * using FAIL_CHECK, so users can determine nature of disagreement
805  * (e.g. numerical precision).
806  * Note FAIL_CHECK accepts << like cout, e.g.
807  * FAIL_CHECK( "Amp at (" << row << ", " << col ") disagreed" );
808  */
809  }
810 
811  // if one node's partition wasn't equal, all-nodes must report not-equal
812  int allAmpsAgree = ampsAgree;
813 #ifdef DISTRIBUTED_MODE
814  MPI_Allreduce(&ampsAgree, &allAmpsAgree, 1, MPI_INT, MPI_LAND, MPI_COMM_WORLD);
815 #endif
816 
817  return allAmpsAgree;
818 }

References Qureg::chunkId, copyStateFromGPU(), DEMAND, Qureg::isDensityMatrix, Qureg::numAmpsPerChunk, Qureg::numAmpsTotal, qreal, QUEST_ENV, Qureg::stateVec, and syncQuESTEnv().

◆ areEqual() [4/9]

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

Definition at line 763 of file utilities.cpp.

763  {
764  return areEqual(qureg, vec, REAL_EPS);
765 }

References areEqual().

◆ areEqual() [5/9]

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

Definition at line 728 of file utilities.cpp.

728  {
729  DEMAND( !qureg.isDensityMatrix );
730  DEMAND( (int) vec.size() == qureg.numAmpsTotal );
731 
732  copyStateFromGPU(qureg);
734 
735  // the starting index in vec of this node's qureg partition.
736  long long int startInd = qureg.chunkId * qureg.numAmpsPerChunk;
737 
738  int ampsAgree = 1;
739  for (long long int i=0; i<qureg.numAmpsPerChunk; i++) {
740  qreal realDif = absReal(qureg.stateVec.real[i] - real(vec[startInd+i]));
741  qreal imagDif = absReal(qureg.stateVec.imag[i] - imag(vec[startInd+i]));
742 
743  if (realDif > precision || imagDif > precision) {
744  ampsAgree = 0;
745 
746  // debug
747  printf("Disagreement at %lld: %g + i(%g) VS %g + i(%g)\n",
748  startInd+i, qureg.stateVec.real[i], qureg.stateVec.imag[i],
749  real(vec[startInd+i]), imag(vec[startInd+i]));
750 
751  break;
752  }
753  }
754 
755  // if one node's partition wasn't equal, all-nodes must report not-equal
756  int allAmpsAgree = ampsAgree;
757 #ifdef DISTRIBUTED_MODE
758  MPI_Allreduce(&ampsAgree, &allAmpsAgree, 1, MPI_INT, MPI_LAND, MPI_COMM_WORLD);
759 #endif
760 
761  return allAmpsAgree;
762 }

References Qureg::chunkId, copyStateFromGPU(), DEMAND, Qureg::isDensityMatrix, Qureg::numAmpsPerChunk, Qureg::numAmpsTotal, qreal, QUEST_ENV, Qureg::stateVec, and syncQuESTEnv().

◆ areEqual() [6/9]

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

Definition at line 724 of file utilities.cpp.

724  {
725  return areEqual(qureg1, qureg2, REAL_EPS);
726 }

References areEqual().

◆ areEqual() [7/9]

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

Definition at line 701 of file utilities.cpp.

701  {
702  DEMAND( qureg1.isDensityMatrix == qureg2.isDensityMatrix );
703  DEMAND( qureg1.numAmpsTotal == qureg2.numAmpsTotal );
704 
705  copyStateFromGPU(qureg1);
706  copyStateFromGPU(qureg2);
708 
709  // loop terminates when areEqual = 0
710  int ampsAgree = 1;
711  for (long long int i=0; ampsAgree && i<qureg1.numAmpsPerChunk; i++)
712  ampsAgree = (
713  absReal(qureg1.stateVec.real[i] - qureg2.stateVec.real[i]) < precision
714  && absReal(qureg1.stateVec.imag[i] - qureg2.stateVec.imag[i]) < precision);
715 
716  // if one node's partition wasn't equal, all-nodes must report not-equal
717  int allAmpsAgree = ampsAgree;
718 #ifdef DISTRIBUTED_MODE
719  MPI_Allreduce(&ampsAgree, &allAmpsAgree, 1, MPI_INT, MPI_LAND, MPI_COMM_WORLD);
720 #endif
721 
722  return allAmpsAgree;
723 }

References copyStateFromGPU(), DEMAND, Qureg::isDensityMatrix, Qureg::numAmpsPerChunk, Qureg::numAmpsTotal, QUEST_ENV, Qureg::stateVec, and syncQuESTEnv().

◆ areEqual() [8/9]

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

Definition at line 387 of file utilities.cpp.

387  {
388  DEMAND( a.size() == b.size() );
389 
390  for (size_t i=0; i<a.size(); i++)
391  if (abs(a[i] - b[i]) > REAL_EPS)
392  return false;
393  return true;
394 }

References DEMAND.

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

◆ areEqual() [9/9]

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

Definition at line 823 of file utilities.cpp.

823  {
824 
825  qreal dif;
826  for (size_t i=0; i<vec.size(); i++) {
827  dif = abs(real(vec[i]) - reals[i]);
828  if (dif > REAL_EPS)
829  return false;
830  dif = abs(imag(vec[i]) - imags[i]);
831  if (dif > REAL_EPS)
832  return false;
833  }
834  return true;
835 }

References qreal.

◆ 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

Definition at line 1268 of file utilities.cpp.

1268  {
1269  return Catch::Generators::GeneratorWrapper<int*>(
1270  std::unique_ptr<Catch::Generators::IGenerator<int*>>(
1271  new SequenceGenerator<int>(1, numBits)));
1272 }

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.

Definition at line 292 of file QuEST_validation.c.

292  {
293  unsigned int l = 0;
294  while (num >>= 1)
295  l++;
296  return l;
297 }

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

◆ getConjugateTranspose()

QMatrix getConjugateTranspose ( QMatrix  a)

Returns the conjugate transpose of the complex square matrix a.

Author
Tyson Jones

Definition at line 179 of file utilities.cpp.

179  {
180  QMatrix b = a;
181  for (size_t r=0; r<a.size(); r++)
182  for (size_t c=0; c<a.size(); c++)
183  b[r][c] = conj(a[c][r]);
184  return b;
185 }

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

◆ 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

Definition at line 187 of file utilities.cpp.

187  {
188 
189  // ensure diagonal
190  for (size_t r=0; r<a.size(); r++)
191  for (size_t c=0; c<a.size(); c++) {
192  if (r == c)
193  continue;
194  DEMAND( a[r][c] == 0. );
195  }
196 
197  // exp(diagonal) = diagonal(exp)
198  QMatrix diag = a;
199  for (size_t i=0; i<a.size(); i++)
200  diag[i][i] = exp(diag[i][i]);
201 
202  return diag;
203 }

References DEMAND.

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

Definition at line 205 of file utilities.cpp.

205  {
206  QMatrix iden = getIdentityMatrix(a.size());
207  QMatrix expo = (cos(angle/2) * iden) + (-1i * sin(angle/2) * a);
208  return expo;
209 }

References getIdentityMatrix().

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

Definition at line 293 of file utilities.cpp.

295  {
296  DEMAND( numCtrls >= 0 );
297  DEMAND( numTargs >= 0 );
298  DEMAND( numQubits >= (numCtrls+numTargs) );
299  DEMAND( op.size() == (1u << numTargs) );
300 
301  // copy {ctrls} and {targs}to restore at end
302  std::vector<int> ctrlsCopy(ctrls, ctrls+numCtrls);
303  std::vector<int> targsCopy(targs, targs+numTargs);
304 
305  // full-state matrix of qubit swaps
306  QMatrix swaps = getIdentityMatrix(1 << numQubits);
307  QMatrix unswaps = getIdentityMatrix(1 << numQubits);
308  QMatrix matr;
309 
310  // swap targs to {0, ..., numTargs-1}
311  for (int i=0; i<numTargs; i++) {
312  if (i != targs[i]) {
313  matr = getSwapMatrix(i, targs[i], numQubits);
314  swaps = matr * swaps;
315  unswaps = unswaps * matr;
316 
317  // even if this is the last targ, ctrls might still need updating
319  i, targs[i], (i < numTargs-1)? &targs[i+1] : NULL,
320  numTargs-i-1, ctrls, numCtrls);
321  }
322  }
323 
324  // swap ctrls to {numTargs, ..., numTargs+numCtrls-1}
325  for (int i=0; i<numCtrls; i++) {
326  int newInd = numTargs+i;
327  if (newInd != ctrls[i]) {
328  matr = getSwapMatrix(newInd, ctrls[i], numQubits);
329  swaps = matr * swaps;
330  unswaps = unswaps * matr;
331 
332  // update remaining ctrls (if any exist)
333  if (i < numCtrls-1)
334  updateIndices(newInd, ctrls[i], NULL, 0, &ctrls[i+1], numCtrls-i-1);
335  }
336  }
337 
338  // construct controlled-op matrix for qubits {0, ..., numCtrls+numTargs-1}
339  size_t dim = 1 << (numCtrls+numTargs);
340  QMatrix fullOp = getIdentityMatrix(dim);
341  setSubMatrix(fullOp, op, dim-op.size(), dim-op.size());
342 
343  // create full-state controlled-op matrix (left-pad identities)
344  if (numQubits > numCtrls+numTargs) {
345  size_t pad = 1 << (numQubits - numCtrls - numTargs);
346  fullOp = getKroneckerProduct(getIdentityMatrix(pad), fullOp);
347  }
348 
349  // apply swap to either side (to swap qubits back and forth)
350  fullOp = unswaps * fullOp * swaps;
351 
352  // restore {ctrls and targs}
353  for (int i=0; i<numCtrls; i++)
354  ctrls[i] = ctrlsCopy[i];
355  for (int i=0; i<numTargs; i++)
356  targs[i] = targsCopy[i];
357 
358  return fullOp;
359 }

References DEMAND, getIdentityMatrix(), getKroneckerProduct(), getSwapMatrix(), setSubMatrix(), and updateIndices().

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

◆ getIdentityMatrix()

QMatrix getIdentityMatrix ( size_t  dim)

Returns a dim-by-dim identity matrix.

Author
Tyson Jones

Definition at line 151 of file utilities.cpp.

151  {
152  DEMAND( dim > 1 );
153  QMatrix matr = getZeroMatrix(dim);
154  for (size_t i=0; i<dim; i++)
155  matr[i][i] = 1;
156  return matr;
157 }

References DEMAND, and getZeroMatrix().

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

Definition at line 159 of file utilities.cpp.

159  {
160  DEMAND( ket.size() == bra.size() );
161  QMatrix mat = getZeroMatrix(ket.size());
162 
163  for (size_t r=0; r<ket.size(); r++)
164  for (size_t c=0; c<ket.size(); c++)
165  mat[r][c] = ket[r] * conj(bra[c]);
166  return mat;
167 }

References DEMAND, and getZeroMatrix().

Referenced by getRandomDensityMatrix(), and TEST_CASE().

◆ getKroneckerProduct()

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

Definition at line 169 of file utilities.cpp.

169  {
170  QMatrix prod = getZeroMatrix(a.size() * b.size());
171  for (size_t r=0; r<b.size(); r++)
172  for (size_t c=0; c<b.size(); c++)
173  for (size_t i=0; i<a.size(); i++)
174  for (size_t j=0; j<a.size(); j++)
175  prod[r+b.size()*i][c+b.size()*j] = a[i][j] * b[r][c];
176  return prod;
177 }

References getZeroMatrix().

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

◆ getNormalised()

QVector getNormalised ( QVector  vec)

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

Author
Tyson Jones

Definition at line 431 of file utilities.cpp.

431  {
432  qreal norm = 0;
433  qreal y, t, c;
434  c = 0;
435 
436  for (size_t i=0; i<vec.size(); i++) {
437  y = real(vec[i])*real(vec[i]) - c;
438  t = norm + y;
439  c = ( t - norm ) - y;
440  norm = t;
441 
442  y = imag(vec[i])*imag(vec[i]) - c;
443  t = norm + y;
444  c = ( t - norm ) - y;
445  norm = t;
446  }
447 
448  for (size_t i=0; i<vec.size(); i++)
449  vec[i] /= sqrt(norm);
450  return vec;
451 }

References qreal.

Referenced by getRandomStateVector().

◆ 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

Definition at line 457 of file utilities.cpp.

457  {
458  DEMAND( numQb > 0 );
459 
460  // generate random probabilities to weight random pure states
461  int dim = 1<<numQb;
462  qreal probs[dim];
463  qreal probNorm = 0;
464  for (int i=0; i<dim; i++) {
465  probs[i] = getRandomReal(0, 1);
466  probNorm += probs[i];
467  }
468  for (int i=0; i<dim; i++)
469  probs[i] /= probNorm;
470 
471  // add random pure states
472  QMatrix dens = getZeroMatrix(dim);
473  for (int i=0; i<dim; i++) {
474  QVector pure = getRandomStateVector(numQb);
475  dens += probs[i] * getKetBra(pure, pure);
476  }
477 
478  return dens;
479 }

References DEMAND, getKetBra(), getRandomReal(), getRandomStateVector(), getZeroMatrix(), and qreal.

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

Definition at line 481 of file utilities.cpp.

481  {
482  return round(getRandomReal(min, max-1));
483 }

References getRandomReal().

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

Definition at line 533 of file utilities.cpp.

533  {
534  DEMAND( numOps >= 1 );
535  DEMAND( numOps <= 4*numQb*numQb );
536 
537  // generate random unitaries
538  std::vector<QMatrix> ops;
539  for (int i=0; i<numOps; i++)
540  ops.push_back(getRandomUnitary(numQb));
541 
542  // generate random weights
543  qreal weights[numOps];
544  for (int i=0; i<numOps; i++)
545  weights[i] = getRandomReal(0, 1);
546 
547  // normalise random weights
548  qreal weightSum = 0;
549  for (int i=0; i<numOps; i++)
550  weightSum += weights[i];
551  for (int i=0; i<numOps; i++)
552  weights[i] = sqrt(weights[i]/weightSum);
553 
554  // normalise ops
555  for (int i=0; i<numOps; i++)
556  ops[i] *= weights[i];
557 
558  // check what we produced was a valid Kraus map
559  QMatrix iden = getIdentityMatrix(1 << numQb);
560  QMatrix prodSum = getZeroMatrix(1 << numQb);
561  for (int i=0; i<numOps; i++)
562  prodSum += getConjugateTranspose(ops[i]) * ops[i];
563  DEMAND( areEqual(prodSum, iden) );
564 
565  return ops;
566 }

References areEqual(), DEMAND, getConjugateTranspose(), getIdentityMatrix(), getRandomReal(), getRandomUnitary(), getZeroMatrix(), and qreal.

Referenced by TEST_CASE().

◆ 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

Definition at line 368 of file utilities.cpp.

368  {
369  DEMAND( dim > 1 );
370 
371  QMatrix matr = getZeroMatrix(dim);
372  for (int i=0; i<dim; i++) {
373  for (int j=0; j<dim; j++) {
374 
375  // generate 2 normally-distributed random numbers via Box-Muller
376  qreal a = rand()/(qreal) RAND_MAX;
377  qreal b = rand()/(qreal) RAND_MAX;
378  qreal r1 = sqrt(-2 * log(a)) * cos(2 * 3.14159265 * b);
379  qreal r2 = sqrt(-2 * log(a)) * sin(2 * 3.14159265 * b);
380 
381  matr[i][j] = r1 + r2*1i;
382  }
383  }
384  return matr;
385 }

References DEMAND, getZeroMatrix(), and qreal.

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

Definition at line 420 of file utilities.cpp.

420  {
421  QVector vec = QVector(dim);
422  for (int i=0; i<dim; i++)
423  vec[i] = getRandomReal(-1,1) + 1i*getRandomReal(-1,1);
424 
425  // check we didn't get the impossibly-unlikely zero-amplitude outcome
426  DEMAND( real(vec[0]) != 0 );
427 
428  return vec;
429 }

References DEMAND, and getRandomReal().

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

Definition at line 410 of file utilities.cpp.

410  {
411  DEMAND( min <= max );
412  qreal r = min + (max - min) * (rand() / (qreal) RAND_MAX);
413 
414  // check bounds satisfied
415  DEMAND( r >= min );
416  DEMAND( r <= max );
417  return r;
418 }

References DEMAND, and qreal.

Referenced by getRandomDensityMatrix(), getRandomInt(), getRandomKrausMap(), getRandomQVector(), 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

Definition at line 453 of file utilities.cpp.

453  {
454  return getNormalised(getRandomQVector(1<<numQb));
455 }

References getNormalised(), and getRandomQVector().

Referenced by getRandomDensityMatrix(), 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

Definition at line 485 of file utilities.cpp.

485  {
486  DEMAND( numQb >= 1 );
487 
488  QMatrix matr = getRandomQMatrix(1 << numQb);
489 
490  for (size_t i=0; i<matr.size(); i++) {
491  QVector row = matr[i];
492 
493  // compute new orthogonal row by subtracting proj row onto prevs
494  for (int k=i-1; k>=0; k--) {
495 
496  // compute row . prev = sum_n row_n conj(prev_n)
497  qcomp prod = 0;
498  for (size_t n=0; n<row.size(); n++)
499  prod += row[n] * conj(matr[k][n]);
500 
501  // subtract (proj row onto prev) = (prod * prev) from final row
502  for (size_t n=0; n<row.size(); n++)
503  matr[i][n] -= prod * matr[k][n];
504  }
505 
506  // compute row magnitude
507  qreal mag = 0;
508  for (size_t j=0; j<row.size(); j++)
509  mag += pow(abs(matr[i][j]), 2);
510  mag = sqrt(mag);
511 
512  // normalise row
513  for (size_t j=0; j<row.size(); j++)
514  matr[i][j] /= mag;
515  }
516 
517  // ensure matrix is indeed unitary
518  QMatrix conjprod = matr * getConjugateTranspose(matr);
519  QMatrix iden = getIdentityMatrix(1 << numQb);
520 
521  // generating big unitary matrices is hard; if we fail, default to identity
522  if ( numQb >= 3 && !areEqual(conjprod, iden) ) {
523 
524  matr = getIdentityMatrix(1 << numQb);
525  conjprod = matr;
526  }
527  DEMAND( areEqual(conjprod, iden) );
528 
529  // return the new orthonormal matrix
530  return matr;
531 }

References areEqual(), DEMAND, getConjugateTranspose(), getIdentityMatrix(), getRandomQMatrix(), qcomp, and qreal.

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

Definition at line 219 of file utilities.cpp.

219  {
220  DEMAND( numQb > 1 );
221  DEMAND( (qb1 >= 0 && qb1 < numQb) );
222  DEMAND( (qb2 >= 0 && qb2 < numQb) );
223 
224  if (qb1 > qb2)
225  std::swap(qb1, qb2);
226 
227  if (qb1 == qb2)
228  return getIdentityMatrix(1 << numQb);
229 
230  QMatrix swap;
231 
232  if (qb2 == qb1 + 1) {
233  // qubits are adjacent
234  swap = QMatrix{{1,0,0,0},{0,0,1,0},{0,1,0,0},{0,0,0,1}};
235 
236  } else {
237  // qubits are distant
238  int block = 1 << (qb2 - qb1);
239  swap = getZeroMatrix(block*2);
240  QMatrix iden = getIdentityMatrix(block/2);
241 
242  // Lemma 3.1 of arxiv.org/pdf/1711.09765.pdf
243  QMatrix p0{{1,0},{0,0}};
244  QMatrix l0{{0,1},{0,0}};
245  QMatrix l1{{0,0},{1,0}};
246  QMatrix p1{{0,0},{0,1}};
247 
248  /* notating a^(n+1) = identity(1<<n) (otimes) a, we construct the matrix
249  * [ p0^(N) l1^N ]
250  * [ l0^(N) p1^N ]
251  * where N = qb2 - qb1 */
252  setSubMatrix(swap, getKroneckerProduct(iden, p0), 0, 0);
253  setSubMatrix(swap, getKroneckerProduct(iden, l0), block, 0);
254  setSubMatrix(swap, getKroneckerProduct(iden, l1), 0, block);
255  setSubMatrix(swap, getKroneckerProduct(iden, p1), block, block);
256  }
257 
258  // pad swap with outer identities
259  if (qb1 > 0)
260  swap = getKroneckerProduct(swap, getIdentityMatrix(1<<qb1));
261  if (qb2 < numQb-1)
262  swap = getKroneckerProduct(getIdentityMatrix(1<<(numQb-qb2-1)), swap);
263 
264  return swap;
265 }

References DEMAND, getIdentityMatrix(), getKroneckerProduct(), getZeroMatrix(), and setSubMatrix().

Referenced by getFullOperatorMatrix().

◆ getZeroMatrix()

QMatrix getZeroMatrix ( size_t  dim)

Returns a dim-by-dim square complex matrix, initialised to all zeroes.

Author
Tyson Jones

Definition at line 143 of file utilities.cpp.

143  {
144  DEMAND( dim > 1 );
145  QMatrix matr = QMatrix(dim);
146  for (size_t i=0; i<dim; i++)
147  matr[i].resize(dim);
148  return matr;
149 }

References DEMAND.

Referenced by getIdentityMatrix(), getKetBra(), getKroneckerProduct(), getRandomDensityMatrix(), getRandomKrausMap(), getRandomQMatrix(), getSwapMatrix(), TEST_CASE(), and toQMatrix().

◆ 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

Definition at line 1278 of file utilities.cpp.

1278  {
1279  return Catch::Generators::GeneratorWrapper<pauliOpType*>(
1280  std::unique_ptr<Catch::Generators::IGenerator<pauliOpType*>>(
1281  new SequenceGenerator<pauliOpType>(PAULI_Z, numPaulis)));
1282 }

References PAULI_Z.

◆ 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

Definition at line 1273 of file utilities.cpp.

1273  {
1274  return Catch::Generators::GeneratorWrapper<int*>(
1275  std::unique_ptr<Catch::Generators::IGenerator<int*>>(
1276  new SequenceGenerator<int>(base-1, numDigits)));
1277 }

◆ setRandomPauliSum() [1/2]

void setRandomPauliSum ( PauliHamil  hamil)

Populates hamil with random coefficients and pauli codes.

Author
Tyson Jones

Definition at line 1062 of file utilities.cpp.

1062  {
1063  setRandomPauliSum(hamil.termCoeffs, hamil.pauliCodes, hamil.numQubits, hamil.numSumTerms);
1064 }

References PauliHamil::numQubits, PauliHamil::numSumTerms, PauliHamil::pauliCodes, setRandomPauliSum(), and PauliHamil::termCoeffs.

◆ 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

Definition at line 1054 of file utilities.cpp.

1054  {
1055  int i=0;
1056  for (int n=0; n<numTerms; n++) {
1057  coeffs[n] = getRandomReal(-5, 5);
1058  for (int q=0; q<numQubits; q++)
1059  codes[i++] = (pauliOpType) getRandomInt(0,4);
1060  }
1061 }

References getRandomInt(), and getRandomReal().

Referenced by setRandomPauliSum(), and 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

Definition at line 211 of file utilities.cpp.

211  {
212  DEMAND( sub.size() + r <= dest.size() );
213  DEMAND( sub.size() + c <= dest.size() );
214  for (size_t i=0; i<sub.size(); i++)
215  for (size_t j=0; j<sub.size(); j++)
216  dest[r+i][c+j] = sub[i][j];
217 }

References DEMAND.

Referenced by getFullOperatorMatrix(), and getSwapMatrix().

◆ 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

Definition at line 1199 of file utilities.cpp.

1201  {
1202  return Catch::Generators::GeneratorWrapper<int*>(
1203  std::unique_ptr<Catch::Generators::IGenerator<int*>>(
1204  new SubListGenerator(list, len, sublen)));
1205 }

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

Definition at line 846 of file utilities.cpp.

846  {
847  DEMAND( qm.size() == 2 );
848  ComplexMatrix2 cm;
849  macro_copyQMatrix(cm, qm);
850  return cm;
851 }

References DEMAND, and macro_copyQMatrix.

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

Definition at line 852 of file utilities.cpp.

852  {
853  DEMAND( qm.size() == 4 );
854  ComplexMatrix4 cm;
855  macro_copyQMatrix(cm, qm);
856  return cm;
857 }

References DEMAND, and macro_copyQMatrix.

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

Definition at line 858 of file utilities.cpp.

858  {
859  DEMAND( qm.size() == (1u<<cm.numQubits) );
860  macro_copyQMatrix(cm, qm);
861 }

References DEMAND, macro_copyQMatrix, and ComplexMatrixN::numQubits.

Referenced by TEST_CASE().

◆ toQMatrix() [1/8]

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

Definition at line 887 of file utilities.cpp.

887  {
888  qcomp a = qcomp(alpha.real, alpha.imag);
889  qcomp b = qcomp(beta.real, beta.imag);
890  QMatrix matr{
891  {a, -conj(b)},
892  {b, conj(a)}};
893  return matr;
894 }

References Complex::imag, qcomp, and Complex::real.

◆ toQMatrix() [2/8]

QMatrix toQMatrix ( ComplexMatrix2  src)

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

Author
Tyson Jones

Definition at line 869 of file utilities.cpp.

869  {
870  QMatrix dest = getZeroMatrix(2);
871  macro_copyComplexMatrix(dest, src);
872  return dest;
873 }

References getZeroMatrix(), and macro_copyComplexMatrix.

Referenced by TEST_CASE(), and toQMatrix().

◆ toQMatrix() [3/8]

QMatrix toQMatrix ( ComplexMatrix4  src)

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

Author
Tyson Jones

Definition at line 874 of file utilities.cpp.

874  {
875  QMatrix dest = getZeroMatrix(4);
876  macro_copyComplexMatrix(dest, src);
877  return dest;
878 }

References getZeroMatrix(), and macro_copyComplexMatrix.

◆ toQMatrix() [4/8]

QMatrix toQMatrix ( ComplexMatrixN  src)

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

Author
Tyson Jones

Definition at line 879 of file utilities.cpp.

879  {
880  DEMAND( src.real != NULL );
881  DEMAND( src.imag != NULL );
882  QMatrix dest = getZeroMatrix(1 << src.numQubits);
883  macro_copyComplexMatrix(dest, src);
884  return dest;
885 }

References DEMAND, getZeroMatrix(), ComplexMatrixN::imag, macro_copyComplexMatrix, ComplexMatrixN::numQubits, and ComplexMatrixN::real.

◆ toQMatrix() [5/8]

QMatrix toQMatrix ( DiagonalOp  op)

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

Author
Tyson Jones

Definition at line 1018 of file utilities.cpp.

1018  {
1019  QVector vec = toQVector(op);
1020  QMatrix mat = getZeroMatrix(1LL << op.numQubits);
1021  for (size_t i=0; i<mat.size(); i++)
1022  mat[i][i] = vec[i];
1023  return mat;
1024 }

References getZeroMatrix(), DiagonalOp::numQubits, and toQVector().

◆ toQMatrix() [6/8]

QMatrix toQMatrix ( PauliHamil  hamil)

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

Author
Tyson Jones

Definition at line 1095 of file utilities.cpp.

1095  {
1096  return toQMatrix(hamil.termCoeffs, hamil.pauliCodes, hamil.numQubits, hamil.numSumTerms);
1097 }

References PauliHamil::numQubits, PauliHamil::numSumTerms, PauliHamil::pauliCodes, PauliHamil::termCoeffs, and toQMatrix().

◆ toQMatrix() [7/8]

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

Definition at line 1066 of file utilities.cpp.

1066  {
1067 
1068  // produce a numTargs-big matrix 'pauliSum' by pauli-matrix tensoring and summing
1069  QMatrix iMatr{{1,0},{0,1}};
1070  QMatrix xMatr{{0,1},{1,0}};
1071  QMatrix yMatr{{0,-1i},{1i,0}};
1072  QMatrix zMatr{{1,0},{0,-1}};
1073  QMatrix pauliSum = getZeroMatrix(1<<NUM_QUBITS);
1074 
1075  for (int t=0; t<numTerms; t++) {
1076  QMatrix pauliProd = QMatrix{{1}};
1077 
1078  for (int q=0; q<numQubits; q++) {
1079  int i = q + t*numQubits;
1080 
1081  QMatrix fac;
1082  pauliOpType code = paulis[i];
1083  if (code == PAULI_I) fac = iMatr;
1084  if (code == PAULI_X) fac = xMatr;
1085  if (code == PAULI_Y) fac = yMatr;
1086  if (code == PAULI_Z) fac = zMatr;
1087  pauliProd = getKroneckerProduct(fac, pauliProd);
1088  }
1089  pauliSum += coeffs[t] * pauliProd;
1090  }
1091 
1092  // a now 2^numQubits by 2^numQubits Hermitian matrix
1093  return pauliSum;
1094 }

References getKroneckerProduct(), getZeroMatrix(), NUM_QUBITS, PAULI_I, PAULI_X, PAULI_Y, and PAULI_Z.

◆ toQMatrix() [8/8]

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

Definition at line 896 of file utilities.cpp.

896  {
897  DEMAND( qureg.isDensityMatrix );
898 #ifdef DISTRIBUTED_MODE
899  DEMAND( qureg.numAmpsTotal < MPI_MAX_AMPS_IN_MSG );
900 #endif
901 
902  // ensure local qureg.stateVec is up to date
903  copyStateFromGPU(qureg);
905 
906  qreal* fullRe;
907  qreal* fullIm;
908 
909  // in distributed mode, give every node the full state vector
910 #ifdef DISTRIBUTED_MODE
911  fullRe = (qreal*) malloc(qureg.numAmpsTotal * sizeof *fullRe);
912  fullIm = (qreal*) malloc(qureg.numAmpsTotal * sizeof *fullIm);
913  MPI_Allgather(
914  qureg.stateVec.real, qureg.numAmpsPerChunk, MPI_QuEST_REAL,
915  fullRe, qureg.numAmpsPerChunk, MPI_QuEST_REAL, MPI_COMM_WORLD);
916  MPI_Allgather(
917  qureg.stateVec.imag, qureg.numAmpsPerChunk, MPI_QuEST_REAL,
918  fullIm, qureg.numAmpsPerChunk, MPI_QuEST_REAL, MPI_COMM_WORLD);
919 #else
920  fullRe = qureg.stateVec.real;
921  fullIm = qureg.stateVec.imag;
922 #endif
923 
924  // copy full state vector into a QVector
925  long long int dim = (1 << qureg.numQubitsRepresented);
926  QMatrix matr = getZeroMatrix(dim);
927  for (long long int n=0; n<qureg.numAmpsTotal; n++)
928  matr[n%dim][n/dim] = qcomp(fullRe[n], fullIm[n]);
929 
930  // clean up if we malloc'd the distributed array
931 #ifdef DISTRIBUTED_MODE
932  free(fullRe);
933  free(fullIm);
934 #endif
935  return matr;
936 }

References copyStateFromGPU(), DEMAND, getZeroMatrix(), Qureg::isDensityMatrix, Qureg::numAmpsPerChunk, Qureg::numAmpsTotal, Qureg::numQubitsRepresented, qcomp, qreal, QUEST_ENV, Qureg::stateVec, and syncQuESTEnv().

◆ 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

Definition at line 1039 of file utilities.cpp.

1039  {
1040  DEMAND( qureg.isDensityMatrix );
1041  DEMAND( (1 << qureg.numQubitsRepresented) == (int) mat.size() );
1042 
1044 
1045  int len = (1 << qureg.numQubitsRepresented);
1046  for (int i=0; i<qureg.numAmpsPerChunk; i++) {
1047  int ind = qureg.chunkId*qureg.numAmpsPerChunk + i;
1048  qureg.stateVec.real[i] = real(mat[ind%len][ind/len]);
1049  qureg.stateVec.imag[i] = imag(mat[ind%len][ind/len]);
1050  }
1051  copyStateToGPU(qureg);
1052 }

References Qureg::chunkId, copyStateToGPU(), DEMAND, Qureg::isDensityMatrix, Qureg::numAmpsPerChunk, Qureg::numQubitsRepresented, QUEST_ENV, Qureg::stateVec, and syncQuESTEnv().

◆ 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

Definition at line 1026 of file utilities.cpp.

1026  {
1027  DEMAND( !qureg.isDensityMatrix );
1028  DEMAND( qureg.numAmpsTotal == (int) vec.size() );
1029 
1031 
1032  for (int i=0; i<qureg.numAmpsPerChunk; i++) {
1033  int ind = qureg.chunkId*qureg.numAmpsPerChunk + i;
1034  qureg.stateVec.real[i] = real(vec[ind]);
1035  qureg.stateVec.imag[i] = imag(vec[ind]);
1036  }
1037  copyStateToGPU(qureg);
1038 }

References Qureg::chunkId, copyStateToGPU(), DEMAND, Qureg::isDensityMatrix, Qureg::numAmpsPerChunk, Qureg::numAmpsTotal, QUEST_ENV, Qureg::stateVec, and syncQuESTEnv().

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

Definition at line 980 of file utilities.cpp.

980  {
981  long long int totalElems = (1LL << op.numQubits);
982 #ifdef DISTRIBUTED_MODE
983  DEMAND( totalElems < MPI_MAX_AMPS_IN_MSG );
984 #endif
985 
986  qreal* fullRe;
987  qreal* fullIm;
988 
989  // in distributed mode, give every node the full diagonal operator
990 #ifdef DISTRIBUTED_MODE
991  fullRe = (qreal*) malloc(totalElems * sizeof *fullRe);
992  fullIm = (qreal*) malloc(totalElems * sizeof *fullIm);
993 
994  MPI_Allgather(
995  op.real, op.numElemsPerChunk, MPI_QuEST_REAL,
996  fullRe, op.numElemsPerChunk, MPI_QuEST_REAL, MPI_COMM_WORLD);
997  MPI_Allgather(
998  op.imag, op.numElemsPerChunk, MPI_QuEST_REAL,
999  fullIm, op.numElemsPerChunk, MPI_QuEST_REAL, MPI_COMM_WORLD);
1000 #else
1001  fullRe = op.real;
1002  fullIm = op.imag;
1003 #endif
1004 
1005  // copy full state vector into a QVector
1006  QVector vec = QVector(totalElems);
1007  for (long long int i=0; i<totalElems; i++)
1008  vec[i] = qcomp(fullRe[i], fullIm[i]);
1009 
1010  // clean up if we malloc'd distrib array
1011 #ifdef DISTRIBUTED_MODE
1012  free(fullRe);
1013  free(fullIm);
1014 #endif
1015  return vec;
1016 }

References DEMAND, DiagonalOp::imag, DiagonalOp::numElemsPerChunk, DiagonalOp::numQubits, qcomp, qreal, and DiagonalOp::real.

◆ 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

Definition at line 938 of file utilities.cpp.

938  {
939  DEMAND( !qureg.isDensityMatrix );
940 #ifdef DISTRIBUTED_MODE
941  DEMAND( qureg.numAmpsTotal < MPI_MAX_AMPS_IN_MSG );
942 #endif
943 
944  // ensure local qureg.stateVec is up to date
945  copyStateFromGPU(qureg);
947 
948  qreal* fullRe;
949  qreal* fullIm;
950 
951  // in distributed mode, give every node the full state vector
952 #ifdef DISTRIBUTED_MODE
953  fullRe = (qreal*) malloc(qureg.numAmpsTotal * sizeof *fullRe);
954  fullIm = (qreal*) malloc(qureg.numAmpsTotal * sizeof *fullIm);
955 
956  MPI_Allgather(
957  qureg.stateVec.real, qureg.numAmpsPerChunk, MPI_QuEST_REAL,
958  fullRe, qureg.numAmpsPerChunk, MPI_QuEST_REAL, MPI_COMM_WORLD);
959  MPI_Allgather(
960  qureg.stateVec.imag, qureg.numAmpsPerChunk, MPI_QuEST_REAL,
961  fullIm, qureg.numAmpsPerChunk, MPI_QuEST_REAL, MPI_COMM_WORLD);
962 #else
963  fullRe = qureg.stateVec.real;
964  fullIm = qureg.stateVec.imag;
965 #endif
966 
967  // copy full state vector into a QVector
968  QVector vec = QVector(qureg.numAmpsTotal);
969  for (long long int i=0; i<qureg.numAmpsTotal; i++)
970  vec[i] = qcomp(fullRe[i], fullIm[i]);
971 
972  // clean up if we malloc'd distrib array
973 #ifdef DISTRIBUTED_MODE
974  free(fullRe);
975  free(fullIm);
976 #endif
977  return vec;
978 }

References copyStateFromGPU(), DEMAND, Qureg::isDensityMatrix, Qureg::numAmpsPerChunk, Qureg::numAmpsTotal, qcomp, qreal, QUEST_ENV, Qureg::stateVec, and syncQuESTEnv().

Referenced by TEST_CASE(), and toQMatrix().

pauliOpType
Codes for specifying Pauli operators.
Definition: QuEST.h:96
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...
Definition: utilities.cpp:293
void copyStateFromGPU(Qureg qureg)
In GPU mode, this copies the state-vector (or density matrix) from GPU memory (qureg....
Definition: QuEST_cpu.c:39
QuESTEnv QUEST_ENV
The global QuESTEnv instance, to be created and destroyed once in this main(), so that the MPI enviro...
Definition: main.cpp:20
void syncQuESTEnv(QuESTEnv env)
Guarantees that all code up to the given point has been executed on all nodes (if running in distribu...
@ PAULI_Z
Definition: QuEST.h:96
#define macro_copyComplexMatrix(dest, src)
Copies ComplexMatrix structures into a QMatrix.
Definition: utilities.cpp:864
@ PAULI_I
Definition: QuEST.h:96
QMatrix getConjugateTranspose(QMatrix a)
Returns the conjugate transpose of the complex square matrix a.
Definition: utilities.cpp:179
int getRandomInt(int min, int max)
Returns a random integer between min (inclusive) and max (exclusive), from the uniform distribution.
Definition: utilities.cpp:481
#define NUM_QUBITS
The default number of qubits in the registers created for unit testing (both statevectors and density...
Definition: utilities.hpp:36
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-ne...
Definition: utilities.cpp:219
QMatrix getRandomUnitary(int numQb)
Returns a uniformly random (under Haar) 2^numQb-by-2^numQb unitary matrix.
Definition: utilities.cpp:485
qreal getRandomReal(qreal min, qreal max)
Returns a random real between min (inclusive) and max (exclusive), from the uniform distribution.
Definition: utilities.cpp:410
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...
Definition: utilities.cpp:159
Represents a 4x4 matrix of complex numbers.
Definition: QuEST.h:125
#define qreal
QMatrix toQMatrix(ComplexMatrix2 src)
Returns a copy of the given 2-by-2 matrix.
Definition: utilities.cpp:869
unsigned int calcLog2(long unsigned int res)
Returns log2 of numbers which must be gauranteed to be 2^n.
Definition: utilities.cpp:361
@ PAULI_X
Definition: QuEST.h:96
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...
Definition: utilities.cpp:1054
int chunkId
The position of the chunk of the state vector held by this process in the full state vector.
Definition: QuEST.h:217
std::vector< qcomp > QVector
A complex vector, which can be zero-initialised with QVector(numAmps).
Definition: utilities.hpp:60
qreal * imag
The imaginary values of the 2^numQubits complex elements.
Definition: QuEST.h:191
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 + d...
Definition: utilities.cpp:211
QVector toQVector(Qureg qureg)
Returns an equal-size copy of the given state-vector qureg.
Definition: utilities.cpp:938
long long int numAmpsPerChunk
Number of probability amplitudes held in stateVec by this process In the non-MPI version,...
Definition: QuEST.h:213
qreal * termCoeffs
The coefficient of each Pauli product. This is a length numSumTerms array.
Definition: QuEST.h:164
#define qcomp
enum pauliOpType * pauliCodes
The Pauli operators acting on each qubit, flattened over every operator.
Definition: QuEST.h:162
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:36
int numQubits
The number of qubits this operator can act on (informing its size)
Definition: QuEST.h:181
int numSumTerms
The number of terms in the weighted sum, or the number of Pauli products.
Definition: QuEST.h:166
@ PAULI_Y
Definition: QuEST.h:96
QVector getRandomStateVector(int numQb)
Returns a random numQb-length L2-normalised state-vector from an undisclosed distribution.
Definition: utilities.cpp:453
qreal ** real
Definition: QuEST.h:139
QMatrix getRandomQMatrix(int dim)
Returns a dim-by-dim complex matrix, where the real and imaginary value of each element are independe...
Definition: utilities.cpp:368
qreal ** imag
Definition: QuEST.h:140
ComplexArray stateVec
Computational state amplitudes - a subset thereof in the MPI version.
Definition: QuEST.h:222
long long int numElemsPerChunk
The number of the 2^numQubits amplitudes stored on each distributed node.
Definition: QuEST.h:183
int isDensityMatrix
Whether this instance is a density-state representation.
Definition: QuEST.h:206
int numQubits
Definition: QuEST.h:138
void updateIndices(int oldEl, int newEl, int *list1, int len1, int *list2, int len2)
Definition: utilities.cpp:276
std::vector< std::vector< qcomp > > QMatrix
A complex square matrix.
Definition: utilities.hpp:49
int numQubits
The number of qubits for which this Hamiltonian is defined.
Definition: QuEST.h:168
QVector getNormalised(QVector vec)
Returns an L2-normalised copy of vec, using Kahan summation for improved accuracy.
Definition: utilities.cpp:431
int numQubitsRepresented
The number of qubits represented in either the state-vector or density matrix.
Definition: QuEST.h:208
long long int numAmpsTotal
Total number of amplitudes, which are possibly distributed among machines.
Definition: QuEST.h:215
qreal * real
The real values of the 2^numQubits complex elements.
Definition: QuEST.h:189
qreal real
Definition: QuEST.h:105
QMatrix getIdentityMatrix(size_t dim)
Returns a dim-by-dim identity matrix.
Definition: utilities.cpp:151
qreal imag
Definition: QuEST.h:106
QMatrix getKroneckerProduct(QMatrix a, QMatrix b)
Returns the kronecker product of a and b, where a and b are square but possibly differently-sized com...
Definition: utilities.cpp:169
QMatrix getZeroMatrix(size_t dim)
Returns a dim-by-dim square complex matrix, initialised to all zeroes.
Definition: utilities.cpp:143
QVector getRandomQVector(int dim)
Returns a dim-length vector with random complex amplitudes in the square joining {-1-i,...
Definition: utilities.cpp:420
bool areEqual(QVector a, QVector b)
Returns true if the absolute value of the difference between every amplitude in vectors a and b is le...
Definition: utilities.cpp:387
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,...
Definition: utilities.cpp:573
#define DEMAND(cond)
Definition: utilities.cpp:24
#define macro_copyQMatrix(dest, src)
Definition: utilities.cpp:838
Represents a 2x2 matrix of complex numbers.
Definition: QuEST.h:114