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 copySubstateFromGPU (Qureg qureg, long long int startInd, long long int numAmps)
 In GPU mode, this copies a substate of the state-vector (or density matrix) from GPU VRAM (qureg.deviceStateVec) into RAM (qureg.stateVec). More...
 
void copySubstateToGPU (Qureg qureg, long long int startInd, long long int numAmps)
 In GPU mode, this copies a substate of the state-vector (or density matrix) from RAM (qureg.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;
void copyStateToGPU(Qureg qureg)
In GPU mode, this copies the state-vector (or density matrix) from RAM (qureg.stateVec) to VRAM / GPU...
Definition: QuEST_cpu.c:43
void copyStateFromGPU(Qureg qureg)
In GPU mode, this copies the state-vector (or density matrix) from GPU memory (qureg....
Definition: QuEST_cpu.c:46

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)

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

Referenced by TEST_CASE(), and toQureg().

◆ copySubstateFromGPU()

void copySubstateFromGPU ( Qureg  qureg,
long long int  startInd,
long long int  numAmps 
)

In GPU mode, this copies a substate of the state-vector (or density matrix) from GPU VRAM (qureg.deviceStateVec) into RAM (qureg.stateVec).

In CPU mode, this function has no effect. In conjunction with copySubstateToGPU(), this allows a user to directly modify a subset of the amplitudes the state-vector in a hardware agnostic way.

For example, to multiply the first amplitude by factor 2, one could do

copySubstateFromGPU(qureg, 0, 1);
qureg.stateVec.real[0] *= 2;
qureg.stateVec.imag[0] *= 2;
copySubstateToGPU(qureg, 0, 1);
void copySubstateToGPU(Qureg qureg, long long int startInd, long long int numAmps)
In GPU mode, this copies a substate of the state-vector (or density matrix) from RAM (qureg....
Definition: QuEST.c:1755
void copySubstateFromGPU(Qureg qureg, long long int startInd, long long int numAmps)
In GPU mode, this copies a substate of the state-vector (or density matrix) from GPU VRAM (qureg....
Definition: QuEST.c:1761

Note users should never access qureg.deviceStateVec directly.

See also
Parameters
[in,out]quregthe qureg of which to copy .deviceStateVec to .stateVec to in GPU mode
[in]startIndthe index of the first amplitude to copy
[in]numAmpsthe number of contiguous amplitudes to copy (starting with startInd)
Exceptions
invalidQuESTInputError()
  • if startInd is an invalid amplitude index
  • if numAmps is greater than the remaining amplitudes in the state, from startInd
Author
Tyson Jones

◆ copySubstateToGPU()

void copySubstateToGPU ( Qureg  qureg,
long long int  startInd,
long long int  numAmps 
)

In GPU mode, this copies a substate of the state-vector (or density matrix) from RAM (qureg.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 copySubstateFromGPU(), this allows a user to directly modify a subset of the amplitudes the state-vector in a harware agnostic way, without having to load the entire state via copyStateFromGPU().

Note though that users should instead use setAmps() if possible.

For example, to multiply the first amplitude by factor 2, one could do

copySubstateFromGPU(qureg, 0, 1);
qureg.stateVec.real[0] *= 2;
qureg.stateVec.imag[0] *= 2;
copySubstateToGPU(qureg, 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
[in]startIndthe index of the first amplitude to copy
[in]numAmpsthe number of contiguous amplitudes to copy (starting with startInd)
Exceptions
invalidQuESTInputError()
  • if startInd is an invalid amplitude index
  • if numAmps is greater than the remaining amplitudes in the state, from startInd
Author
Tyson Jones

◆ 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

◆ 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);
void getQuESTSeeds(QuESTEnv env, unsigned long int **seeds, int *numSeeds)
Obtain the seeds presently used in random number generation.
Definition: QuEST.c:1750

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

◆ 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)

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);
}
void invalidQuESTInputError(const char *errMsg, const char *errFunc)
An internal function called when invalid arguments are passed to a QuEST API call,...
Definition: QuEST_validation.c:231

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

◆ 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

◆ 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

◆ 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

◆ 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

◆ 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

◆ 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)

Referenced by 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)

Referenced by createQuESTEnv().

◆ 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

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