Debugging

Utilities for seeding and debugging, such as state-logging. More...

Functions

void copyStateFromGPU (Qureg qureg)
 In GPU mode, this copies the state-vector (or density matrix) from GPU memory (qureg.deviceStateVec) to RAM (qureg.stateVec), where it can be accessed/modified by the user. More...
 
void copyStateToGPU (Qureg qureg)
 In GPU mode, this copies the state-vector (or density matrix) from RAM (qureg.stateVec) to VRAM / GPU-memory (qureg.deviceStateVec), which is the version operated upon by other calls to the API. More...
 
void getEnvironmentString (QuESTEnv env, char str[200])
 Sets str to a string containing information about the runtime environment, including whether simulation is using CUDA (for GPU), OpenMP (for multithreading) and/or MPI (for distribution). More...
 
void getQuESTSeeds (QuESTEnv env, unsigned long int **seeds, int *numSeeds)
 Obtain the seeds presently used in random number generation. More...
 
void initDebugState (Qureg qureg)
 Initialises qureg to be in the un-normalised, non-physical state with with $n$-th complex amplitude given by $2n/10 + i(2n+1)/10$. More...
 
void invalidQuESTInputError (const char *errMsg, const char *errFunc)
 An internal function called when invalid arguments are passed to a QuEST API call, which the user can optionally override by redefining. More...
 
void reportPauliHamil (PauliHamil hamil)
 Print the PauliHamil to screen. More...
 
void reportQuESTEnv (QuESTEnv env)
 Report information about the QuEST environment. More...
 
void reportQuregParams (Qureg qureg)
 Report metainformation about a set of qubits: number of qubits, number of probability amplitudes. More...
 
void reportState (Qureg qureg)
 Print the current state vector of probability amplitudes for a set of qubits to file. More...
 
void reportStateToScreen (Qureg qureg, QuESTEnv env, int reportRank)
 Print the current state vector of probability amplitudes for a set of qubits to standard out. More...
 
void seedQuEST (QuESTEnv *env, unsigned long int *seedArray, int numSeeds)
 Seeds the random number generator with a custom array of key(s), overriding the default keys. More...
 
void seedQuESTDefault (QuESTEnv *env)
 Seeds the random number generator with the (master node) current time and process ID. More...
 
void syncQuESTEnv (QuESTEnv env)
 Guarantees that all code up to the given point has been executed on all nodes (if running in distributed mode) More...
 
int syncQuESTSuccess (int successCode)
 Performs a logical AND on all successCodes held by all processes. More...
 

Detailed Description

Utilities for seeding and debugging, such as state-logging.

Author
Ania Brown
Tyson Jones
Balint Koczor
Nicolas Vogt of HQS (one-qubit damping)

Function Documentation

◆ copyStateFromGPU()

void copyStateFromGPU ( Qureg  qureg)

In GPU mode, this copies the state-vector (or density matrix) from GPU memory (qureg.deviceStateVec) to RAM (qureg.stateVec), where it can be accessed/modified by the user.

In CPU mode, this function has no effect. In conjunction with copyStateToGPU(), this allows a user to directly modify the state-vector in a harware agnostic way. Note though that users should instead use setAmps() if possible.

For example, to set the first real element to 1, one could do:

qureg.stateVec.real[0] = 1;

Note users should never access qureg.deviceStateVec directly.

See also
Parameters
[in,out]quregthe qureg of which to copy .deviceStateVec to .stateVec in GPU mode
Author
Ania Brown
Tyson Jones (doc)

Definition at line 45 of file QuEST_cpu.c.

45  {
46 }

References DEBUG, Qureg::deviceStateVec, Qureg::numAmpsPerChunk, and Qureg::stateVec.

Referenced by areEqual(), densmatr_calcTotalProb(), statevec_calcTotalProb(), statevec_compareStates(), statevec_reportStateToScreen(), TEST_CASE(), toQMatrix(), and toQVector().

◆ copyStateToGPU()

void copyStateToGPU ( Qureg  qureg)

In GPU mode, this copies the state-vector (or density matrix) from RAM (qureg.stateVec) to VRAM / GPU-memory (qureg.deviceStateVec), which is the version operated upon by other calls to the API.

In CPU mode, this function has no effect. In conjunction with copyStateFromGPU() (which should be called first), this allows a user to directly modify the state-vector in a harware agnostic way. Note though that users should instead use setAmps() if possible.

For example, to set the first real element to 1, one could do:

qureg.stateVec.real[0] = 1;

Note users should never access qureg.deviceStateVec directly.

See also
Parameters
[in,out]quregthe qureg of which to copy .stateVec to .deviceStateVec in GPU mode
Author
Ania Brown
Tyson Jones (doc)

Definition at line 42 of file QuEST_cpu.c.

42  {
43 }

References DEBUG, Qureg::deviceStateVec, Qureg::numAmpsPerChunk, and Qureg::stateVec.

Referenced by statevec_initStateFromSingleFile(), TEST_CASE(), and toQureg().

◆ getEnvironmentString()

void getEnvironmentString ( QuESTEnv  env,
char  str[200] 
)

Sets str to a string containing information about the runtime environment, including whether simulation is using CUDA (for GPU), OpenMP (for multithreading) and/or MPI (for distribution).

The number of CPU threads and distributed ranks is also reported. Note there is currently no reporting of the number of GPU cores used.

The string format is:

"CUDA=b OpenMP=b MPI=b threads=n ranks=n"

where b is 0 or 1, and n are integers.

Parameters
[in]envobject representing the execution environment. A single instance is used for each program
[out]strto be populated with the output string
Author
Ania Brown
Tyson Jones

Definition at line 198 of file QuEST_cpu_distributed.c.

198  {
199  int ompStatus=0;
200  int numThreads=1;
201 # ifdef _OPENMP
202  ompStatus=1;
203  numThreads=omp_get_max_threads();
204 # endif
205  sprintf(str, "CUDA=0 OpenMP=%d MPI=1 threads=%d ranks=%d", ompStatus, numThreads, env.numRanks);
206 }

References QuESTEnv::numRanks.

◆ getQuESTSeeds()

void getQuESTSeeds ( QuESTEnv  env,
unsigned long int **  seeds,
int *  numSeeds 
)

Obtain the seeds presently used in random number generation.

This function sets argument seeds to the address of the array of keys which have seeded QuEST's Mersenne Twister random number generator. numSeeds is set to the length of seeds. These are the seeds which inform the outcomes of random functions like measure(), and are set using seedQuEST() and seedQuESTDefault().

The output seeds array must not be freed, and should not be modified.

Obtaining QuEST's seeds is useful for seeding your own random number generators, so that a simulation (with random QuEST measurements, and your own random decisions) can be precisely repeated later, just by calling seedQuEST().

Note this function merely sets the arguments to the attributes for env. I.e.

unsigned long int* seeds;
int numSeeds;
getQuESTSeeds(env, &seeds, &numSeeds);
func(seeds, numSeeds);

is equivalent to

func(env.seeds, env.numSeeds);

However, one should not rely upon their local pointer from getQuESTSeeds() to be automatically updated after a subsequent call to seedQuEST() or seedQuESTDefault(). Instead, getQuESTSeeds() should be recalled.

See also
Parameters
[in]envthe QuESTEnv runtime environment
[in]seedsa pointer to an unitialised array to be modified
[in]numSeedsa pointer to an integer to be modified
Author
Tyson Jones

Definition at line 1622 of file QuEST.c.

1622  {
1623  *seeds = env.seeds;
1624  *numSeeds = env.numSeeds;
1625 }

References QuESTEnv::numSeeds, and QuESTEnv::seeds.

Referenced by main().

◆ initDebugState()

void initDebugState ( Qureg  qureg)

Initialises qureg to be in the un-normalised, non-physical state with with $n$-th complex amplitude given by $2n/10 + i(2n+1)/10$.

This is used internally for debugging and testing.

Parameters
[in,out]quregthe register to have its amplitudes overwritten
Author
Ania Brown
Tyson Jones (doc)

Definition at line 1578 of file QuEST.c.

1578  {
1579  statevec_initDebugState(qureg);
1580 }

References statevec_initDebugState().

Referenced by TEST_CASE().

◆ invalidQuESTInputError()

void invalidQuESTInputError ( const char *  errMsg,
const char *  errFunc 
)

An internal function called when invalid arguments are passed to a QuEST API call, which the user can optionally override by redefining.

This function is a weak symbol, so that users can choose how input errors are handled, by redefining it in their own code. Users must ensure that the triggered API call does not continue (e.g. the user exits or throws an exception), else QuEST will continue with the valid input and likely trigger a seg-fault. This function is triggered before any internal state-change, hence it is safe to interrupt with exceptions.

E.g. in C

void invalidQuESTInputError(const char* errMsg, const char* errFunc) {
// log to file
printf("ERROR! Writing to file...\n");
FILE *fp = fopen("errorlog.txt", "w");
fprintf(fp, "incorrect usage of function '%s': %s", errFunc, errMsg);
fclose(fp);
// exit
exit(1);
}

This function is compatible with C++ exceptions, though note the user of extern "C":

extern "C" void invalidQuESTInputError(const char* errMsg, const char* errFunc) {
string err = "in function " + string(errFunc) + ": " + string(errMsg);
throw std::invalid_argument(err);
}
Parameters
[in]errMsga string describing the nature of the argument error
[in]errFuncthe name of the invalidly-called API function
Exceptions
invalidQuESTInputError()
  • unless overriden by the user
Author
Tyson Jones

An internal function called when invalid arguments are passed to a QuEST API call, which the user can optionally override by redefining.

an negative qubit index). This is redefined here to, in lieu of printing and exiting, throw a C++ exception which can be caught (and hence unit tested for) by Catch2

Definition at line 209 of file QuEST_validation.c.

209  {
210  default_invalidQuESTInputError(errMsg, errFunc);
211 }

References default_invalidQuESTInputError().

Referenced by QuESTAssert(), validateFileOpened(), validateHamilFileCoeffParsed(), validateHamilFileParams(), validateHamilFilePauliCode(), and validateHamilFilePauliParsed().

◆ reportPauliHamil()

void reportPauliHamil ( PauliHamil  hamil)

Print the PauliHamil to screen.

The output features a new line for each term, each with format

c p1 p2 p3 ... pN

where c is the real coefficient of the term, and p1 ... pN are numbers 0, 1, 2, 3 to indicate identity, pauliX, pauliY and pauliZ operators respectively, acting on qubits 0 through N-1 (all qubits). A tab character separates c and p1, but single spaces separate the Pauli operators.

See also
Parameters
[in]hamilan instantiated PauliHamil
Exceptions
invalidQuESTInputErrorif the parameters of hamil are invalid, i.e. if numQubits <= 0, or if numSumTerms <= 0, or if pauliCodes contains an invalid Pauli code.
Author
Tyson Jones

Definition at line 1598 of file QuEST.c.

1598  {
1599  validatePauliHamil(hamil, __func__);
1600 
1601  for (int t=0; t<hamil.numSumTerms; t++) {
1602  printf(REAL_QASM_FORMAT, hamil.termCoeffs[t]);
1603  printf("\t");
1604  for (int q=0; q<hamil.numQubits; q++)
1605  printf("%d ", (int) hamil.pauliCodes[q+t*hamil.numQubits]);
1606  printf("\n");
1607  }
1608 }

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

◆ reportQuESTEnv()

void reportQuESTEnv ( QuESTEnv  env)

Report information about the QuEST environment.

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

Definition at line 183 of file QuEST_cpu_distributed.c.

183  {
184  if (env.rank==0){
185  printf("EXECUTION ENVIRONMENT:\n");
186  printf("Running distributed (MPI) version\n");
187  printf("Number of ranks is %d\n", env.numRanks);
188 # ifdef _OPENMP
189  printf("OpenMP enabled\n");
190  printf("Number of threads available is %d\n", omp_get_max_threads());
191 # else
192  printf("OpenMP disabled\n");
193 # endif
194  printf("Precision: size of qreal is %ld bytes\n", sizeof(qreal) );
195  }
196 }

References QuESTEnv::numRanks, qreal, and QuESTEnv::rank.

◆ reportQuregParams()

void reportQuregParams ( Qureg  qureg)

Report metainformation about a set of qubits: number of qubits, number of probability amplitudes.

Parameters
[in]quregobject representing the set of qubits
Author
Ania Brown

Definition at line 237 of file QuEST_common.c.

237  {
238  long long int numAmps = 1LL << qureg.numQubitsInStateVec;
239  long long int numAmpsPerRank = numAmps/qureg.numChunks;
240  if (qureg.chunkId==0){
241  printf("QUBITS:\n");
242  printf("Number of qubits is %d.\n", qureg.numQubitsInStateVec);
243  printf("Number of amps is %lld.\n", numAmps);
244  printf("Number of amps per rank is %lld.\n", numAmpsPerRank);
245  }
246 }

References Qureg::chunkId, Qureg::numChunks, and Qureg::numQubitsInStateVec.

◆ reportState()

void reportState ( Qureg  qureg)

Print the current state vector of probability amplitudes for a set of qubits to file.

File format:

real, imag
realComponent1, imagComponent1
realComponent2, imagComponent2
...
realComponentN, imagComponentN

File naming convention:

For each node that the program runs on, a file 'state_rank_[node_rank].csv' is generated. If there is more than one node, ranks after the first do not include the header

real, imag

so that files are easier to combine.

Parameters
[in,out]quregobject representing the set of qubits
Author
Ania Brown

Definition at line 219 of file QuEST_common.c.

219  {
220  FILE *state;
221  char filename[100];
222  long long int index;
223  sprintf(filename, "state_rank_%d.csv", qureg.chunkId);
224  state = fopen(filename, "w");
225  if (qureg.chunkId==0) fprintf(state, "real, imag\n");
226 
227  for(index=0; index<qureg.numAmpsPerChunk; index++){
228  # if QuEST_PREC==1 || QuEST_PREC==2
229  fprintf(state, "%.12f, %.12f\n", qureg.stateVec.real[index], qureg.stateVec.imag[index]);
230  # elif QuEST_PREC == 4
231  fprintf(state, "%.12Lf, %.12Lf\n", qureg.stateVec.real[index], qureg.stateVec.imag[index]);
232  #endif
233  }
234  fclose(state);
235 }

References Qureg::chunkId, Qureg::numAmpsPerChunk, and Qureg::stateVec.

◆ reportStateToScreen()

void reportStateToScreen ( Qureg  qureg,
QuESTEnv  env,
int  reportRank 
)

Print the current state vector of probability amplitudes for a set of qubits to standard out.

For debugging purposes. Each rank should print output serially. Only print output for systems <= 5 qubits

Author
Ania Brown

Definition at line 1594 of file QuEST.c.

1594  {
1595  statevec_reportStateToScreen(qureg, env, reportRank);
1596 }

References statevec_reportStateToScreen().

◆ seedQuEST()

void seedQuEST ( QuESTEnv env,
unsigned long int *  seedArray,
int  numSeeds 
)

Seeds the random number generator with a custom array of key(s), overriding the default keys.

This determines the sequence of outcomes in functions like measure() and measureWithStats().

In distributed mode, the key(s) passed to the master node will be broadcast to all other nodes, such that every node generates the same sequence of pseudorandom numbers.

This function will copy the contents of seedArray into a permanent array env.seeds, so seedArray is afterward safe to free.

QuEST uses the Mersenne Twister for random number generation.

See also
Parameters
[in]enva pointer to the QuESTEnv runtime environment
[in]seedArrayArray of integers to use as seed. This allows the MT to be initialised with more than a 32-bit integer if required
[in]numSeedsLength of seedArray
Author
Ania Brown
Tyson Jones (doc)

Definition at line 1388 of file QuEST_cpu_distributed.c.

1388  {
1389 
1390  // it is imperative every node agrees on the seed, so that random decisions
1391  // agree on every node. Hence we use only the master node keys.
1392  MPI_Bcast(seedArray, numSeeds, MPI_UNSIGNED_LONG, 0, MPI_COMM_WORLD);
1393 
1394  // free existing seed array, if exists
1395  if (env->seeds != NULL)
1396  free(env->seeds);
1397 
1398  // record keys in permanent heap
1399  env->seeds = malloc(numSeeds * sizeof *(env->seeds));
1400  for (int i=0; i<numSeeds; i++)
1401  (env->seeds)[i] = seedArray[i];
1402  env->numSeeds = numSeeds;
1403 
1404  // pass keys to Mersenne Twister seeder
1405  init_by_array(seedArray, numSeeds);
1406 }

References init_by_array(), QuESTEnv::numSeeds, and QuESTEnv::seeds.

Referenced by main(), and seedQuESTDefault().

◆ seedQuESTDefault()

void seedQuESTDefault ( QuESTEnv env)

Seeds the random number generator with the (master node) current time and process ID.

This is the default seeding used by createQuESTEnv(), and determines the outcomes in functions like measure() and measureWithStats().

In distributed mode, every node agrees on the seed (nominated by the master node) such that every node generates the same sequence of pseudorandom numbers.

QuEST uses the Mersenne Twister for random number generation.

See also
  • Use seedQuEST() to provide a custom seed, overriding the default.
  • Use getQuESTSeeds() to obtain the seeds currently being used for RNG.
Parameters
[in]enva pointer to the QuESTEnv runtime environment
Author
Ania Brown
Balint Koczor (Windows compatibility)
Tyson Jones (doc)

Definition at line 1614 of file QuEST.c.

1614  {
1615 
1616  // seed Mersenne Twister random number generator with two keys -- time and pid
1617  unsigned long int keys[2];
1618  getQuESTDefaultSeedKey(keys);
1619  seedQuEST(env, keys, 2);
1620 }

References getQuESTDefaultSeedKey(), and seedQuEST().

Referenced by createQuESTEnv(), and main().

◆ syncQuESTEnv()

void syncQuESTEnv ( QuESTEnv  env)

Guarantees that all code up to the given point has been executed on all nodes (if running in distributed mode)

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

Definition at line 164 of file QuEST_cpu_distributed.c.

164  {
165  MPI_Barrier(MPI_COMM_WORLD);
166 }

Referenced by areEqual(), deleteFilesWithPrefixSynch(), statevec_initStateFromSingleFile(), statevec_reportStateToScreen(), TEST_CASE(), toQMatrix(), toQureg(), toQVector(), and writeToFileSynch().

◆ syncQuESTSuccess()

int syncQuESTSuccess ( int  successCode)

Performs a logical AND on all successCodes held by all processes.

If any one process has a zero successCode all processes will return a zero success code.

Parameters
[in]successCode1 if process task succeeded, 0 if process task failed
Returns
1 if all processes succeeded, 0 if any one process failed
Author
Ania Brown

Definition at line 168 of file QuEST_cpu_distributed.c.

168  {
169  int totalSuccess;
170  MPI_Allreduce(&successCode, &totalSuccess, 1, MPI_INT, MPI_LAND, MPI_COMM_WORLD);
171  return totalSuccess;
172 }
void init_by_array(unsigned long init_key[], int key_length)
Definition: mt19937ar.c:80
void copyStateFromGPU(Qureg qureg)
In GPU mode, this copies the state-vector (or density matrix) from GPU memory (qureg....
Definition: QuEST_cpu.c:45
int rank
Definition: QuEST.h:364
int numSeeds
Definition: QuEST.h:367
int numChunks
Number of chunks the state vector is broken up into – the number of MPI processes used.
Definition: QuEST.h:338
void getQuESTDefaultSeedKey(unsigned long int *key)
Definition: QuEST_common.c:195
void getQuESTSeeds(QuESTEnv env, unsigned long int **seeds, int *numSeeds)
Obtain the seeds presently used in random number generation.
Definition: QuEST.c:1622
#define qreal
int numQubitsInStateVec
Number of qubits in the state-vector - this is double the number represented for mixed states.
Definition: QuEST.h:329
void invalidQuESTInputError(const char *errMsg, const char *errFunc)
An internal function called when invalid arguments are passed to a QuEST API call,...
int chunkId
The position of the chunk of the state vector held by this process in the full state vector.
Definition: QuEST.h:336
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
enum pauliOpType * pauliCodes
The Pauli operators acting on each qubit, flattened over every operator.
Definition: QuEST.h:281
void statevec_reportStateToScreen(Qureg qureg, QuESTEnv env, int reportRank)
Print the current state vector of probability amplitudes for a set of qubits to standard out.
Definition: QuEST_cpu.c:1439
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
int numRanks
Definition: QuEST.h:365
int numSumTerms
The number of terms in the weighted sum, or the number of Pauli products.
Definition: QuEST.h:285
unsigned long int * seeds
Definition: QuEST.h:366
void statevec_initDebugState(Qureg qureg)
Initialise the state vector of probability amplitudes to an (unphysical) state with each component of...
Definition: QuEST_cpu.c:1657
void validatePauliHamil(PauliHamil hamil, const char *caller)
ComplexArray stateVec
Computational state amplitudes - a subset thereof in the MPI version.
Definition: QuEST.h:341
int numQubits
The number of qubits informing the Hilbert dimension of the Hamiltonian.
Definition: QuEST.h:287
void seedQuEST(QuESTEnv *env, unsigned long int *seedArray, int numSeeds)
Seeds the random number generator with a custom array of key(s), overriding the default keys.
void default_invalidQuESTInputError(const char *errMsg, const char *errFunc)