Examples#
In this section, we show a basic example on how to use the cuPauliProp API for an end-to-end quantum circuit simulation.
Building code#
Assuming your cuPauliprop header files are located in CUPAULIPROP_INCLUDE_DIR and the shared library is located in CUPAULIPROP_LIB_DIR, the example code can be built via the following command:
nvcc kicked_ising_example.cpp -I${CUPAULIPROP_INCLUDE_DIR} -L${CUPAULIPROP_LIB_DIR} -lcupauliprop -o kicked_ising_example
When building statically, use the following command:
nvcc kicked_ising_example.cpp -I${CUPAULIPROP_INCLUDE_DIR} ${CUPAULIPROP_LIB_DIR}/libcupauliprop_static.a -o kicked_ising_example
When using cuQuantum downloaded from the NVIDIA devzone website and extracted in CUQUANTUM_ROOT, the header files are located in ${CUQUANTUM_ROOT}/include and the shared library is located in ${CUQUANTUM_ROOT}/lib.
For convenience, we also provide a Makefile and CMakeLists.txt file in the same folder as the example code.
Kicked Ising example#
In this example, we demonstrate use of the cuPauliProp API to perform an end-to-end simulation of IBM’s
127 qubit kicked Ising experiment, as presented in Nature volume 618, pages 500–505 (2023). We will
specifically simulate a single datum of Figure 4 b, which is the expectation value of observable
\(Z_{62}\) after a 20-repetition Trotter circuit. The full demonstration is available at
kicked_ising_example.cpp.
Our example code will make use of the below helper macros which handle errors thrown by the CUDA or cuPauliProp APIs.
27// ========================================================================
28// CUDA and cuPauliProp error handling
29// ========================================================================
30
31#define HANDLE_CUPP_ERROR(x) \
32{ \
33 const auto err = x; \
34 if (err != CUPAULIPROP_STATUS_SUCCESS) \
35 { \
36 printf("cuPauliProp error in line %d\n", __LINE__); \
37 fflush(stdout); \
38 std::abort(); \
39 } \
40};
41
42
43#define HANDLE_CUDA_ERROR(x) \
44{ \
45 const auto err = x; \
46 if (err != cudaSuccess) \
47 { \
48 const char * error = cudaGetErrorString(err); \
49 printf("CUDA Error: %s in line %d\n", error, __LINE__); \
50 fflush(stdout); \
51 std::abort(); \
52 } \
53};
We next decide the memory budget for our simulation. It will turn out that for our chosen combination of
simulated circuit, output observable and truncation parameters, only around 60 MiB of total device memory
is needed. We will demonstrate how to dedicate all available device memory for simulation however, which
this example optionally performs when USE_MAX_VRAM = true.
57// ========================================================================
58// Memory usage
59// ========================================================================
60
61// Each Pauli expansion has two pre-allocated GPU buffers, storing packed
62// integers (which encode Pauli strings) and corresponding coefficients.
63// As much memory can be dedicated as your hardware allows, while the min-
64// imum required is specific and very sensitive to the simulated circuit,
65// studied observable, and the chosen truncation hyperparameters.
66// Some operations also require additional workspace memory which is also
67// ideally pre-allocated, and can be established using the API 'Prepare'
68// functions (e.g. cupaulipropPauliExpansionViewPrepareTraceWithZeroState).
69// In this demo, we dedicate either (a percentage of) the entirety of GPU
70// memory uniformly between the required memory buffers, or instead use a
71// fixed hardcoded amount which has been prior tested to be consistent with
72// our other simulation parameters (like truncations); these choices are
73// toggled via USE_MAX_VRAM below.
74
75// =true to use MAX_VRAM_PERCENT of VRAM, and =false to use fixed memories below
76bool USE_MAX_VRAM = false;
77double MAX_VRAM_PERCENT = 90; // 0-100%
78
79size_t FIXED_EXPANSION_PAULI_MEM = 16 * (1LLU << 20); // bytes = 16 MiB
80size_t FIXED_EXPANSION_COEF_MEM = 4 * (1LLU << 20); // bytes = 4 MiB
81size_t FIXED_WORKSPACE_MEM = 20 * (1LLU << 20); // bytes = 20 MiB
We here scaffold the quantum circuit we will subsequently instantiate. In short, it is a Trotter circuit encoding the dynamics of a kicked Ising system, with a heavy-hex topology matching the connectivity of IBM’s 127 qubit Eagle device.
85// ========================================================================
86// Circuit preparation (Trotterised Ising on IBM heavy hex topology)
87// ========================================================================
88
89// This demo simulates the circuits experimentally executed by IBM in article
90// 'Nature volume 618, pages 500–505 (2023)'. This is a circuit Trotterising
91// the evolution operator of a 2D transverse-field Ising model, but where the
92// prescribed ZZ rotations have a fixed angle of -pi/2, and where the X angles
93// are arbitrarily set/swept; later, we will fix the X angle to be pi/4. The
94// Hamiltonian ZZ interactions are confined to a heavy-hex topology, matching
95// the connectivity of the IBM Eagle processor 'ibm_kyiv', as we fix below.
96
97const int NUM_CIRCUIT_QUBITS = 127;
98const int NUM_ROTATIONS_PER_LAYER = 48;
99const int NUM_PAULIS_PER_X_ROTATION = 1;
100const int NUM_PAULIS_PER_Z_ROTATION = 2;
101
102const double PI = 3.14159265358979323846;
103const double ZZ_ROTATION_ANGLE = - PI / 2.0;
104
105// Indices of ZZ-interacting qubits which undergo the first (red) Trotter round
106const int32_t ZZ_QUBITS_RED[NUM_ROTATIONS_PER_LAYER][NUM_PAULIS_PER_Z_ROTATION] = {
107 { 2, 1}, { 33, 39}, { 59, 60}, { 66, 67}, { 72, 81}, {118, 119},
108 { 21, 20}, { 26, 25}, { 13, 12}, { 31, 32}, { 70, 74}, {122, 123},
109 { 96, 97}, { 57, 56}, { 63, 64}, {107, 108}, {103, 104}, { 46, 45},
110 { 28, 35}, { 7, 6}, { 79, 78}, { 5, 4}, {109, 114}, { 62, 61},
111 { 58, 71}, { 37, 52}, { 76, 77}, { 0, 14}, { 36, 51}, {106, 105},
112 { 73, 85}, { 88, 87}, { 68, 55}, {116, 115}, { 94, 95}, {100, 110},
113 { 17, 30}, { 92, 102}, { 50, 49}, { 83, 84}, { 48, 47}, { 98, 99},
114 { 8, 9}, {121, 120}, { 23, 24}, { 44, 43}, { 22, 15}, { 53, 41}
115};
116
117// Indices of ZZ-interacting qubits which undergo the second (blue) Trotter round
118const int32_t ZZ_QUBITS_BLUE[NUM_ROTATIONS_PER_LAYER][NUM_PAULIS_PER_Z_ROTATION] = {
119 { 53, 60}, {123, 124}, { 21, 22}, { 11, 12}, { 67, 68}, { 2, 3},
120 { 66, 65}, {122, 121}, {110, 118}, { 6, 5}, { 94, 90}, { 28, 29},
121 { 14, 18}, { 63, 62}, {111, 104}, {100, 99}, { 45, 44}, { 4, 15},
122 { 20, 19}, { 57, 58}, { 77, 71}, { 76, 75}, { 26, 27}, { 16, 8},
123 { 35, 47}, { 31, 30}, { 48, 49}, { 69, 70}, {125, 126}, { 89, 74},
124 { 80, 79}, {116, 117}, {114, 113}, { 10, 9}, {106, 93}, {101, 102},
125 { 92, 83}, { 98, 91}, { 82, 81}, { 54, 64}, { 96, 109}, { 85, 84},
126 { 87, 86}, {108, 112}, { 34, 24}, { 42, 43}, { 40, 41}, { 39, 38}
127};
128
129// Indices of ZZ-interacting qubits which undergo the third (green) Trotter round
130const int32_t ZZ_QUBITS_GREEN[NUM_ROTATIONS_PER_LAYER][NUM_PAULIS_PER_Z_ROTATION] = {
131 { 10, 11}, { 54, 45}, {111, 122}, { 64, 65}, { 60, 61}, {103, 102},
132 { 72, 62}, { 4, 3}, { 33, 20}, { 58, 59}, { 26, 16}, { 28, 27},
133 { 8, 7}, {104, 105}, { 73, 66}, { 87, 93}, { 85, 86}, { 55, 49},
134 { 68, 69}, { 89, 88}, { 80, 81}, {117, 118}, {101, 100}, {114, 115},
135 { 96, 95}, { 29, 30}, {106, 107}, { 83, 82}, { 91, 79}, { 0, 1},
136 { 56, 52}, { 90, 75}, {126, 112}, { 36, 32}, { 46, 47}, { 77, 78},
137 { 97, 98}, { 17, 12}, {119, 120}, { 22, 23}, { 24, 25}, { 43, 34},
138 { 42, 41}, { 40, 39}, { 37, 38}, {125, 124}, { 50, 51}, { 18, 19}
139};
We instantiate the circuit as a std::vector of cupaulipropQuantumOperator_t types. In this demonstration,
all circuit gates are Pauli rotations of 1 or 2 qubits, chiefly
We fix \(\theta_Z\) to be \(-\pi/2\), and will later substitute a value for the X-rotation strength \(\theta_X\).
143// ========================================================================
144// Circuit construction
145// ========================================================================
146
147// Each 'step' of the Trotter circuit alternates a layer of single-qubit X
148// rotations on every qubit, then a sequence of two-qubit Z rotations on the
149// heavy-hex topology, upon qubit pairs in the red, blue and green lists
150// above. Note that ZZ rotations about -pi/2 are actually Clifford, though
151// we still here treat them like a generic Pauli rotation. The functions
152// below construct a Trotter circuit with a variable number of steps.
153
154
155std::vector<cupaulipropQuantumOperator_t> getXRotationLayer(
156 cupaulipropHandle_t handle, double xRotationAngle
157) {
158 std::vector<cupaulipropQuantumOperator_t> layer(NUM_CIRCUIT_QUBITS);
159
160 const cupaulipropPauliKind_t paulis[NUM_PAULIS_PER_X_ROTATION] = {CUPAULIPROP_PAULI_X};
161
162 for (int32_t i=0; i<NUM_CIRCUIT_QUBITS; i++) {
163 HANDLE_CUPP_ERROR(cupaulipropCreatePauliRotationGateOperator(
164 handle, xRotationAngle, NUM_PAULIS_PER_X_ROTATION, &i, paulis, &layer[i]));
165 }
166
167 return layer;
168}
169
170
171std::vector<cupaulipropQuantumOperator_t> getZZRotationLayer(
172 cupaulipropHandle_t handle,
173 const int32_t topology[NUM_ROTATIONS_PER_LAYER][NUM_PAULIS_PER_Z_ROTATION]
174) {
175 std::vector<cupaulipropQuantumOperator_t> layer(NUM_ROTATIONS_PER_LAYER);
176
177 const cupaulipropPauliKind_t paulis[NUM_PAULIS_PER_Z_ROTATION] = {
178 CUPAULIPROP_PAULI_Z, CUPAULIPROP_PAULI_Z};
179
180 for (uint32_t i=0; i<NUM_ROTATIONS_PER_LAYER; i++) {
181 HANDLE_CUPP_ERROR(cupaulipropCreatePauliRotationGateOperator(
182 handle, ZZ_ROTATION_ANGLE, NUM_PAULIS_PER_Z_ROTATION, topology[i], paulis, &layer[i]));
183 }
184
185 return layer;
186}
187
188
189std::vector<cupaulipropQuantumOperator_t> getIBMHeavyHexIsingCircuit(
190 cupaulipropHandle_t handle, double xRotationAngle, int numTrotterSteps
191) {
192 std::vector<cupaulipropQuantumOperator_t> circuit;
193
194 for (int n=0; n<numTrotterSteps; n++) {
195 auto layerX = getXRotationLayer (handle, xRotationAngle);
196 auto layerRedZZ = getZZRotationLayer(handle, ZZ_QUBITS_RED);
197 auto layerBlueZZ = getZZRotationLayer(handle, ZZ_QUBITS_BLUE);
198 auto layerGreenZZ = getZZRotationLayer(handle, ZZ_QUBITS_GREEN);
199
200 circuit.insert(circuit.end(), layerX.begin(), layerX.end());
201 circuit.insert(circuit.end(), layerRedZZ.begin(), layerRedZZ.end());
202 circuit.insert(circuit.end(), layerBlueZZ.begin(), layerBlueZZ.end());
203 circuit.insert(circuit.end(), layerGreenZZ.begin(), layerGreenZZ.end());
204 }
205
206 return circuit;
207}
This demonstration will simulate the circuit in the Heisenberg picture. As such, our input Pauli expansion will be initialized to an experimentally measured observable operator. We’ll later choose observable \(Z_{62}\), but we here define a function which encodes any Pauli string as a sequence of “packed integers” as accepted by the cuPauliProp API.
211// ========================================================================
212// Observable preparation
213// ========================================================================
214
215// This demo simulates the IBM circuit via back-propagating the measurement
216// observable through the adjoint circuit. As such, we encode the measured
217// observable into our initial Pauli expansion, in the format recognised by
218// cuPauliProp. Pauli strings are represented with "packed integers" wherein
219// every bit encodes a Pauli operator upon a corresponding qubit. We maintain
220// two masks which separately encode the position of X and Z Pauli operators,
221// indicated by a set bit at the qubit index, with a common set bit encoding
222// a Y Pauli operator. Simulating more qubits than exist bits in the packed
223// integer type (64) requires using multiple packed integers for each X and Z
224// mask. We store a Pauli string's constituent X and Z masks contiguously in
225// a single array, where the final mask of each per-string is padded with zero
226// bits to be an integer multiple of the packed integer size (64 bits).
227
228// The below function accepts a single Pauli string (i.e. a tensor product of
229// the given Pauli operators at the specified qubit indices) and returns the
230// sequence of packed integers which encode it as per the cuPauliProp API;
231// this sequence can be copied directly to the GPU buffer of a Pauli expansion.
232
233std::vector<cupaulipropPackedIntegerType_t> getPauliStringAsPackedIntegers(
234 std::vector<cupaulipropPauliKind_t> paulis,
235 std::vector<uint32_t> qubits
236) {
237 assert(paulis.size() == qubits.size());
238 assert(*std::max_element(qubits.begin(), qubits.end()) < NUM_CIRCUIT_QUBITS);
239
240 int32_t numPackedInts;
241 HANDLE_CUPP_ERROR(cupaulipropGetNumPackedIntegers(NUM_CIRCUIT_QUBITS, &numPackedInts));
242
243 // A single Pauli string is composed of separate X and Z masks, one after the other
244 std::vector<cupaulipropPackedIntegerType_t> out(numPackedInts * 2, 0);
245 auto xPtr = &out[0];
246 auto zPtr = &out[numPackedInts];
247
248 // Process one input (pauli, qubit) pair at a time
249 for (auto i=0; i<qubits.size(); i++) {
250
251 // The qubit corresponds to a specific bit of a specific packed integer
252 auto numBitsPerPackedInt = 8 * sizeof(cupaulipropPackedIntegerType_t);
253 auto intInd = qubits[i] / numBitsPerPackedInt;
254 auto bitInd = qubits[i] % numBitsPerPackedInt;
255
256 // Overwrite a bit of either the X or Z masks (or both when pauli==Y)
257 if (paulis[i] == CUPAULIPROP_PAULI_X || paulis[i] == CUPAULIPROP_PAULI_Y)
258 xPtr[intInd] = xPtr[intInd] | (1ULL << bitInd);
259 if (paulis[i] == CUPAULIPROP_PAULI_Z || paulis[i] == CUPAULIPROP_PAULI_Y)
260 zPtr[intInd] = zPtr[intInd] | (1ULL << bitInd);
261 }
262
263 return out;
264}
With all convenience functions defined, we are ready to begin the main control flow.
268// ========================================================================
269// Main
270// ========================================================================
271
272// Simulation of the IBM utility experiment proceeds as follows. We setup the
273// cuPauliProp library, attaching a new stream, then proceed to creating two
274// Pauli expansions (since the API is out-of-place, as elaborated upon below).
275// One expansion is initialised to the measured observable of the IBM circuit.
276// We prepare workspace memory, fix truncation hyperparameters, then create
277// the circuit as a list of cuPauliProp operators. We process the circuit in
278// reverse, adjointing each operation, applied upon a newly prepared view of
279// the input expansion, each time checking our dynamically growing memory costs
280// have not exceeded our budgets. Thereafter we compute the overlap between the
281// final back-propagated observable and the experimental initial state (the all-
282// zero state), producing an estimate of the experimental expectation value.
283// Finally, we free all allocated memory like good citizens.
284
285
286int main(int argc, char** argv) {
287 std::cout << "cuPauliProp IBM Heavy-hex Ising Example" << std::endl;
288 std::cout << "=======================================" << std::endl << std::endl;
Our first chore is to setup the cuPauliProp library. For simplicity, we choose not to attach any stream, and therefore use the default stream, such that all cuPauliProp API functions are synchronous.
292 // ========================================================================
293 // Library setup
294 // ========================================================================
295
296 int deviceId = 0;
297 HANDLE_CUDA_ERROR(cudaSetDevice(deviceId));
298
299 // Work on the default stream; use cupaulipropSetStream() to perform all
300 // cuPauliProp operations on a specific stream for asynchronous usage
301 cupaulipropHandle_t handle;
302 HANDLE_CUPP_ERROR(cupaulipropCreate(&handle));
We here choose to use either our pre-decided memory budget, or optionally consult the available device
memory and distribute it between the needed data structures. In total, we will reserve memory for two Pauli
expansions (each containing a separate Pauli string and coefficient buffer), and a workspace buffer. The needed
workspace buffer size is informed by both the capacity of the expansions and the API operations to be later
performed. While it can be precomputed using the API “prepare” functions, such as
cupaulipropPauliExpansionViewPrepareOperatorApplication(), we here merely reserve as much memory for
the workspace as needed by a Pauli expansion, for simplicity.
306 // ========================================================================
307 // Decide memory usage
308 // ========================================================================
309
310 // As outlined in the 'Memory usage' section above, we here either uniformly
311 // allocate all (or a high percentage of) available memory between the needed
312 // memory buffers, or use the pre-decided fixed values. This demo will create
313 // a total of two Pauli expansions (each of which accepts two separate buffers
314 // to store Pauli strings and their corresponding coefficients; these have
315 // different sizes) and one workspace, hence we arrange for an allocation of
316 // five buffers in total.
317
318 size_t expansionPauliMem;
319 size_t expansionCoefMem;
320 size_t workspaceMem;
321 size_t totalUsedMem;
322
323 if (USE_MAX_VRAM) {
324
325 // Find usable device memory
326 size_t totalFreeMem, totalGlobalMem;
327 HANDLE_CUDA_ERROR(cudaMemGetInfo(&totalFreeMem, &totalGlobalMem));
328 size_t totalUsableMem = static_cast<size_t>(totalFreeMem * MAX_VRAM_PERCENT/100);
329
330 // Divide it between the three instances (two expansions, one workspace)
331 size_t instanceMem = totalUsableMem / 3;
332
333 // Determine the ideal ratio between an expansion's Pauli and coef buffers
334 int32_t numPackedInts;
335 HANDLE_CUPP_ERROR(cupaulipropGetNumPackedIntegers(NUM_CIRCUIT_QUBITS, &numPackedInts));
336 size_t pauliMemPerTerm = 2 * numPackedInts * sizeof(cupaulipropPackedIntegerType_t);
337 size_t coefMemPerTerm = sizeof(double);
338 size_t totalMemPerTerm = pauliMemPerTerm + coefMemPerTerm;
339
340 expansionPauliMem = (instanceMem * pauliMemPerTerm) / totalMemPerTerm;
341 expansionCoefMem = (instanceMem * coefMemPerTerm ) / totalMemPerTerm;
342 workspaceMem = instanceMem;
343
344 totalUsedMem = 2*expansionPauliMem + 2*expansionCoefMem + workspaceMem;
345 std::cout << "Dedicated memory: " << MAX_VRAM_PERCENT << "% of " << totalFreeMem;
346 std::cout << " B free = " << totalUsedMem << " B, divided into..." << std::endl;
347
348 } else {
349
350 // Use pre-decided buffer sizes
351 expansionPauliMem = FIXED_EXPANSION_PAULI_MEM;
352 expansionCoefMem = FIXED_EXPANSION_COEF_MEM;
353 workspaceMem = FIXED_WORKSPACE_MEM;
354
355 totalUsedMem = 2*expansionPauliMem + 2*expansionCoefMem + workspaceMem;
356 std::cout << "Dedicated memory: " << totalUsedMem << " B = 60 MiB, divided into..." << std::endl;
357 }
358
359 std::cout << " expansion Pauli buffer: " << expansionPauliMem << " B" << std::endl;
360 std::cout << " expansion coef buffer: " << expansionCoefMem << " B" << std::endl;
361 std::cout << " workspace buffer: " << workspaceMem << " B\n" << std::endl;
It is our responsibility to allocate and initialise the device memory before passing it to the cuPauliProp
Pauli expansion constructor. We create an “input” Pauli expansion initialised to the single Pauli string \(Z_{62}\),
and an “output” Pauli expansion which is initially empty. Notice that we indicate to cupaulipropCreatePauliExpansion()
that the pre-initialised input strings (of which there is one) happen to be sorted and unique, which cuPauliProp
can internally later leverage for potential speedup.
365 // ========================================================================
366 // Pauli expansion preparation
367 // ========================================================================
368
369 // Create buffers for two Pauli expansions, which will serve as 'input' and
370 // 'output' to the out-of-place cuPauliProp API. Note that the capacities of
371 // these buffers constrain the maximum number of Pauli strings maintained
372 // during simulation, and ergo inform the accuracy of the simulation. The
373 // sufficient buffer sizes are specific to the simulated system, and we here
374 // choose a surprisingly small capacity as admitted by the studied circuit.
375
376 void * d_inExpansionPauliBuffer;
377 void * d_outExpansionPauliBuffer;
378 void * d_inExpansionCoefBuffer;
379 void * d_outExpansionCoefBuffer;
380 HANDLE_CUDA_ERROR(cudaMalloc(&d_inExpansionPauliBuffer, expansionPauliMem));
381 HANDLE_CUDA_ERROR(cudaMalloc(&d_inExpansionCoefBuffer, expansionCoefMem));
382 HANDLE_CUDA_ERROR(cudaMalloc(&d_outExpansionPauliBuffer, expansionPauliMem));
383 HANDLE_CUDA_ERROR(cudaMalloc(&d_outExpansionCoefBuffer, expansionCoefMem));
384
385 // Prepare the X and Z masks which encode the experimental observable Z_62,
386 // which has a coefficient of unity, as seen in Figure 4. b) of the IBM work.
387 std::cout << "Observable: Z_62\n" << std::endl;
388 int64_t numObservableTerms = 1;
389 double observableCoef = 1.0;
390 std::vector<cupaulipropPauliKind_t> observablePaulis = {CUPAULIPROP_PAULI_Z};
391 std::vector<uint32_t> observableQubits = {62};
392
393 // Overwrite the 'input' Pauli expansion buffers with the observable data
394 auto observablePackedInts = getPauliStringAsPackedIntegers(observablePaulis, observableQubits);
395 size_t numObservableBytes = observablePackedInts.size() * sizeof(observablePackedInts[0]);
396 HANDLE_CUDA_ERROR(cudaMemcpy(
397 d_inExpansionPauliBuffer, observablePackedInts.data(), numObservableBytes, cudaMemcpyHostToDevice));
398 HANDLE_CUDA_ERROR(cudaMemcpy(
399 d_inExpansionCoefBuffer, &observableCoef, sizeof(observableCoef), cudaMemcpyHostToDevice));
400
401 // Create two Pauli expansions, which will serve as 'input and 'output' to the API.
402 // Because we begin from a real observable coefficient, and our circuit is completely
403 // positive and trace preserving, it is sufficient to use strictly real coefficients
404 // in our expansions, informing dataType below. We indicate that the single prepared
405 // term in the input expansion is technically unique, and the terms sorted, which
406 // permits cuPauliProp to use automatic optimisations during simulation.
407
408 cupaulipropPauliExpansion_t inExpansion;
409 cupaulipropPauliExpansion_t outExpansion;
410
411 int32_t isSorted = 1;
412 int32_t isUnique = 1;
413 cudaDataType_t dataType = CUDA_R_64F;
414
415 HANDLE_CUPP_ERROR(cupaulipropCreatePauliExpansion( // init to above
416 handle, NUM_CIRCUIT_QUBITS,
417 d_inExpansionPauliBuffer, expansionPauliMem,
418 d_inExpansionCoefBuffer, expansionCoefMem,
419 dataType, numObservableTerms, isSorted, isUnique,
420 &inExpansion));
421 HANDLE_CUPP_ERROR(cupaulipropCreatePauliExpansion( // init to empty
422 handle, NUM_CIRCUIT_QUBITS,
423 d_outExpansionPauliBuffer, expansionPauliMem,
424 d_outExpansionCoefBuffer, expansionCoefMem,
425 dataType, 0, 0, 0,
426 &outExpansion));
427
428
We must also ourselves allocate device memory for the workspace.
430 // ========================================================================
431 // Workspace preparation
432 // ========================================================================
433
434 // Some API functions require additional workspace memory which we bind to a
435 // workspace descriptor. Ordinarily we use the 'Prepare' functions to precisely
436 // bound upfront the needed workspace memory, but in this simple demo, we
437 // instead use a workspace memory which we prior know to be sufficient.
438
439 // Create a workspace
440 cupaulipropWorkspaceDescriptor_t workspace;
441 HANDLE_CUPP_ERROR(cupaulipropCreateWorkspaceDescriptor(handle, &workspace));
442
443 void* d_workspaceBuffer;
444 HANDLE_CUDA_ERROR(cudaMalloc(&d_workspaceBuffer, workspaceMem));
445 HANDLE_CUPP_ERROR(cupaulipropWorkspaceSetMemory(
446 handle, workspace,
447 CUPAULIPROP_MEMSPACE_DEVICE, CUPAULIPROP_WORKSPACE_SCRATCH,
448 d_workspaceBuffer, workspaceMem));
449
450 // Note that the 'prepare' functions which check the required workspace memory
451 // will detach this buffer, requiring we re-call SetMemory above. We can avoid
452 // this repeated re-attachment by use of a second, bufferless workspace descri-
453 // ptor which we pass to the 'prepare' functions in lieu of this one.
During subsequent Heisenberg evolution of the observable operator \(Z_{62}\), the number of terms in the Pauli expansion will grow exponentially quickly. To mitigate this growth and keep simulation tractable through the entire circuit, we employ truncation strategies. For the simulated system, it is sufficient to perform truncation after every \(10\) gates, wherein we discard all Pauli strings with corresponding coefficient magnitudes smaller or equal to \(10^{-4}\), or which have a Pauli “weight” (the number of non-identity Paulis therein) exceeding \(8\). This strategy keeps total memory costs below 60 MiB (and so within the capacities of our Pauli expansion and workspace buffers), while ultimately producing an expectation value within the error bars of the error-mitigated experimental value.
457 // ========================================================================
458 // Truncation parameter preparation
459 // ========================================================================
460
461 // The Pauli propagation simulation technique has memory and runtime costs which
462 // (for generic circuits) grow exponentially with the circuit length, which we
463 // curtail through "truncation"; dynamic discarding of Pauli strings in our
464 // expansion which are predicted to contribute negligibly to the final output
465 // expectation value. This file demonstrates simultaneous usage of two truncation
466 // techniques; discarding of Pauli strings with an absolute coefficient less than
467 // 0.0001, or a "Pauli weight" (the number of non-identity operators in the string)
468 // exceeding eight. For example, given a Pauli expansion containing:
469 // 0.1 XYZXYZII + 1E-5 ZIIIIIII + 0.2 XXXXYYYY,
470 // our truncation parameters below would see the latter two strings discarded due
471 // to coefficient and weight truncation respectively.
472
473 cupaulipropCoefficientTruncationParams_t coefTruncParams;
474 coefTruncParams.cutoff = 1E-4;
475
476 cupaulipropPauliWeightTruncationParams_t weightTruncParams;
477 weightTruncParams.cutoff = 8;
478
479 const uint32_t numTruncStrats = 2;
480 cupaulipropTruncationStrategy_t truncStrats[] = {
481 {
482 CUPAULIPROP_TRUNCATION_STRATEGY_COEFFICIENT_BASED,
483 &coefTruncParams
484 },
485 {
486 CUPAULIPROP_TRUNCATION_STRATEGY_PAULI_WEIGHT_BASED,
487 &weightTruncParams
488 }
489 };
490
491 // It is not necessary to perform truncation after every gate, since the
492 // Pauli expansion size may not have grown substantially, and attempting
493 // to truncate may incur superfluous memory enumeration costs. In this
494 // demo, we choose to truncate only after every tenth applied gate. Note
495 // deferring truncation requires additional expansion memory; choosing to
496 // truncate after every gate shrinks this demo's costs to 20 MiB total.
497 const int numGatesBetweenTruncations = 10;
498
499 std::cout << "Coefficient truncation threshold: " << coefTruncParams.cutoff << std::endl;
500 std::cout << "Pauli weight truncation threshold: " << weightTruncParams.cutoff << std::endl;
501 std::cout << "Truncation performed after every: " << numGatesBetweenTruncations << " gates\n" << std::endl;
With the stage set, we are ready to begin simulation! We fix the angle of all X-rotation gates in our circuit to \(\theta_X = \pi/4\), and use \(20\) Trotter repetitions. Our simulated system therefore corresponds to the antepenultimate experimental datum in Figure 4. b of Nature volume 618, pages 500–505 (2023). Let operator \(\hat{U}\) represent the action of this full, forward experimental circuit. The original experimentalists performed \(\hat{U}\) upon the fiducial zero-state \(|0\rangle^{\otimes 127}\), and estimated the expectation value of observable \(Z_{62}\) through repeated sampling.
In contrast, we will here work in the Heisenberg picture, and simulate the adjoint (ergo reversed) circuit upon the observable operator \(Z_{62}\). For each gate, we create a view of the entire input Pauli expansion; we confirm that the output Pauli expansion has sufficient capacity to store the (worst case) result of applying the gate upon the view; we reattach the workspace buffer which is detached by that query; we apply the gate and overwrite the output Pauli expansion; and finally, we swap the input and output Pauli expansions so that the updated expansion becomes the next input.
505 // ========================================================================
506 // Back-propagation of the observable through the circuit
507 // ========================================================================
508
509 // We now simulate the observable operator being back-propagated through the
510 // adjoint circuit, mapping the input expansion (initialised to Z_62) to a
511 // final output expansion containing many weighted Pauli strings. We use the
512 // heavy-hex fixed-angle Ising circuit with 20 total repetitions, fixing the
513 // angle of the X rotation gates to PI/4. Our simulation therefore corresponds
514 // to the middle datum of Fig. 4 b) of the IBM manuscript, for which MPS and
515 // isoTNS siulation techniques showed the greatest divergence from experiment.
516
517 double xRotationAngle = PI / 4.;
518 int numTrotterSteps = 20;
519 auto circuit = getIBMHeavyHexIsingCircuit(handle, xRotationAngle, numTrotterSteps);
520
521 std::cout << "Circuit: 127 qubit IBM heavy-hex Ising circuit, with..." << std::endl;
522 std::cout << " Trotter steps: " << numTrotterSteps << std::endl;
523 std::cout << " Total gates: " << circuit.size() << std::endl;
524 std::cout << " Rx angle: " << xRotationAngle << " (i.e. PI/4)\n" << std::endl;
525
526 // Constrain that every intermediate output expansion contains unique Pauli
527 // strings (forbidding duplicates), but permit the retained strings to be
528 // unsorted. This combination gives cuPauliProp the best chance of automatically
529 // selecting optimal internal functions and postconditions for the simulation.
530 uint32_t adjoint = true;
531 uint32_t makeSorted = false;
532 uint32_t keepDuplicates = false;
533
534 std::cout << "Imposed postconditions:" << std::endl;
535 if (makeSorted) {
536 std::cout << " Pauli strings will be sorted." << std::endl;
537 }
538 if (!keepDuplicates) {
539 std::cout << " Pauli strings will be unique." << std::endl;
540 }
541 if (keepDuplicates && !makeSorted) {
542 std::cout << "No postconditions imposed on Pauli strings." << std::endl;
543 }
544 std::cout << std::endl;
545
546 // Begin timing before any gates are applied
547 auto startTime = std::chrono::high_resolution_clock::now();
548 int64_t maxNumTerms = 0;
549
550 // Iterate the circuit in reverse to effect the adjoint of the total circuit
551 for (int gateInd=circuit.size()-1; gateInd >= 0; --gateInd) {
552 cupaulipropQuantumOperator_t gate = circuit[gateInd];
553
554 // Create a view of the current input expansion, selecting all currently
555 // contained terms. For very large systems, we may have alternatively
556 // chosen a smaller view of the partial state to work around memory limits.
557 cupaulipropPauliExpansionView_t inView;
558 int64_t numExpansionTerms;
559 HANDLE_CUPP_ERROR(cupaulipropPauliExpansionGetNumTerms(handle, inExpansion, &numExpansionTerms));
560 HANDLE_CUPP_ERROR(cupaulipropPauliExpansionGetContiguousRange(
561 handle, inExpansion, 0, numExpansionTerms, &inView));
562
563 // Track the intermediate expansion size, for our curiousity
564 if (numExpansionTerms > maxNumTerms)
565 maxNumTerms = numExpansionTerms;
566
567 // Choose whether or not to perform truncations after this gate
568 int numPassedTruncStrats = (gateInd % numGatesBetweenTruncations == 0)? numTruncStrats : 0;
569
570 // Check the expansion and workspace memories needed to apply the current gate
571 int64_t reqExpansionPauliMem;
572 int64_t reqExpansionCoefMem;
573 int64_t reqWorkspaceMem;
574 HANDLE_CUPP_ERROR(cupaulipropPauliExpansionViewPrepareOperatorApplication(
575 handle, inView, gate, makeSorted, keepDuplicates,
576 numPassedTruncStrats, numPassedTruncStrats > 0 ? truncStrats : nullptr,
577 workspaceMem,
578 &reqExpansionPauliMem, &reqExpansionCoefMem, workspace));
579 HANDLE_CUPP_ERROR(cupaulipropWorkspaceGetMemorySize(
580 handle, workspace,
581 CUPAULIPROP_MEMSPACE_DEVICE, CUPAULIPROP_WORKSPACE_SCRATCH,
582 &reqWorkspaceMem));
583
584 // Verify that our existing buffers and workspace have sufficient memory
585 assert(reqExpansionPauliMem <= expansionPauliMem);
586 assert(reqExpansionCoefMem <= expansionCoefMem);
587 assert(reqWorkspaceMem <= workspaceMem);
588
589 // Beware that cupaulipropPauliExpansionViewPrepareOperatorApplication() above
590 // detaches the memory buffer from the workspace, which we here re-attach.
591 HANDLE_CUPP_ERROR(cupaulipropWorkspaceSetMemory(
592 handle, workspace,
593 CUPAULIPROP_MEMSPACE_DEVICE, CUPAULIPROP_WORKSPACE_SCRATCH,
594 d_workspaceBuffer, workspaceMem));
595
596 // Apply the gate upon the prepared view of the input expansion, evolving the
597 // Pauli strings pointed to within, truncating the result. The input expansion
598 // is unchanged while the output expansion is entirely overwritten.
599 HANDLE_CUPP_ERROR(cupaulipropPauliExpansionViewComputeOperatorApplication(
600 handle, inView, outExpansion, gate,
601 adjoint, makeSorted, keepDuplicates,
602 numPassedTruncStrats, numPassedTruncStrats > 0 ? truncStrats : nullptr,
603 workspace));
604
605 // Free the temporary view since it points to the old input expansion, whereas
606 // we will subsequently treat the modified output expansion as the next input
607 HANDLE_CUPP_ERROR(cupaulipropDestroyPauliExpansionView(inView));
608
609 // Treat outExpansion as the input in the next gate application
610 std::swap(inExpansion, outExpansion);
611 }
612
613 // Restore outExpansion to being the final output for clarity
614 std::swap(inExpansion, outExpansion);
Our desired output quantity is the expectation value of \(Z_{62}\), as given by
This is expressible in terms of our final output Pauli expansion \(s_{\text{out}}\) (the Heisenberg-evolved observable) as
Therefore we calculate the trace of our output Pauli expansion with the zero state, producing our estimate to the experimental expectation value of \(Z_{62}\).
618 // ========================================================================
619 // Evaluation of the the expectation value of observable
620 // ========================================================================
621
622 // The output expansion is now a proxy for the observable back-propagated
623 // through to the front of the circuit (though having discarded components
624 // which negligibly influence the subsequent overlap). The expectation value
625 // of the IBM experiment is the overlap of the output expansion with the
626 // zero state, i.e. Tr(outExpansion * |0><0|), as we now compute.
627
628 // Obtain a view of the full output expansion (we'll free it in 'Clean up')
629 cupaulipropPauliExpansionView_t outView;
630 int64_t numOutTerms;
631 HANDLE_CUPP_ERROR(cupaulipropPauliExpansionGetNumTerms(handle, inExpansion, &numOutTerms));
632 HANDLE_CUPP_ERROR(cupaulipropPauliExpansionGetContiguousRange(
633 handle, inExpansion, 0, numOutTerms, &outView));
634
635 // Check that the existing workspace memory is sufficient to compute the trace
636 int64_t reqWorkspaceMem;
637 HANDLE_CUPP_ERROR(cupaulipropPauliExpansionViewPrepareTraceWithZeroState(
638 handle, outView, workspaceMem, workspace));
639 HANDLE_CUPP_ERROR(cupaulipropWorkspaceGetMemorySize(
640 handle, workspace,
641 CUPAULIPROP_MEMSPACE_DEVICE, CUPAULIPROP_WORKSPACE_SCRATCH,
642 &reqWorkspaceMem));
643 assert(reqWorkspaceMem <= workspaceMem);
644
645 // Beware that we must now reattach the buffer to the workspace
646 HANDLE_CUPP_ERROR(cupaulipropWorkspaceSetMemory(
647 handle, workspace,
648 CUPAULIPROP_MEMSPACE_DEVICE, CUPAULIPROP_WORKSPACE_SCRATCH,
649 d_workspaceBuffer, workspaceMem));
650
651 // Compute the trace; the main and final output of this simulation!
652 double expec;
653 HANDLE_CUPP_ERROR(cupaulipropPauliExpansionViewComputeTraceWithZeroState(
654 handle, outView, &expec, workspace));
655
656 // End timing after trace is evaluated
657 auto endTime = std::chrono::high_resolution_clock::now();
658 auto duration = std::chrono::duration_cast<std::chrono::microseconds>(endTime - startTime);
659 auto durationSecs = (duration.count() / 1e6);
660
661 std::cout << "Expectation value: " << expec << std::endl;
662 std::cout << "Final number of terms: " << numOutTerms << std::endl;
663 std::cout << "Maximum number of terms: " << maxNumTerms << std::endl;
664 std::cout << "Runtime: " << durationSecs << " seconds\n" << std::endl;
Finally, we clean up resources and free all device memory reserved during our simulation. This includes any outstanding views, all quantum operators, all expansions and their underlying buffers, and similarly for the workspace. Our final step is to finalise the cuPauliProp library, concluding simulation.
668 // ========================================================================
669 // Clean up
670 // ========================================================================
671
672 HANDLE_CUPP_ERROR(cupaulipropDestroyPauliExpansionView(outView));
673
674 for (auto & gate : circuit) {
675 HANDLE_CUPP_ERROR(cupaulipropDestroyOperator(gate));
676 }
677
678 HANDLE_CUPP_ERROR(cupaulipropDestroyWorkspaceDescriptor(workspace));
679 HANDLE_CUDA_ERROR(cudaFree(d_workspaceBuffer));
680
681 HANDLE_CUPP_ERROR(cupaulipropDestroyPauliExpansion(inExpansion));
682 HANDLE_CUPP_ERROR(cupaulipropDestroyPauliExpansion(outExpansion));
683
684 HANDLE_CUDA_ERROR(cudaFree(d_inExpansionPauliBuffer));
685 HANDLE_CUDA_ERROR(cudaFree(d_outExpansionPauliBuffer));
686 HANDLE_CUDA_ERROR(cudaFree(d_inExpansionCoefBuffer));
687 HANDLE_CUDA_ERROR(cudaFree(d_outExpansionCoefBuffer));
688
689 HANDLE_CUPP_ERROR(cupaulipropDestroy(handle));
690
691 return EXIT_SUCCESS;
692}
Useful tips#
For debugging, one can set the environment variable
CUPAULIPROP_LOG_LEVEL=n. The leveln= 0, 1, …, 5 corresponds to the logger level as described in the table below. The environment variableCUPAULIPROP_LOG_FILE=<filepath>can be used to redirect the log output to a custom file at<filepath>instead ofstdout.
Level |
Summary |
Long Description |
0 |
Off |
Logging is disabled (default) |
1 |
Errors |
Only errors will be logged |
2 |
Performance Trace |
API calls that launch CUDA kernels will log their parameters and important information |
3 |
Performance Hints |
Hints that can potentially improve the application’s performance |
4 |
Heuristics Trace |
Provides general information about the library execution, may contain details about heuristic status |
5 |
API Trace |
API calls will log their parameter and important information |