Precision-agnostic types and data structures for representing numbers, quantum states, and operators. More...

Data Structures

struct  Complex
 Represents one complex number. More...
 
struct  ComplexMatrix2
 Represents a 2x2 matrix of complex numbers. More...
 
struct  ComplexMatrix4
 Represents a 4x4 matrix of complex numbers. More...
 
struct  ComplexMatrixN
 Represents a general 2^N by 2^N matrix of complex numbers. More...
 
struct  DiagonalOp
 Represents a diagonal complex operator on the full Hilbert state of a Qureg. More...
 
struct  PauliHamil
 A Pauli Hamiltonian, expressed as a real-weighted sum of pauli products, and which can hence represent any Hermitian operator. More...
 
struct  QuESTEnv
 Information about the environment the program is running in. More...
 
struct  Qureg
 Represents a system of qubits. More...
 
struct  Vector
 Represents a 3-vector of real numbers. More...
 

Macros

#define fromComplex(comp)   qcomp(comp.real, comp.imag)
 
#define getStaticComplexMatrixN(numQubits, re, im)
 Creates a ComplexMatrixN struct which lives in the stack and so does not need freeing, but cannot be returned beyond the calling scope. More...
 
#define qcomp
 
#define qreal   double
 
#define QuEST_PREC   2
 
#define toComplex(scalar)   ((Complex) {.real = creal(scalar), .imag = cimag(scalar)})
 

Enumerations

enum  bitEncoding { UNSIGNED =0, TWOS_COMPLEMENT =1 }
 Flags for specifying how the bits in sub-register computational basis states are mapped to indices in functions like applyPhaseFunc(). More...
 
enum  pauliOpType { PAULI_I =0, PAULI_X =1, PAULI_Y =2, PAULI_Z =3 }
 Codes for specifying Pauli operators. More...
 
enum  phaseFunc {
  NORM =0, SCALED_NORM =1, INVERSE_NORM =2, SCALED_INVERSE_NORM =3,
  SCALED_INVERSE_SHIFTED_NORM =4, PRODUCT =5, SCALED_PRODUCT =6, INVERSE_PRODUCT =7,
  SCALED_INVERSE_PRODUCT =8, DISTANCE =9, SCALED_DISTANCE =10, INVERSE_DISTANCE =11,
  SCALED_INVERSE_DISTANCE =12, SCALED_INVERSE_SHIFTED_DISTANCE =13
}
 Flags for specifying named phase functions. More...
 

Functions

Qureg createCloneQureg (Qureg qureg, QuESTEnv env)
 Create a new Qureg which is an exact clone of the passed qureg, which can be either a state-vector or a density matrix. More...
 
ComplexMatrixN createComplexMatrixN (int numQubits)
 Allocate dynamic memory for a square complex matrix of any size, which can be passed to functions like multiQubitUnitary() and applyMatrixN(). More...
 
Qureg createDensityQureg (int numQubits, QuESTEnv env)
 Creates a density matrix Qureg object representing a set of qubits which can enter noisy and mixed states. More...
 
DiagonalOp createDiagonalOp (int numQubits, QuESTEnv env)
 Creates a DiagonalOp representing a diagonal operator on the full Hilbert space of a Qureg. More...
 
DiagonalOp createDiagonalOpFromPauliHamilFile (char *fn, QuESTEnv env)
 Creates and initialiases a diagonal operator from the Z Pauli Hamiltonian encoded in file with filename fn. More...
 
PauliHamil createPauliHamil (int numQubits, int numSumTerms)
 Dynamically allocates a Hamiltonian expressed as a real-weighted sum of products of Pauli operators. More...
 
PauliHamil createPauliHamilFromFile (char *fn)
 Creates a PauliHamil instance, a real-weighted sum of products of Pauli operators, populated with the data in filename fn. More...
 
QuESTEnv createQuESTEnv (void)
 Create the QuEST execution environment. More...
 
Qureg createQureg (int numQubits, QuESTEnv env)
 Creates a state-vector Qureg object representing a set of qubits which will remain in a pure state. More...
 
void destroyComplexMatrixN (ComplexMatrixN matr)
 Destroy a ComplexMatrixN instance created with createComplexMatrixN() More...
 
void destroyDiagonalOp (DiagonalOp op, QuESTEnv env)
 Destroys a DiagonalOp created with createDiagonalOp(), freeing its memory. More...
 
void destroyPauliHamil (PauliHamil hamil)
 Destroy a PauliHamil instance, created with either createPauliHamil() or createPauliHamilFromFile(). More...
 
void destroyQuESTEnv (QuESTEnv env)
 Destroy the QuEST environment. More...
 
void destroyQureg (Qureg qureg, QuESTEnv env)
 Deallocate a Qureg, freeing its memory. More...
 
void initComplexMatrixN (ComplexMatrixN m, qreal real[][1<< m.numQubits], qreal imag[][1<< m.numQubits])
 Initialises a ComplexMatrixN instance to have the passed real and imag values. More...
 
void initDiagonalOp (DiagonalOp op, qreal *real, qreal *imag)
 Overwrites the entire DiagonalOp op with the given real and imag complex elements. More...
 
void initDiagonalOpFromPauliHamil (DiagonalOp op, PauliHamil hamil)
 Populates the diagonal operator op to be equivalent to the given Pauli Hamiltonian hamil, assuming hamil contains only PAULI_Z operators. More...
 
void initPauliHamil (PauliHamil hamil, qreal *coeffs, enum pauliOpType *codes)
 Initialise PauliHamil instance hamil with the given term coefficients and Pauli codes (one for every qubit in every term). More...
 
void setDiagonalOpElems (DiagonalOp op, long long int startInd, qreal *real, qreal *imag, long long int numElems)
 Modifies a subset (starting at index startInd, and ending at index startInd + numElems) of the elements in DiagonalOp op with the given complex numbers (passed as real and imag components). More...
 
void syncDiagonalOp (DiagonalOp op)
 Update the GPU memory with the current values in op.real and op.imag. More...
 

Detailed Description

Precision-agnostic types and data structures for representing numbers, quantum states, and operators.

Macro Definition Documentation

◆ fromComplex

#define fromComplex (   comp)    qcomp(comp.real, comp.imag)

Converts a Complex struct to a qcomp native type

Author
Tyson Jones

◆ getStaticComplexMatrixN

#define getStaticComplexMatrixN (   numQubits,
  re,
  im 
)
Value:
numQubits, \
(qreal[1<<numQubits][1<<numQubits]) UNPACK_ARR re, \
(qreal[1<<numQubits][1<<numQubits]) UNPACK_ARR im, \
(double*[1<<numQubits]) {NULL}, (double*[1<<numQubits]) {NULL} \
)

Creates a ComplexMatrixN struct which lives in the stack and so does not need freeing, but cannot be returned beyond the calling scope.

That is, the .real and .imag arrays of the returned ComplexMatrixN live in the stack as opposed to that returned by createComplexMatrixN() (which live in the heap). Note the real and imag components must be wrapped in paranthesis, e.g.

ComplexMatrixN u = getStaticComplexMatrixN(1, ({{1,2},{3,4}}), ({{0}}));

Here is an example of an incorrect usage, since a 'local' ComplexMatrixN cannot leave the calling scope (otherwise inducing dangling pointers):

ComplexMatrixN getMyMatrix(void) {
return getStaticComplexMatrixN(1, ({{1,2},{3,4}}), ({{0}}));
}

This function is actually a single-line anonymous macro, so can be safely invoked within arguments to other functions, e.g.

qureg, (int[]) {0}, 1,
getStaticComplexMatrixN(1, ({{1,0},{0,1}}), ({{0}}))
);

The returned ComplexMatrixN can be accessed and modified in the same way as that returned by createComplexMatrixN(), e.g.

ComplexMatrixN u = getStaticComplexMatrixN(3, ({{0}}), ({{0}}));
for (int i=0; i<8; i++)
for (int j=0; j<8; j++)
u.real[i][j] = .1;


Note that the first argument numQubits must be a literal.

This macro is only callable in C, since it invokes the function bindArraysToStackComplexMatrixN() which is only callable in C.

See also
Author
Tyson Jones

Definition at line 5521 of file QuEST.h.

◆ qcomp

#define qcomp

A precision-agnostic operator-overloaded complex number type.
This is a complex analog of qreal and is of single, double or quad precision depending on the value of QuEST_PREC. It resolves to the native complex type provided by <complex.h> for both C99 and C++11, so can be used with operators. It can be constructed with qcomp(real, imag).

For example, in C,

qcomp x = 2 + 3i;
x -= 3.2*x;

and in C++,

qcomp x = qcomp(2, 3);
x -= 3*x;

Assuming QuEST_PREC=4, qcomp will be 'complex long double' in C and 'complex<long double>' in C++.

Can be converted to/from Complex, the struct accepted by the QuEST interface, using toComplex and fromComplex.

Authors
Randy Meyers and Dr. Thomas Plum (created C & C++ agnosticism)
Author
Tyson Jones (created precision agnosticism)

◆ qreal

#define qreal   double

A precision-agnostic floating point number, as determined by QuEST_PREC. Is a single, double or quad precision float when QuEST_PREC is 1, 2 or 4 respectively.

Author
Ania Brown
Tyson Jones (doc)

◆ QuEST_PREC

#define QuEST_PREC   2

Sets the precision of qreal and qcomp floating-point numbers, and hence the numerical precision of Qureg.

QuEST_PREC can be set as a preprocessor macro during compilation, or by editing its definition in QuEST_precision.h.
The possible values are:

QuEST_PREC qreal sizeof(qreal)
1 float 4 bytes
2 double 8 bytes
4 long double 16 bytes

Note that quad precision (QuEST_PREC = 4 ) is not compatible with most GPUs.

See also
Author
Ania Brown
Tyson Jones (doc)

◆ toComplex

#define toComplex (   scalar)    ((Complex) {.real = creal(scalar), .imag = cimag(scalar)})

Creates a Complex struct, which can be passed to the QuEST API, from a qcomp

Author
Tyson Jones

Enumeration Type Documentation

◆ bitEncoding

Flags for specifying how the bits in sub-register computational basis states are mapped to indices in functions like applyPhaseFunc().

  • UNSIGNED means the bits encode an unsigned integer, hence

    \[ \begin{aligned} |00\rangle & \rightarrow \, 0 \\ |01\rangle & \rightarrow \, 1 \\ |10\rangle & \rightarrow \, 2 \\ |11\rangle & \rightarrow \, 3 \end{aligned} \]

  • TWOS_COMPLEMENT means the bits encode a signed integer through two's complement, such that

    \[ \begin{aligned} |000\rangle & \rightarrow \, 0 \\ |001\rangle & \rightarrow \, 1 \\ |010\rangle & \rightarrow \, 2 \\ |011\rangle & \rightarrow \, 3 \\ |100\rangle & \rightarrow \,-4 \\ |101\rangle & \rightarrow \,-3 \\ |110\rangle & \rightarrow \,-2 \\ |111\rangle & \rightarrow \,-1 \end{aligned} \]

    Remember that the qubits specified within a sub-register, and their ordering (least to most significant) determine the bits of a computational basis state, before intrepretation as an encoding of an integer.

Author
Tyson Jones
Enumerator
UNSIGNED 
TWOS_COMPLEMENT 

Definition at line 269 of file QuEST.h.

◆ pauliOpType

Codes for specifying Pauli operators.

Author
Tyson Jones
Enumerator
PAULI_I 
PAULI_X 
PAULI_Y 
PAULI_Z 

Definition at line 96 of file QuEST.h.

96 {PAULI_I=0, PAULI_X=1, PAULI_Y=2, PAULI_Z=3};

◆ phaseFunc

enum phaseFunc

Flags for specifying named phase functions.

These can be passed to functions applyNamedPhaseFunc(), applyNamedPhaseFuncOverrides(), applyParamNamedPhaseFunc(), and applyParamNamedPhaseFuncOverrides().

Norm based phase functions:

  • NORM maps state $|x\rangle|y\rangle\dots$ to $\sqrt{x^2 + y^2 + \dots}$
  • SCALED_NORM maps state $|x\rangle|y\rangle\dots$ to $\text{coeff} \sqrt{x^2 + y^2 + \dots}$
  • INVERSE_NORM maps state $|x\rangle|y\rangle\dots$ to $1/\sqrt{x^2 + y^2 + \dots}$
  • SCALED_INVERSE_NORM maps state $|x\rangle|y\rangle\dots$ to $\text{coeff}/\sqrt{x^2 + y^2 + \dots}$
  • SCALED_INVERSE_SHIFTED_NORM maps state $|x\rangle|y\rangle\dots$ to $\text{coeff}/\sqrt{(x-\Delta_x)^2 + (y-\Delta_y)^2 + \dots}$

Product based phase functions:

  • PRODUCT maps state $|x\rangle|y\rangle|z\rangle\dots$ to $x \; y \; z \dots$
  • SCALED_PRODUCT maps state $|x\rangle|y\rangle|z\rangle\dots$ to $\text{coeff} \; x \; y \; z \dots$
  • INVERSE_PRODUCT maps state $|x\rangle|y\rangle|z\rangle\dots$ to $1/(x \; y \; z \dots)$
  • SCALED_INVERSE_PRODUCT maps state $|x\rangle|y\rangle|z\rangle\dots$ to $\text{coeff}/(x \; y \; z \dots)$

Euclidean distance based phase functions:

  • DISTANCE maps state $|x_1\rangle|x_2\rangle|y_1\rangle|y_2\rangle\dots$ to $\sqrt{(x_1-x_2)^2 + (y_1-y_2)^2 + \dots}$
  • SCALED_DISTANCE maps state $|x_1\rangle|x_2\rangle|y_1\rangle|y_2\rangle\dots$ to $\text{coeff}\sqrt{(x_1-x_2)^2 + (y_1-y_2)^2 + \dots}$
  • INVERSE_DISTANCE maps state $|x_1\rangle|x_2\rangle|y_1\rangle|y_2\rangle\dots$ to $1/\sqrt{(x_1-x_2)^2 + (y_1-y_2)^2 + \dots}$
  • SCALED_INVERSE_DISTANCE maps state $|x_1\rangle|x_2\rangle|y_1\rangle|y_2\rangle\dots$ to $\text{coeff}/\sqrt{(x_1-x_2)^2 + (y_1-y_2)^2 + \dots}$
  • SCALED_INVERSE_SHIFTED_DISTANCE maps state $|x_1\rangle|x_2\rangle|y_1\rangle|y_2\rangle\dots$ to $\text{coeff}/\sqrt{(x_1-x_2-\Delta_x)^2 + (y_1-y_2-\Delta_y)^2 + \dots}$
Author
Tyson Jones
Richard Meister (shifted functions)
Enumerator
NORM 
SCALED_NORM 
INVERSE_NORM 
SCALED_INVERSE_NORM 
SCALED_INVERSE_SHIFTED_NORM 
PRODUCT 
SCALED_PRODUCT 
INVERSE_PRODUCT 
SCALED_INVERSE_PRODUCT 
DISTANCE 
SCALED_DISTANCE 
INVERSE_DISTANCE 
SCALED_INVERSE_DISTANCE 
SCALED_INVERSE_SHIFTED_DISTANCE 

Definition at line 231 of file QuEST.h.

Function Documentation

◆ createCloneQureg()

Qureg createCloneQureg ( Qureg  qureg,
QuESTEnv  env 
)

Create a new Qureg which is an exact clone of the passed qureg, which can be either a state-vector or a density matrix.

The returned Qureg will have the same dimensions as the passed qureg and begin in an identical quantum state. This must be destroyed by the user later with destroyQureg().

See also
Returns
an object representing the set of qubits
Parameters
[in]quregan existing Qureg to be cloned
[in]envthe QuESTEnv
Author
Tyson Jones

Definition at line 64 of file QuEST.c.

64  {
65 
66  Qureg newQureg;
67  statevec_createQureg(&newQureg, qureg.numQubitsInStateVec, env);
68  newQureg.isDensityMatrix = qureg.isDensityMatrix;
71 
72  qasm_setup(&newQureg);
73  statevec_cloneQureg(newQureg, qureg);
74  return newQureg;
75 }

References Qureg::isDensityMatrix, Qureg::numQubitsInStateVec, Qureg::numQubitsRepresented, qasm_setup(), statevec_cloneQureg(), and statevec_createQureg().

Referenced by TEST_CASE().

◆ createComplexMatrixN()

ComplexMatrixN createComplexMatrixN ( int  numQubits)

Allocate dynamic memory for a square complex matrix of any size, which can be passed to functions like multiQubitUnitary() and applyMatrixN().

The returned matrix will have dimensions

\[ 2^{\text{numQubits}} \times 2^{\text{numQubits}}, \]

stored as nested arrays ComplexMatrixN.real and ComplexMatrixN.imag, initialised to zero.

Unlike a Qureg, the memory of a ComplexMatrixN is always stored in RAM, and non-distributed. Hence, elements can be directly accessed and modified:

int numQubits = 5;
int dim = (1 << numQubits);
for (int r=0; r<dim; r++) {
for (int c=0; c<dim; c++) {
m.real[r][c] = rand();
m.imag[r][c] = rand();
}
}


A ComplexMatrixN can be initialised in bulk using initComplexMatrixN(), though this is not C++ compatible.

Like ComplexMatrix2 and ComplexMatrix4 (which are incidentally stored in the stack), the returned ComplexMatrixN is safe to return from functions.

The ComplexMatrixN must eventually be freed using destroyComplexMatrixN(), since it is created in the dynamic heap. One can instead use getStaticComplexMatrixN() to create a ComplexMatrixN struct in the stack (which doesn't need to be later destroyed), though this may cause a stack overflow if the matrix is too large (approx 10+ qubits).

See also
Parameters
[in]numQubitsthe number of qubits of which the returned ComplexMatrixN will correspond
Returns
a dynamic ComplexMatrixN struct, that is one where the .real and .imag fields are arrays kept in the heap and must be later destroyed.
Exceptions
invalidQuESTInputError()
  • if numQubits <= 0
  • if the memory was not allocated successfully
Author
Tyson Jones

Definition at line 1348 of file QuEST.c.

1348  {
1349  validateNumQubitsInMatrix(numQubits, __func__);
1350 
1351  int numRows = 1 << numQubits;
1352 
1353  ComplexMatrixN m = {
1354  .numQubits = numQubits,
1355  .real = malloc(numRows * sizeof *m.real),
1356  .imag = malloc(numRows * sizeof *m.imag)};
1357 
1358  for (int n=0; n < 1<<numQubits; n++) {
1359  m.real[n] = calloc(numRows, sizeof **m.real);
1360  m.imag[n] = calloc(numRows, sizeof **m.imag);
1361  }
1362 
1363  // error if the ComplexMatrixN was not successfully malloc'ds
1364  validateMatrixInit(m, __func__);
1365 
1366  return m;
1367  }

References ComplexMatrixN::imag, ComplexMatrixN::numQubits, ComplexMatrixN::real, validateMatrixInit(), and validateNumQubitsInMatrix().

Referenced by densmatr_mixMultiQubitKrausMap(), densmatr_mixTwoQubitKrausMap(), and TEST_CASE().

◆ createDensityQureg()

Qureg createDensityQureg ( int  numQubits,
QuESTEnv  env 
)

Creates a density matrix Qureg object representing a set of qubits which can enter noisy and mixed states.

Allocates space for a matrix of complex amplitudes, which assuming a single qreal floating-point number requires qrealBytes, requires memory

\[ \text{qrealBytes} \times 2 \times 2^{2 \times\text{numQubits}}\;\;\text{(bytes)}, \]

though there are additional memory costs in GPU and distributed modes. Notice this is the memory cost of a state-vector created with createQureg() of twice as many qubits.

The returned Qureg begins in the zero state, as produced by initZeroState().

Once created, the following Qureg fields are relevant in all backends:

Behind the scenes, density matrice are stored as state-vectors, flattened column-wise. As such, individual amplitudes should be fetched with getDensityAmp(), in lieu of direct access.

QuESTEnv env must be prior created with createQuESTEnv().

Serial

In serial and local (non-distributed) multithreaded modes, a density matrix Qureg costs only the memory above. For example, at double precision (QuEST_PREC = 2, qrealBytes = 8), the memory costs are:

numQubits memory
10 16 MiB
12 256 MiB
14 4 GiB
16 64 GiB
18 1 TiB
20 16 TiB

GPU

In GPU-accelerated mode, an additional density matrix is created in GPU memory. Therefore both RAM and VRAM must be of sufficient memory to store the state-vector, each of the size indicated in the Serial table above.

Note that many GPUs do not support quad precision qreal.

Distributed

In distributed mode, the density matrix is uniformly partitioned between the N distributed nodes (column-wise).

Only a power-of-2 number of nodes N may be used (e.g. N = 1, 2, 4, 8, ...). There must additionally be at least 1 amplitude of a density matrix stored on each node. This means one cannot create a density matrix Qureg with fewer than $\log_2(\text{N})/2$ qubits.

Additional memory is allocated on each node for communication buffers, of size equal to the density matrix partition. Hence the total memory per-node required is:

\[ 2 \times \text{qrealBytes} \times 2 \times 2^{2\times\text{numQubits}}/N \;\;\text{(bytes)}, \]

For example, at double precision (QuEST_PREC = 2, qrealBytes = 8), the memory costs are:

numQubits memory per node
N = 2 N = 4 N = 8 N = 16 N = 32
10 16 MiB 8 MiB 4 MiB 2 MiB 1 MiB
15 16 GiB 8 GiB 4 GiB 2 GiB 1 GiB
20 16 TiB 8 TiB 4 TiB 2 TiB 1 TiB
See also
Returns
an object representing the set of qubits
Parameters
[in]numQubitsnumber of qubits in the system
[in]envobject representing the execution environment (local, multinode etc)
Exceptions
invalidQuESTInputError()
  • if numQubits <= 0
  • if numQubits is so large that the number of amplitudes cannot fit in a long long int type,
  • if in distributed mode, there are more nodes than elements in the would-be state-vector
exit
  • if in GPU mode, but GPU memory cannot be allocated.
Author
Tyson Jones

Definition at line 50 of file QuEST.c.

50  {
51  validateNumQubitsInQureg(2*numQubits, env.numRanks, __func__);
52 
53  Qureg qureg;
54  statevec_createQureg(&qureg, 2*numQubits, env);
55  qureg.isDensityMatrix = 1;
56  qureg.numQubitsRepresented = numQubits;
57  qureg.numQubitsInStateVec = 2*numQubits;
58 
59  qasm_setup(&qureg);
60  initZeroState(qureg); // safe call to public function
61  return qureg;
62 }

References initZeroState(), Qureg::isDensityMatrix, Qureg::numQubitsInStateVec, Qureg::numQubitsRepresented, QuESTEnv::numRanks, qasm_setup(), statevec_createQureg(), and validateNumQubitsInQureg().

Referenced by TEST_CASE().

◆ createDiagonalOp()

DiagonalOp createDiagonalOp ( int  numQubits,
QuESTEnv  env 
)

Creates a DiagonalOp representing a diagonal operator on the full Hilbert space of a Qureg.

The resulting operator need not be unitary nor Hermitian, and can be applied to any Qureg of a compatible number of qubits.

This function allocates space for $2^{\text{numQubits}}$ complex amplitudes, which are initially zero. This is the same cost as a local state-vector of equal number of qubits; see the Serial section of createQureg(). Note that this is a paralell data-type, so its ultimate memory costs depend on the hardware backends, as elaborated below.

The operator elements should be modified with initDiagonalOp() and setDiagonalOpElems(), and must be later freed with destroyDiagonalOp().

GPU

In GPU-accelerated mode, this function also creates additional equally-sized persistent memory on the GPU. If you wish to modify the operator elements directly (in lieu of setDiagonalOpElems()), you must thereafter call syncDiagonalOp() to update the operator stored in VRAM.

For example,

for (long long int i=0; i<op.numElemsPerChunk; i++) {
op.real[i] = rand();
op.imag[i] = rand();
}

Distribution

In distributed mode, the memory for the diagonal operator is divided evenly between the $N$ available nodes, such that each node contains only $2^{\text{numQubits}}/N$ complex values. This is assigned to DiagonalOp.numElemsPerChunk.

Users must therefore exercise care in modifying DiagonalOp.real and DiagonalOp.imag directly.

For example, the following is valid code when when distributed between N = 2 nodes:

// create diag({1,2,3,4,5,6,7,8, 9,10,11,12,13,14,15,16})
int numQubits = 4;
DiagonalOp op = createDiagonalOp(numQubits4, env);
for (int i=0; i<8; i++) {
if (env.rank == 0)
op.real[i] = (i+1);
if (env.rank == 1)
op.real[i] = (i+1+8);
}


See also
Returns
a dynamic DiagonalOp instance initialised to diag(0,0,...).
Parameters
[in]numQubitsnumber of qubits which inform the Hilbert dimension of the operator.
[in]envthe QuESTEnv
Exceptions
invalidQuESTInputError()
  • if numQubits <= 0
  • if numQubits is so large that the number of elements cannot fit in a long long int type,
  • if in distributed mode, there are more nodes than elements in the operator
exit
  • if the memory could not be allocated
Author
Tyson Jones

Definition at line 1518 of file QuEST.c.

1518  {
1519  validateNumQubitsInDiagOp(numQubits, env.numRanks, __func__);
1520 
1521  return agnostic_createDiagonalOp(numQubits, env);
1522 }

References agnostic_createDiagonalOp(), QuESTEnv::numRanks, and validateNumQubitsInDiagOp().

Referenced by TEST_CASE().

◆ createDiagonalOpFromPauliHamilFile()

DiagonalOp createDiagonalOpFromPauliHamilFile ( char *  fn,
QuESTEnv  env 
)

Creates and initialiases a diagonal operator from the Z Pauli Hamiltonian encoded in file with filename fn.

This is a convenience function to prepare a diagonal operator from a plaintext description of an all-Z Pauli Hamiltonian. The returned DiagonalOp is a distributed data structure, and significantly faster to use (through functions like calcExpecDiagonalOp()) than PauliHamil functions (like calcExpecPauliHamil()).

The returned DiagonalOp must be later freed with destroyDiagonalOp().
Note a PauliHamil from fn is temporarily internally created.

This function is equivalent to


See also
Parameters
[in]fnfilename of a plaintext file encoding an all-Z Pauli Hamiltonian
[in]envthe session QuESTEnv
Returns
a created DiagonalOp equivalent to the Hamiltonian in fn
Exceptions
invalidQuESTInputError()
  • if file fn cannot be read
  • if file fn does not encode a valid PauliHamil
  • if the encoded PauliHamil consists of operators other than PAULI_Z and PAULI_I
Author
Tyson Jones
Milos Prokop (serial prototype)

Definition at line 1558 of file QuEST.c.

1558  {
1559  PauliHamil h = createPauliHamilFromFile(fn); // validates fn
1560  validateDiagPauliHamilFromFile(h, env.numRanks, __func__); // destroys h if invalid
1561 
1564 
1565  destroyPauliHamil(h);
1566  return op;
1567 }

References agnostic_createDiagonalOp(), agnostic_initDiagonalOpFromPauliHamil(), createPauliHamilFromFile(), destroyPauliHamil(), PauliHamil::numQubits, QuESTEnv::numRanks, and validateDiagPauliHamilFromFile().

Referenced by TEST_CASE().

◆ createPauliHamil()

PauliHamil createPauliHamil ( int  numQubits,
int  numSumTerms 
)

Dynamically allocates a Hamiltonian expressed as a real-weighted sum of products of Pauli operators.

A PauliHamil is merely an encapsulation of the multiple parameters of functions like applyPauliSum().

The Pauli operators (PauliHamil.pauliCodes) are all initialised to identity (PAULI_I), but the coefficients (PauliHamil.termCoeffs) are not initialised.

The Hamiltonian can be used in functions like applyPauliHamil() and applyTrotterCircuit(), with Qureg instances of the same number of qubits.

A PauliHamil can be modified directly (see PauliHamil), or through initPauliHamil(). It can furthermore be created and initialised from a file description directly with createPauliHamilFromFile().

The returned dynamic PauliHamil instance must later be freed via destroyPauliHamil().

See also
Parameters
[in]numQubitsthe number of qubits on which this Hamiltonian acts
[in]numSumTermsthe number of weighted terms in the sum, or the number of Pauli products
Returns
a dynamic PauliHamil struct, with fields pauliCodes and termCoeffs stored in the heap
Exceptions
invalidQuESTInputError()
  • if numQubits <= 0
  • if numSumTerms <= 0
Author
Tyson Jones

Definition at line 1398 of file QuEST.c.

1398  {
1399  validateHamilParams(numQubits, numSumTerms, __func__);
1400 
1401  PauliHamil h;
1402  h.numQubits = numQubits;
1403  h.numSumTerms = numSumTerms;
1404  h.termCoeffs = malloc(numSumTerms * sizeof *h.termCoeffs);
1405  h.pauliCodes = malloc(numQubits*numSumTerms * sizeof *h.pauliCodes);
1406 
1407  // initialise pauli codes to identity
1408  for (int i=0; i<numQubits*numSumTerms; i++)
1409  h.pauliCodes[i] = PAULI_I;
1410 
1411  return h;
1412 }

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

Referenced by createPauliHamilFromFile(), and TEST_CASE().

◆ createPauliHamilFromFile()

PauliHamil createPauliHamilFromFile ( char *  fn)

Creates a PauliHamil instance, a real-weighted sum of products of Pauli operators, populated with the data in filename fn.

Each line in the plaintext file is interpreted as a separate product of Pauli operators in the total sum. Each line must be a space-separated list with format

c p1 p2 p3 ... pN

where c is the real coefficient of the term, and p1 ... pN are numbers in {0,1,2,3} to indicate PAULI_I, PAULI_X, PAULI_Y, PAULI_Z operators respectively, which act on qubits 0 through N-1 (all qubits).

For example, the file containing

0.31 1 0 1 2
-0.2 3 2 0 0

encodes a two-term four-qubit Hamiltonian

\[ 0.31 \, X_0 \, X_2 \, Y_3 -0.2 \, Z_0 \, Y_1 \,. \]

The initialised PauliHamil can be previewed with reportPauliHamil().

The number of qubits and terms are inferred from the file. The created Hamiltonian can be used just like one created via createPauliHamil(). It can be modified directly (see PauliHamil), or through initPauliHamil().

The returned dynamic PauliHamil instance must later be freed via destroyPauliHamil().

See also
Parameters
[in]fnfilename of the plaintext file specifying the pauli operators and coefficients
Returns
a dynamic PauliHamil struct
Exceptions
invalidQuESTInputError()
  • if the file with name fn cannot be read
  • if the file is not correctly formatted as described above
Author
Tyson Jones

Definition at line 1420 of file QuEST.c.

1420  {
1421 
1422  /* The validation in this function must close the file handle and free
1423  * allocated memory before raising an error (whether that's a C exit, or
1424  * an overriden C++ exception).
1425  */
1426 
1427  FILE* file = fopen(fn, "r");
1428  int success = (file != NULL);
1429  validateFileOpened(success, fn, __func__);
1430 
1431  /* file format: coeff {term} \n where {term} is #numQubits values of
1432  * 0 1 2 3 signifying I X Y Z acting on that qubit index
1433  */
1434 
1435  // count the number of qubits (ignore trailing whitespace)
1436  int numQubits = -1; // to exclude coeff at start
1437  char ch = getc(file);
1438  char prev = '0'; // anything not space
1439  while (ch != '\n' && ch != EOF) {
1440  if (ch == ' ' && prev != ' ') // skip multiple spaces
1441  numQubits++;
1442  prev = ch;
1443  ch = getc(file);
1444  }
1445  // edge-case: if we hit EOF/newline without a space
1446  if (prev != ' ')
1447  numQubits++;
1448 
1449  /* TODO:
1450  * The below code may break on Windows where newlines are multiple characters
1451  */
1452 
1453  // count the number of terms (being cautious of trailing newlines)
1454  int numTerms = 0;
1455  prev = '\n';
1456  rewind(file);
1457  while ((ch=getc(file)) != EOF) {
1458  if (ch == '\n' && prev != '\n')
1459  numTerms++;
1460  prev = ch;
1461  }
1462  // edge-case: if we hit EOF without a newline, count that line
1463  if (prev != '\n')
1464  numTerms++;
1465 
1466  // validate the inferred number of terms and qubits (closes file if error)
1467  validateHamilFileParams(numQubits, numTerms, file, fn, __func__);
1468 
1469  // allocate space for PauliHamil data
1470  PauliHamil h = createPauliHamil(numQubits, numTerms);
1471 
1472  // specifier for a qreal number then a space
1473  char strSpec[50];
1474  strcpy(strSpec, REAL_SPECIFIER);
1475  strcat(strSpec, " ");
1476 
1477  // collect coefficients and terms
1478  rewind(file);
1479  for (int t=0; t<numTerms; t++) {
1480 
1481  // record coefficient, and validate (closes file and frees h if error)
1482  success = fscanf(file, strSpec, &(h.termCoeffs[t])) == 1;
1483  validateHamilFileCoeffParsed(success, h, file, fn, __func__);
1484 
1485  // record Pauli operations, and validate (closes file and frees h if error)
1486  for (int q=0; q<numQubits; q++) {
1487  int i = t*numQubits + q;
1488 
1489  // verbose, to avoid type warnings
1490  int code;
1491  success = fscanf(file, "%d ", &code) == 1;
1492  h.pauliCodes[i] = (enum pauliOpType) code;
1493  validateHamilFilePauliParsed(success, h, file, fn, __func__);
1494  validateHamilFilePauliCode(h.pauliCodes[i], h, file, fn, __func__);
1495  }
1496 
1497  // the trailing newline is magically eaten
1498  }
1499 
1500  fclose(file);
1501  return h;
1502 }

References createPauliHamil(), PauliHamil::pauliCodes, PauliHamil::termCoeffs, validateFileOpened(), validateHamilFileCoeffParsed(), validateHamilFileParams(), validateHamilFilePauliCode(), and validateHamilFilePauliParsed().

Referenced by createDiagonalOpFromPauliHamilFile(), and TEST_CASE().

◆ createQuESTEnv()

QuESTEnv createQuESTEnv ( void  )

Create the QuEST execution environment.

This should be called only once, and the environment should be freed with destroyQuESTEnv at the end of the user's code. If something needs to be done to set up the execution environment, such as initializing MPI when running in distributed mode, it is handled here.

See also
Returns
object representing the execution environment. A single instance is used for each program
Author
Ania Brown

Definition at line 129 of file QuEST_cpu_distributed.c.

129  {
130 
131  QuESTEnv env;
132 
133  // init MPI environment
134  int rank, numRanks, initialized;
135  MPI_Initialized(&initialized);
136  if (!initialized){
137  MPI_Init(NULL, NULL);
138  MPI_Comm_size(MPI_COMM_WORLD, &numRanks);
139  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
140 
141  env.rank=rank;
142  env.numRanks=numRanks;
143 
144  } else {
145 
146  printf("ERROR: Trying to initialize QuESTEnv multiple times. Ignoring...\n");
147 
148  // ensure env is initialised anyway, so the compiler is happy
149  MPI_Comm_size(MPI_COMM_WORLD, &numRanks);
150  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
151  env.rank=rank;
152  env.numRanks=numRanks;
153  }
154 
155  validateNumRanks(env.numRanks, __func__);
156 
157  env.seeds = NULL;
158  env.numSeeds = 0;
159  seedQuESTDefault(&env);
160 
161  return env;
162 }

References GPUExists(), QuESTEnv::numRanks, QuESTEnv::numSeeds, QuESTEnv::rank, seedQuESTDefault(), QuESTEnv::seeds, and validateNumRanks().

Referenced by main().

◆ createQureg()

Qureg createQureg ( int  numQubits,
QuESTEnv  env 
)

Creates a state-vector Qureg object representing a set of qubits which will remain in a pure state.

Allocates space for a state-vector of complex amplitudes, which assuming a single qreal floating-point number requires qrealBytes, requires memory

\[ \text{qrealBytes} \times 2 \times 2^\text{numQubits}\;\;\text{(bytes)}, \]

though there are additional memory costs in GPU and distributed modes.

The returned Qureg begins in the zero state, as produced by initZeroState().

Once created, the following Qureg fields are relevant in all backends:

QuESTEnv env must be prior created with createQuESTEnv().

Serial

In serial and local (non-distributed) multithreaded modes, a state-vector Qureg costs only the memory above. For example, at double precision (QuEST_PREC = 2, qrealBytes = 8), the memory costs are:

numQubits memory
10 16 KiB
16 1 MiB
20 16 MiB
26 1 GiB
30 16 GiB

Individual amplitudes should be fetched and modified with functions like getAmp() and setAmps(). However, it is sometimes useful to access the state-vector directly, for example to create your own low-level (high performance) multithreaded functions. In those instants, Qureg.stateVec can be accessed directly, storing the real and imaginary components of the state-vector amplitudes in:

  • Qureg.stateVec.real
  • Qureg.stateVec.imag

The total number of amplitudes in the state-vector is

For example,

Qureg qureg = createQureg(10, env);
// ruin a perfectly good state-vector
for (long long int i=0; i<qureg.numAmpsTotal; i++) {
qureg.stateVec.real[i] = rand();
qureg.stateVec.imag[i] = rand();
}


GPU

In GPU-accelerated mode, an additional state-vector is created in GPU memory. Therefore both RAM and VRAM must be of sufficient memory to store the state-vector, each of the size indicated in the Serial table above.

Note that many GPUs do not support quad precision qreal.

Individual amplitudes of the created Qureg should be fetched and modified with functions like getAmp() and setAmps(). This is especially important since the GPU state-vector can be accessed directly, and changes to Qureg.stateVec will be ignored and overwritten. To modify the state-vector "directly", one must use copyStateFromGPU() and copyStateToGPU() before and after.

For example,

Qureg qureg = createQureg(10, env);
// ruin a perfectly good state-vector
for (long long int i=0; i<qureg.numAmpsTotal; i++) {
qureg.stateVec.real[i] = rand();
qureg.stateVec.imag[i] = rand();
}


Distributed

In distributed mode, the state-vector is uniformly partitioned between the N distributed nodes.

Only a power-of-2 number of nodes N may be used (e.g. N = 1, 2, 4, 8, ...). There must additionally be at least 1 amplitude of a state-vector stored on each node. This means one cannot create a state-vector Qureg with fewer than $\log_2(\text{N})$ qubits.

In addition to Qureg.stateVec, additional memory is allocated on each node for communication buffers, of size equal to the state-vector partition. Hence the total memory per-node required is:

\[ 2 \times \text{qrealBytes} \times 2 \times 2^\text{numQubits}/N \;\;\text{(bytes)}, \]

For example, at double precision (QuEST_PREC = 2, qrealBytes = 8), the memory costs are:

numQubits memory per node
N = 2 N = 4 N = 8 N = 16 N = 32
10 16 KiB 8 KiB 4 KiB 2 KiB 1 KiB
20 16 MiB 8 MiB 4 MiB 2 MiB 1 MiB
30 16 GiB 8 GiB 4 GiB 2 GiB 1 GiB
40 16 TiB 8 TiB 4 TiB 2 TiB 1 TiB

State-vector amplitudes should be set and modified using getAmp() and setAmps(). Direct modification is possible, but should be done extremely carefully, since each node only stores a partition of the full state-vector, which itself mightn't fit on any single node. Furthermore, an asynchronous MPI process may may unexpectedly modify local amplitudes; avoid this with syncQuESTEnv().

The fields relevant to distribution are:

Therefore, this code is valid

// set state |i> to have amplitude i
for (long long int i=0; i<qureg.numAmpsPerChunk; i++)
qureg.stateVec.real[i] = i + qureg.chunkId * qureg.numAmpsPerChunk;

while the following erroneous code would cause a segmentation fault:

// incorrectly attempt to set state |i> to have amplitude i
for (long long int i=0; i<qureg.numAmpsTotal; i++)
qureg.stateVec.real[i] = i;


See also
Returns
an object representing the set of qubits
Parameters
[in]numQubitsnumber of qubits in the system
[in]envobject representing the execution environment (local, multinode etc)
Exceptions
invalidQuESTInputError()
  • if numQubits <= 0
  • if numQubits is so large that the number of amplitudes cannot fit in a long long int type,
  • if in distributed mode, there are more nodes than elements in the would-be state-vector
exit
  • if in GPU mode, but GPU memory cannot be allocated.
Author
Ania Brown
Tyson Jones (validation, doc)

Definition at line 36 of file QuEST.c.

36  {
37  validateNumQubitsInQureg(numQubits, env.numRanks, __func__);
38 
39  Qureg qureg;
40  statevec_createQureg(&qureg, numQubits, env);
41  qureg.isDensityMatrix = 0;
42  qureg.numQubitsRepresented = numQubits;
43  qureg.numQubitsInStateVec = numQubits;
44 
45  qasm_setup(&qureg);
46  initZeroState(qureg); // safe call to public function
47  return qureg;
48 }

References initZeroState(), Qureg::isDensityMatrix, Qureg::numQubitsInStateVec, Qureg::numQubitsRepresented, QuESTEnv::numRanks, qasm_setup(), statevec_createQureg(), and validateNumQubitsInQureg().

Referenced by TEST_CASE().

◆ destroyComplexMatrixN()

void destroyComplexMatrixN ( ComplexMatrixN  matr)

Destroy a ComplexMatrixN instance created with createComplexMatrixN()

It is invalid to attempt to destroy a matrix created with getStaticComplexMatrixN().

See also
Parameters
[in]matrthe dynamic matrix (created with createComplexMatrixN()) to deallocate
Exceptions
invalidQuESTInputError()
  • if matr was not yet allocated.
malloc_error
Author
Tyson Jones

Definition at line 1369 of file QuEST.c.

1369  {
1370  /* this checks m.real/imag != NULL, which is only ever set when the mallocs
1371  * in createComplexMatrixN fail, which already prompts an error. Hence
1372  * this check if useless
1373  */
1374  validateMatrixInit(m, __func__);
1375 
1376  int numRows = 1 << m.numQubits;
1377  for (int r=0; r < numRows; r++) {
1378  free(m.real[r]);
1379  free(m.imag[r]);
1380  }
1381  free(m.real);
1382  free(m.imag);
1383 }

References ComplexMatrixN::imag, ComplexMatrixN::numQubits, ComplexMatrixN::real, and validateMatrixInit().

Referenced by densmatr_mixMultiQubitKrausMap(), densmatr_mixTwoQubitKrausMap(), and TEST_CASE().

◆ destroyDiagonalOp()

void destroyDiagonalOp ( DiagonalOp  op,
QuESTEnv  env 
)

Destroys a DiagonalOp created with createDiagonalOp(), freeing its memory.

See also
Parameters
[in]opthe DiagonalOp to destroy
[in]envthe QuESTEnv
Exceptions
invalidQuESTInputError()
  • if op was not previously created
Author
Tyson Jones

Definition at line 1524 of file QuEST.c.

1524  {
1525  // env accepted for API consistency
1526  validateDiagOpInit(op, __func__);
1527 
1529 }

References agnostic_destroyDiagonalOp(), and validateDiagOpInit().

Referenced by TEST_CASE().

◆ destroyPauliHamil()

void destroyPauliHamil ( PauliHamil  hamil)

Destroy a PauliHamil instance, created with either createPauliHamil() or createPauliHamilFromFile().

Parameters
[in]hamila dynamic PauliHamil instantiation
Author
Tyson Jones

Definition at line 1414 of file QuEST.c.

1414  {
1415 
1416  free(h.termCoeffs);
1417  free(h.pauliCodes);
1418 }

References PauliHamil::pauliCodes, and PauliHamil::termCoeffs.

Referenced by createDiagonalOpFromPauliHamilFile(), TEST_CASE(), validateDiagPauliHamilFromFile(), validateHamilFileCoeffParsed(), validateHamilFilePauliCode(), and validateHamilFilePauliParsed().

◆ destroyQuESTEnv()

void destroyQuESTEnv ( QuESTEnv  env)

Destroy the QuEST environment.

If something needs to be done to clean up the execution environment, such as finalizing MPI when running in distributed mode, it is handled here

See also
Parameters
[in]envobject representing the execution environment. A single instance is used for each program
Author
Ania Brown

Definition at line 174 of file QuEST_cpu_distributed.c.

174  {
175  free(env.seeds);
176 
177  int finalized;
178  MPI_Finalized(&finalized);
179  if (!finalized) MPI_Finalize();
180  else printf("ERROR: Trying to close QuESTEnv multiple times. Ignoring\n");
181 }

References QuESTEnv::seeds.

Referenced by main().

◆ destroyQureg()

void destroyQureg ( Qureg  qureg,
QuESTEnv  env 
)

Deallocate a Qureg, freeing its memory.

This frees all memory bound to qureg, including its state-vector or density matrix in RAM, in VRAM (in GPU mode), and communication buffers (in distributed mode).

The qureg must have been previously created with createQureg(), createDensityQureg() or createCloneQureg().

See also
Parameters
[in,out]quregthe Qureg to be destroyed
[in]envthe QuESTEnv
Author
Ania Brown
Tyson Jones (improved doc)

Definition at line 77 of file QuEST.c.

77  {
78  statevec_destroyQureg(qureg, env);
79  qasm_free(qureg);
80 }

References qasm_free(), and statevec_destroyQureg().

Referenced by TEST_CASE().

◆ initComplexMatrixN()

void initComplexMatrixN ( ComplexMatrixN  m,
qreal  real[][1<< m.numQubits],
qreal  imag[][1<< m.numQubits] 
)

Initialises a ComplexMatrixN instance to have the passed real and imag values.

This allows succint population of any-sized ComplexMatrixN, e.g. through 2D arrays:

ComplexMatrixN m = createComplexMatrixN(3);
initComplexMatrixN(m, 
    (qreal[8][8]) {{1,2,3,4,5,6,7,8}, {0}},
    (qreal[8][8]) {{0}});

m can be created by either createComplexMatrixN() or getStaticComplexMatrixN().

This function is only callable in C, since C++ signatures cannot contain variable-length 2D arrays

Parameters
[in]mthe matrix to initialise
[in]realmatrix of real values; can be 2D array of array of pointers
[in]imagmatrix of imaginary values; can be 2D array of array of pointers
Exceptions
invalidQuESTInputError()
Author
Tyson Jones

Definition at line 1386 of file QuEST.c.

1386  {
1387  validateMatrixInit(m, __func__);
1388 
1389  int dim = 1 << m.numQubits;
1390  for (int i=0; i<dim; i++)
1391  for (int j=0; j<dim; j++) {
1392  m.real[i][j] = re[i][j];
1393  m.imag[i][j] = im[i][j];
1394  }
1395 }

References ComplexMatrixN::imag, ComplexMatrixN::numQubits, ComplexMatrixN::real, and validateMatrixInit().

◆ initDiagonalOp()

void initDiagonalOp ( DiagonalOp  op,
qreal real,
qreal imag 
)

Overwrites the entire DiagonalOp op with the given real and imag complex elements.

Both real and imag must have length equal to pow(2, op.numQubits).

In GPU mode, this updates both the RAM (op.real and op.imag) and persistent GPU memory; there is no need to call syncDiagonalOp() afterward.

In distributed mode, this function assumes real and imag exist fully on every node. For DiagonalOp which are too large to fit into a single node, use setDiagonalOpElems() or syncDiagonalOp().

See also
Parameters
[in,out]opthe diagonal operator to modify
[in]realthe real components of the full set of new elements
[in]imagthe imaginary components of the full set of new elements
Exceptions
invalidQuESTInputError()
  • if op was not created
segmentation-fault
  • if either real or imag have length smaller than pow(2, op.numQubits)
Author
Tyson Jones

Definition at line 1537 of file QuEST.c.

1537  {
1538  validateDiagOpInit(op, __func__);
1539 
1540  agnostic_setDiagonalOpElems(op, 0, real, imag, 1LL << op.numQubits);
1541 }

References agnostic_setDiagonalOpElems(), DiagonalOp::numQubits, and validateDiagOpInit().

Referenced by TEST_CASE().

◆ initDiagonalOpFromPauliHamil()

void initDiagonalOpFromPauliHamil ( DiagonalOp  op,
PauliHamil  hamil 
)

Populates the diagonal operator op to be equivalent to the given Pauli Hamiltonian hamil, assuming hamil contains only PAULI_Z operators.

Given a PauliHamil hamil featuring only PAULI_Z and PAULI_I, with term coefficients $\{\lambda_j\}$, which hence has form

\[ \begin{aligned} \text{hamil} &= \sum\limits_j^{\text{numSumTerms}} \lambda_j \bigotimes\limits_{k_j} \hat{Z}_k \\ &\equiv \begin{pmatrix} r_1 \\ & r_2 \\ & & r_3 \\ & & & \ddots \\ & & & & r_{2^{\,\text{numQubits}}} \end{pmatrix}, \end{aligned} \]

this function modifies op to

\[ \text{op} \; \rightarrow \; \text{diag} \big( \; r_1, \; r_2, \; r_3, \; \dots, \; r_{2^{\,\text{numQubits}}} \, \big), \]

where the real amplitudes have form

\[ r_i = \sum\limits_j \, s_{ij} \, \lambda_j, \;\;\;\; s_{ij} = \pm 1 \,. \]

This is useful since calculations with DiagonalOp are significantly faster than the equivalent calculations with a general PauliHamil. For example, applyDiagonalOp() requires a factor numSumTerms * numQubits fewer operations than applyPauliHamil().

In distributed mode, each node will contain only a sub-partition of the full diagonal matrix.
In GPU mode, both the CPU and GPU memory copies of op will be updated, so there is no need to call syncDiagonalOp() afterward.

See also
Parameters
[in,out]opan existing DiagonalOp (e.g. created with createDiagonalOp()) to modify
[in]hamila PauliHamil of equal dimension to op, containing only PAULI_Z and PAULI_I operators
Exceptions
invalidQuESTInputError()
  • if hamil has invalid parameters (numQubits <= 0, numSumTerms <= 0)
  • if op and hamil have unequal dimensions
  • if hamil contains any operator other than PAULI_Z and PAULI_I
segmentation-fault
  • if either op or hamil have not been already created
Author
Tyson Jones
Milos Prokop (serial prototype)

Definition at line 1550 of file QuEST.c.

1550  {
1551  validateHamilParams(hamil.numQubits, hamil.numSumTerms, __func__);
1552  validateDiagOpInit(op, __func__);
1553  validateDiagPauliHamil(op, hamil, __func__);
1554 
1556 }

References agnostic_initDiagonalOpFromPauliHamil(), PauliHamil::numQubits, PauliHamil::numSumTerms, validateDiagOpInit(), validateDiagPauliHamil(), and validateHamilParams().

Referenced by TEST_CASE().

◆ initPauliHamil()

void initPauliHamil ( PauliHamil  hamil,
qreal coeffs,
enum pauliOpType codes 
)

Initialise PauliHamil instance hamil with the given term coefficients and Pauli codes (one for every qubit in every term).

Arguments coeffs and codes encode a weighted sum of Pauli operators, with the same format as other QuEST functions (like calcExpecPauliSum()).

This is useful to set the elements of the PauliHamil in batch.
For example

int numQubits = 3;
int numTerms = 2;
PauliHamil hamil = createPauliHamil(numQubits, numTerms);
// hamil = 0.5 X0 Y1 - 0.5 Z1 X3
(qreal[]) {0.5, -0.5},

The initialised PauliHamil can be previewed with reportPauliHamil().

hamil must be already created with createPauliHamil(), or createPauliHamilFromFile().

See also
Parameters
[in,out]hamilan existing PauliHamil instance to be modified
[in]coeffsan array of sum term coefficients, which must have length hamil.numSumTerms
[in]codesa flat array of Pauli codes, of length hamil.numSumTerms*hamil.numQubits
Exceptions
invalidQuESTInputError()
  • if hamil has invalid parameters (numQubits <= 0, numSumTerms <= 0)
  • if any code in codes is not a valid Pauli code (pauliOpType)
Author
Tyson Jones

Definition at line 1504 of file QuEST.c.

1504  {
1505  validateHamilParams(hamil.numQubits, hamil.numSumTerms, __func__);
1506  validatePauliCodes(codes, hamil.numSumTerms*hamil.numQubits, __func__);
1507 
1508  int i=0;
1509  for (int t=0; t<hamil.numSumTerms; t++) {
1510  hamil.termCoeffs[t] = coeffs[t];
1511  for (int q=0; q<hamil.numQubits; q++) {
1512  hamil.pauliCodes[i] = codes[i];
1513  i++;
1514  }
1515  }
1516 }

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

Referenced by TEST_CASE().

◆ setDiagonalOpElems()

void setDiagonalOpElems ( DiagonalOp  op,
long long int  startInd,
qreal real,
qreal imag,
long long int  numElems 
)

Modifies a subset (starting at index startInd, and ending at index startInd + numElems) of the elements in DiagonalOp op with the given complex numbers (passed as real and imag components).

In GPU mode, this updates both the RAM (op.real and op.imag), and the persistent GPU memory.

In distributed mode, this function assumes the subset real and imag exist (at least) on the node containing the ultimately updated elements.
For example, below is the correct way to modify the full 8 elements of op when split between 2 nodes.

int numElems = 4;
qreal re[] = {1,2,3,4};
qreal im[] = {1,2,3,4};
setDiagonalOpElems(op, 0, re, im, numElems);
// modify re and im to the next set of elements
for (int i=0; i<4; i++) {
re[i] += 4;
im[i] += 4;
}
setDiagonalOpElems(op, 4, re, im, numElems);

In this way, one can avoid a single node containing all new elements which might not fit. If more elements are passed than exist on an individual node, each node merely ignores the additional elements.

Parameters
[in,out]opthe DiagonalOp to modify
[in]startIndthe starting index (globally) of the subset of elements to modify
[in]realthe real components of the new elements
[in]imagthe imaginary components of the new elements
[in]numElemsthe number of new elements (the length of real and imag)
Exceptions
invalidQuESTInputError()
  • if op was not created
  • if startInd is an invalid index (<0 or >=pow(2,op.numQubits))
  • if numElems is an invalid number of elements (<=0 or >pow(2,op.numQubits))
  • if there are fewer than numElems elements in op after startInd
segmentation-fault
  • if either real or imag have fewer elements than numElems
Author
Tyson Jones

Definition at line 1543 of file QuEST.c.

1543  {
1544  validateDiagOpInit(op, __func__);
1545  validateNumElems(op, startInd, numElems, __func__);
1546 
1547  agnostic_setDiagonalOpElems(op, startInd, real, imag, numElems);
1548 }

References agnostic_setDiagonalOpElems(), validateDiagOpInit(), and validateNumElems().

Referenced by TEST_CASE().

◆ syncDiagonalOp()

void syncDiagonalOp ( DiagonalOp  op)

Update the GPU memory with the current values in op.real and op.imag.

This is required after making manual changes to op, but is not required after functions initDiagonalOp() and setDiagonalOpElems().

This function has no effect in other modes besides GPU mode.

See also
Parameters
[in,out]opthe DiagonalOp to synch to GPU
Exceptions
invalidQuESTInputError()
  • if op was not created
Author
Tyson Jones

Definition at line 1531 of file QuEST.c.

1531  {
1532  validateDiagOpInit(op, __func__);
1533 
1535 }

References agnostic_syncDiagonalOp(), and validateDiagOpInit().

Referenced by TEST_CASE().

void agnostic_destroyDiagonalOp(DiagonalOp op)
Definition: QuEST_cpu.c:1368
@ INVERSE_PRODUCT
Definition: QuEST.h:233
pauliOpType
Codes for specifying Pauli operators.
Definition: QuEST.h:96
void initPauliHamil(PauliHamil hamil, qreal *coeffs, enum pauliOpType *codes)
Initialise PauliHamil instance hamil with the given term coefficients and Pauli codes (one for every ...
Definition: QuEST.c:1504
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
@ DISTANCE
Definition: QuEST.h:234
int rank
Definition: QuEST.h:364
void validateHamilFileCoeffParsed(int parsed, PauliHamil h, FILE *file, char *fn, const char *caller)
void seedQuESTDefault(QuESTEnv *env)
Seeds the random number generator with the (master node) current time and process ID.
Definition: QuEST.c:1614
void validateHamilParams(int numQubits, int numTerms, const char *caller)
void qasm_free(Qureg qureg)
Definition: QuEST_qasm.c:887
@ PAULI_I
Definition: QuEST.h:96
ComplexMatrixN createComplexMatrixN(int numQubits)
Allocate dynamic memory for a square complex matrix of any size, which can be passed to functions lik...
Definition: QuEST.c:1348
void syncDiagonalOp(DiagonalOp op)
Update the GPU memory with the current values in op.real and op.imag.
Definition: QuEST.c:1531
int numSeeds
Definition: QuEST.h:367
void validateNumQubitsInQureg(int numQubits, int numRanks, const char *caller)
@ TWOS_COMPLEMENT
Definition: QuEST.h:269
void validateHamilFilePauliParsed(int parsed, PauliHamil h, FILE *file, char *fn, const char *caller)
void agnostic_syncDiagonalOp(DiagonalOp op)
Definition: QuEST_cpu.c:1373
void statevec_destroyQureg(Qureg qureg, QuESTEnv env)
Definition: QuEST_cpu.c:1328
void validateHamilFilePauliCode(enum pauliOpType code, PauliHamil h, FILE *file, char *fn, const char *caller)
@ NORM
Definition: QuEST.h:232
@ SCALED_INVERSE_DISTANCE
Definition: QuEST.h:234
Information about the environment the program is running in.
Definition: QuEST.h:362
@ UNSIGNED
Definition: QuEST.h:269
PauliHamil createPauliHamilFromFile(char *fn)
Creates a PauliHamil instance, a real-weighted sum of products of Pauli operators,...
Definition: QuEST.c:1420
void setDiagonalOpElems(DiagonalOp op, long long int startInd, qreal *real, qreal *imag, long long int numElems)
Modifies a subset (starting at index startInd, and ending at index startInd + numElems) of the elemen...
Definition: QuEST.c:1543
ComplexMatrixN bindArraysToStackComplexMatrixN(int numQubits, qreal re[][1<< numQubits], qreal im[][1<< numQubits], qreal **reStorage, qreal **imStorage)
Definition: QuEST_common.c:652
Represents a general 2^N by 2^N matrix of complex numbers.
Definition: QuEST.h:186
@ INVERSE_DISTANCE
Definition: QuEST.h:234
#define qreal
void validateMatrixInit(ComplexMatrixN matr, const char *caller)
void validateFileOpened(int opened, char *fn, const char *caller)
void multiQubitUnitary(Qureg qureg, int *targs, int numTargs, ComplexMatrixN u)
Apply a general multi-qubit unitary (including a global phase factor) with any number of target qubit...
Definition: QuEST.c:296
@ PAULI_X
Definition: QuEST.h:96
DiagonalOp agnostic_createDiagonalOp(int numQubits, QuESTEnv env)
Definition: QuEST_cpu.c:1346
int numQubitsInStateVec
Number of qubits in the state-vector - this is double the number represented for mixed states.
Definition: QuEST.h:329
int chunkId
The position of the chunk of the state vector held by this process in the full state vector.
Definition: QuEST.h:336
qreal * imag
The imaginary values of the 2^numQubits complex elements.
Definition: QuEST.h:310
long long int numAmpsPerChunk
Number of probability amplitudes held in stateVec by this process In the non-MPI version,...
Definition: QuEST.h:332
qreal * termCoeffs
The real coefficient of each Pauli product. This is an array of length PauliHamil....
Definition: QuEST.h:283
#define getStaticComplexMatrixN(numQubits, re, im)
Creates a ComplexMatrixN struct which lives in the stack and so does not need freeing,...
Definition: QuEST.h:5521
#define qcomp
enum pauliOpType * pauliCodes
The Pauli operators acting on each qubit, flattened over every operator.
Definition: QuEST.h:281
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:42
@ SCALED_PRODUCT
Definition: QuEST.h:233
void statevec_cloneQureg(Qureg targetQureg, Qureg copyQureg)
works for both statevectors and density matrices
Definition: QuEST_cpu.c:1572
int numRanks
Definition: QuEST.h:365
int numQubits
The number of qubits this operator can act on (informing its size)
Definition: QuEST.h:300
int numSumTerms
The number of terms in the weighted sum, or the number of Pauli products.
Definition: QuEST.h:285
Represents a diagonal complex operator on the full Hilbert state of a Qureg.
Definition: QuEST.h:297
@ PAULI_Y
Definition: QuEST.h:96
A Pauli Hamiltonian, expressed as a real-weighted sum of pauli products, and which can hence represen...
Definition: QuEST.h:277
void validateNumElems(DiagonalOp op, long long int startInd, long long int numElems, const char *caller)
@ SCALED_INVERSE_SHIFTED_NORM
Definition: QuEST.h:232
qreal ** real
Definition: QuEST.h:189
void validateDiagPauliHamil(DiagonalOp op, PauliHamil hamil, const char *caller)
unsigned long int * seeds
Definition: QuEST.h:366
@ INVERSE_NORM
Definition: QuEST.h:232
Represents a system of qubits.
Definition: QuEST.h:322
void validateDiagPauliHamilFromFile(PauliHamil hamil, int numRanks, const char *caller)
void validateNumQubitsInMatrix(int numQubits, const char *caller)
qreal ** imag
Definition: QuEST.h:190
void agnostic_initDiagonalOpFromPauliHamil(DiagonalOp op, PauliHamil hamil)
Definition: QuEST_cpu.c:1377
@ PRODUCT
Definition: QuEST.h:233
ComplexArray stateVec
Computational state amplitudes - a subset thereof in the MPI version.
Definition: QuEST.h:341
void validateDiagOpInit(DiagonalOp op, const char *caller)
long long int numElemsPerChunk
The number of the 2^numQubits amplitudes stored on each distributed node.
Definition: QuEST.h:302
void initDiagonalOpFromPauliHamil(DiagonalOp op, PauliHamil hamil)
Populates the diagonal operator op to be equivalent to the given Pauli Hamiltonian hamil,...
Definition: QuEST.c:1550
int isDensityMatrix
Whether this instance is a density-state representation.
Definition: QuEST.h:325
int numQubits
Definition: QuEST.h:188
void validateHamilFileParams(int numQubits, int numTerms, FILE *file, char *fn, const char *caller)
int numQubits
The number of qubits informing the Hilbert dimension of the Hamiltonian.
Definition: QuEST.h:287
int numQubitsRepresented
The number of qubits represented in either the state-vector or density matrix.
Definition: QuEST.h:327
long long int numAmpsTotal
Total number of amplitudes, which are possibly distributed among machines.
Definition: QuEST.h:334
@ SCALED_DISTANCE
Definition: QuEST.h:234
qreal * real
The real values of the 2^numQubits complex elements.
Definition: QuEST.h:308
Qureg createQureg(int numQubits, QuESTEnv env)
Creates a state-vector Qureg object representing a set of qubits which will remain in a pure state.
Definition: QuEST.c:36
void destroyPauliHamil(PauliHamil h)
Destroy a PauliHamil instance, created with either createPauliHamil() or createPauliHamilFromFile().
Definition: QuEST.c:1414
void agnostic_setDiagonalOpElems(DiagonalOp op, long long int startInd, qreal *real, qreal *imag, long long int numElems)
Definition: QuEST_cpu.c:4228
void validateNumQubitsInDiagOp(int numQubits, int numRanks, const char *caller)
void validateNumRanks(int numRanks, const char *caller)
@ SCALED_INVERSE_SHIFTED_DISTANCE
Definition: QuEST.h:234
@ SCALED_NORM
Definition: QuEST.h:232
@ SCALED_INVERSE_PRODUCT
Definition: QuEST.h:233
void initZeroState(Qureg qureg)
Initialise qureg into the zero state.
Definition: QuEST.c:113
void qasm_setup(Qureg *qureg)
Definition: QuEST_qasm.c:61
void validatePauliCodes(enum pauliOpType *pauliCodes, int numPauliCodes, const char *caller)
void statevec_createQureg(Qureg *qureg, int numQubits, QuESTEnv env)
Definition: QuEST_cpu.c:1290
PauliHamil createPauliHamil(int numQubits, int numSumTerms)
Dynamically allocates a Hamiltonian expressed as a real-weighted sum of products of Pauli operators.
Definition: QuEST.c:1398
DiagonalOp createDiagonalOp(int numQubits, QuESTEnv env)
Creates a DiagonalOp representing a diagonal operator on the full Hilbert space of a Qureg.
Definition: QuEST.c:1518
@ SCALED_INVERSE_NORM
Definition: QuEST.h:232