xref: /libCEED/examples/petsc/bps.c (revision 023b8a51676743c24bb03c40af89971dbec6e8fb)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 //                        libCEED + PETSc Example: CEED BPs
9 //
10 // This example demonstrates a simple usage of libCEED with PETSc to solve the CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps.
11 //
12 // The code uses higher level communication protocols in DMPlex.
13 //
14 // Build with:
15 //
16 //     make bps [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
17 //
18 // Sample runs:
19 //
20 //     ./bps -problem bp1 -degree 3
21 //     ./bps -problem bp2 -degree 3
22 //     ./bps -problem bp3 -degree 3
23 //     ./bps -problem bp4 -degree 3
24 //     ./bps -problem bp5 -degree 3 -ceed /cpu/self
25 //     ./bps -problem bp6 -degree 3 -ceed /gpu/cuda
26 //
27 //TESTARGS -ceed {ceed_resource} -test -problem bp5 -degree 3 -ksp_max_it_clip 15,15
28 //TESTARGS -ceed {ceed_resource} -test -problem bp3 -degree 3 -ksp_max_it_clip 50,50 -simplex
29 
30 /// @file
31 /// CEED BPs example using PETSc with DMPlex
32 /// See bpsraw.c for a "raw" implementation using a structured grid.
33 const char help[] = "Solve CEED BPs using PETSc with DMPlex\n";
34 
35 #include "bps.h"
36 
37 #include <ceed.h>
38 #include <petsc.h>
39 #include <petscdmplex.h>
40 #include <petscksp.h>
41 #include <petscsys.h>
42 #include <stdbool.h>
43 #include <string.h>
44 
45 #include "include/bpsproblemdata.h"
46 #include "include/libceedsetup.h"
47 #include "include/matops.h"
48 #include "include/petscutils.h"
49 #include "include/petscversion.h"
50 #include "include/structs.h"
51 
52 #if PETSC_VERSION_LT(3, 12, 0)
53 #ifdef PETSC_HAVE_CUDA
54 #include <petsccuda.h>
55 // Note: With PETSc prior to version 3.12.0, providing the source path to include 'cublas_v2.h' will be needed to use 'petsccuda.h'.
56 #endif
57 #endif
58 
59 // -----------------------------------------------------------------------------
60 // Main body of program, called in a loop for performance benchmarking purposes
61 // -----------------------------------------------------------------------------
62 static PetscErrorCode RunWithDM(RunParams rp, DM dm, const char *ceed_resource) {
63   double               my_rt_start, my_rt, rt_min, rt_max;
64   PetscInt             xl_size, l_size, g_size;
65   PetscScalar         *r;
66   Vec                  X, X_loc, rhs, rhs_loc;
67   Mat                  mat_O;
68   KSP                  ksp;
69   OperatorApplyContext op_apply_ctx, op_error_ctx;
70   Ceed                 ceed;
71   CeedData             ceed_data;
72   CeedQFunction        qf_error;
73   CeedOperator         op_error;
74   CeedVector           rhs_ceed, target;
75   VecType              vec_type;
76   PetscMemType         mem_type;
77 
78   PetscFunctionBeginUser;
79   // Set up libCEED
80   CeedInit(ceed_resource, &ceed);
81   CeedMemType mem_type_backend;
82   CeedGetPreferredMemType(ceed, &mem_type_backend);
83 
84   PetscCall(DMGetVecType(dm, &vec_type));
85   if (!vec_type) {  // Not yet set by user -dm_vec_type
86     switch (mem_type_backend) {
87       case CEED_MEM_HOST:
88         vec_type = VECSTANDARD;
89         break;
90       case CEED_MEM_DEVICE: {
91         const char *resolved;
92         CeedGetResource(ceed, &resolved);
93         if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA;
94         else if (strstr(resolved, "/gpu/hip/occa")) vec_type = VECSTANDARD;  // https://github.com/CEED/libCEED/issues/678
95         else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP;
96         else vec_type = VECSTANDARD;
97       }
98     }
99     PetscCall(DMSetVecType(dm, vec_type));
100   }
101 
102   // Create global and local solution vectors
103   PetscCall(DMCreateGlobalVector(dm, &X));
104   PetscCall(VecGetLocalSize(X, &l_size));
105   PetscCall(VecGetSize(X, &g_size));
106   PetscCall(DMCreateLocalVector(dm, &X_loc));
107   PetscCall(VecGetSize(X_loc, &xl_size));
108   PetscCall(VecDuplicate(X, &rhs));
109 
110   // Operator
111   PetscCall(PetscMalloc1(1, &op_apply_ctx));
112   PetscCall(PetscMalloc1(1, &op_error_ctx));
113   PetscCall(MatCreateShell(rp->comm, l_size, l_size, g_size, g_size, op_apply_ctx, &mat_O));
114   PetscCall(MatShellSetOperation(mat_O, MATOP_MULT, (void (*)(void))MatMult_Ceed));
115   PetscCall(MatShellSetOperation(mat_O, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag));
116   PetscCall(MatShellSetVecType(mat_O, vec_type));
117 
118   // Print summary
119   if (!rp->test_mode) {
120     PetscInt P = rp->degree + 1, Q = P + rp->q_extra;
121 
122     const char *used_resource;
123     CeedGetResource(ceed, &used_resource);
124 
125     VecType vec_type;
126     PetscCall(VecGetType(X, &vec_type));
127 
128     PetscInt c_start, c_end;
129     PetscCall(DMPlexGetHeightStratum(dm, 0, &c_start, &c_end));
130     DMPolytopeType cell_type;
131     PetscCall(DMPlexGetCellType(dm, c_start, &cell_type));
132     CeedElemTopology elem_topo = ElemTopologyP2C(cell_type);
133     PetscMPIInt      comm_size;
134     PetscCall(MPI_Comm_size(rp->comm, &comm_size));
135     PetscCall(PetscPrintf(rp->comm,
136                           "\n-- CEED Benchmark Problem %" CeedInt_FMT " -- libCEED + PETSc --\n"
137                           "  MPI:\n"
138                           "    Hostname                                : %s\n"
139                           "    Total ranks                             : %d\n"
140                           "    Ranks per compute node                  : %d\n"
141                           "  PETSc:\n"
142                           "    PETSc Vec Type                          : %s\n"
143                           "  libCEED:\n"
144                           "    libCEED Backend                         : %s\n"
145                           "    libCEED Backend MemType                 : %s\n"
146                           "  Mesh:\n"
147                           "    Solution Order (P)                      : %" CeedInt_FMT "\n"
148                           "    Quadrature  Order (Q)                   : %" CeedInt_FMT "\n"
149                           "    Additional quadrature points (q_extra)  : %" CeedInt_FMT "\n"
150                           "    Global nodes                            : %" PetscInt_FMT "\n"
151                           "    Local Elements                          : %" PetscInt_FMT "\n"
152                           "    Element topology                        : %s\n"
153                           "    Owned nodes                             : %" PetscInt_FMT "\n"
154                           "    DoF per node                            : %" PetscInt_FMT "\n",
155                           rp->bp_choice + 1, rp->hostname, comm_size, rp->ranks_per_node, vec_type, used_resource, CeedMemTypes[mem_type_backend], P,
156                           Q, rp->q_extra, g_size / rp->num_comp_u, c_end - c_start, CeedElemTopologies[elem_topo], l_size / rp->num_comp_u,
157                           rp->num_comp_u));
158   }
159 
160   // Create RHS vector
161   PetscCall(VecDuplicate(X_loc, &rhs_loc));
162   PetscCall(VecZeroEntries(rhs_loc));
163   PetscCall(VecGetArrayAndMemType(rhs_loc, &r, &mem_type));
164   CeedVectorCreate(ceed, xl_size, &rhs_ceed);
165   CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r);
166 
167   PetscCall(PetscMalloc1(1, &ceed_data));
168   PetscCall(SetupLibceedByDegree(dm, ceed, rp->degree, rp->dim, rp->q_extra, rp->dim, rp->num_comp_u, g_size, xl_size, bp_options[rp->bp_choice],
169                                  ceed_data, true, rhs_ceed, &target));
170 
171   // Gather RHS
172   CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL);
173   PetscCall(VecRestoreArrayAndMemType(rhs_loc, &r));
174   PetscCall(VecZeroEntries(rhs));
175   PetscCall(DMLocalToGlobal(dm, rhs_loc, ADD_VALUES, rhs));
176   CeedVectorDestroy(&rhs_ceed);
177 
178   // Create the error QFunction
179   CeedQFunctionCreateInterior(ceed, 1, bp_options[rp->bp_choice].error, bp_options[rp->bp_choice].error_loc, &qf_error);
180   CeedQFunctionAddInput(qf_error, "u", rp->num_comp_u, CEED_EVAL_INTERP);
181   CeedQFunctionAddInput(qf_error, "true_soln", rp->num_comp_u, CEED_EVAL_NONE);
182   CeedQFunctionAddInput(qf_error, "qdata", ceed_data->q_data_size, CEED_EVAL_NONE);
183   CeedQFunctionAddOutput(qf_error, "error", rp->num_comp_u, CEED_EVAL_INTERP);
184 
185   // Create the error operator
186   CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_error);
187   CeedOperatorSetField(op_error, "u", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE);
188   CeedOperatorSetField(op_error, "true_soln", ceed_data->elem_restr_u_i, CEED_BASIS_COLLOCATED, target);
189   CeedOperatorSetField(op_error, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
190   CeedOperatorSetField(op_error, "error", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE);
191 
192   // Set up apply operator context
193   PetscCall(SetupApplyOperatorCtx(rp->comm, dm, ceed, ceed_data, X_loc, op_apply_ctx));
194   PetscCall(KSPCreate(rp->comm, &ksp));
195   {
196     PC pc;
197     PetscCall(KSPGetPC(ksp, &pc));
198     if (rp->bp_choice == CEED_BP1 || rp->bp_choice == CEED_BP2) {
199       PetscCall(PCSetType(pc, PCJACOBI));
200       if (rp->simplex) {
201         PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL));
202       } else {
203         PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
204       }
205     } else {
206       PetscCall(PCSetType(pc, PCNONE));
207     }
208     PetscCall(KSPSetType(ksp, KSPCG));
209     PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
210     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
211   }
212   PetscCall(KSPSetOperators(ksp, mat_O, mat_O));
213 
214   // First run's performance log is not considered for benchmarking purposes
215   PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1));
216   my_rt_start = MPI_Wtime();
217   PetscCall(KSPSolve(ksp, rhs, X));
218   my_rt = MPI_Wtime() - my_rt_start;
219   PetscCall(MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, rp->comm));
220   // Set maxits based on first iteration timing
221   if (my_rt > 0.02) {
222     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, rp->ksp_max_it_clip[0]));
223   } else {
224     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, rp->ksp_max_it_clip[1]));
225   }
226   PetscCall(KSPSetFromOptions(ksp));
227 
228   // Timed solve
229   PetscCall(VecZeroEntries(X));
230   PetscCall(PetscBarrier((PetscObject)ksp));
231 
232   // -- Performance logging
233   PetscCall(PetscLogStagePush(rp->solve_stage));
234 
235   // -- Solve
236   my_rt_start = MPI_Wtime();
237   PetscCall(KSPSolve(ksp, rhs, X));
238   my_rt = MPI_Wtime() - my_rt_start;
239 
240   // -- Performance logging
241   PetscCall(PetscLogStagePop());
242 
243   // Output results
244   {
245     KSPType            ksp_type;
246     KSPConvergedReason reason;
247     PetscReal          rnorm;
248     PetscInt           its;
249     PetscCall(KSPGetType(ksp, &ksp_type));
250     PetscCall(KSPGetConvergedReason(ksp, &reason));
251     PetscCall(KSPGetIterationNumber(ksp, &its));
252     PetscCall(KSPGetResidualNorm(ksp, &rnorm));
253     if (!rp->test_mode || reason < 0 || rnorm > 1e-8) {
254       PetscCall(PetscPrintf(rp->comm,
255                             "  KSP:\n"
256                             "    KSP Type                                : %s\n"
257                             "    KSP Convergence                         : %s\n"
258                             "    Total KSP Iterations                    : %" PetscInt_FMT "\n"
259                             "    Final rnorm                             : %e\n",
260                             ksp_type, KSPConvergedReasons[reason], its, (double)rnorm));
261     }
262     if (!rp->test_mode) {
263       PetscCall(PetscPrintf(rp->comm, "  Performance:\n"));
264     }
265     {
266       // Set up error operator context
267       PetscCall(SetupErrorOperatorCtx(rp->comm, dm, ceed, ceed_data, X_loc, op_error, op_error_ctx));
268       PetscScalar l2_error;
269       PetscCall(ComputeL2Error(X, &l2_error, op_error_ctx));
270       PetscReal tol = 5e-2;
271       if (!rp->test_mode || l2_error > tol) {
272         PetscCall(MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, rp->comm));
273         PetscCall(MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, rp->comm));
274         PetscCall(PetscPrintf(rp->comm,
275                               "    L2 Error                                : %e\n"
276                               "    CG Solve Time                           : %g (%g) sec\n",
277                               (double)l2_error, rt_max, rt_min));
278       }
279     }
280     if (!rp->test_mode) {
281       PetscCall(PetscPrintf(rp->comm, "    DoFs/Sec in CG                          : %g (%g) million\n", 1e-6 * g_size * its / rt_max,
282                             1e-6 * g_size * its / rt_min));
283     }
284   }
285 
286   if (rp->write_solution) {
287     PetscViewer vtk_viewer_soln;
288 
289     PetscCall(PetscViewerCreate(rp->comm, &vtk_viewer_soln));
290     PetscCall(PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK));
291     PetscCall(PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"));
292     PetscCall(VecView(X, vtk_viewer_soln));
293     PetscCall(PetscViewerDestroy(&vtk_viewer_soln));
294   }
295 
296   // Cleanup
297   PetscCall(VecDestroy(&X));
298   PetscCall(VecDestroy(&X_loc));
299   PetscCall(VecDestroy(&op_apply_ctx->Y_loc));
300   PetscCall(VecDestroy(&op_error_ctx->Y_loc));
301   PetscCall(MatDestroy(&mat_O));
302   PetscCall(PetscFree(op_apply_ctx));
303   PetscCall(PetscFree(op_error_ctx));
304   PetscCall(CeedDataDestroy(0, ceed_data));
305 
306   PetscCall(VecDestroy(&rhs));
307   PetscCall(VecDestroy(&rhs_loc));
308   PetscCall(KSPDestroy(&ksp));
309   CeedVectorDestroy(&target);
310   CeedQFunctionDestroy(&qf_error);
311   CeedOperatorDestroy(&op_error);
312   CeedDestroy(&ceed);
313   PetscFunctionReturn(0);
314 }
315 
316 static PetscErrorCode Run(RunParams rp, PetscInt num_resources, char *const *ceed_resources, PetscInt num_bp_choices, const BPType *bp_choices) {
317   DM dm;
318 
319   PetscFunctionBeginUser;
320   // Setup DM
321   PetscCall(CreateDistributedDM(rp, &dm));
322 
323   for (PetscInt b = 0; b < num_bp_choices; b++) {
324     DM       dm_deg;
325     VecType  vec_type;
326     PetscInt q_extra = rp->q_extra;
327     rp->bp_choice    = bp_choices[b];
328     rp->num_comp_u   = bp_options[rp->bp_choice].num_comp_u;
329     rp->q_extra      = q_extra < 0 ? bp_options[rp->bp_choice].q_extra : q_extra;
330     PetscCall(DMClone(dm, &dm_deg));
331     PetscCall(DMGetVecType(dm, &vec_type));
332     PetscCall(DMSetVecType(dm_deg, vec_type));
333     // Create DM
334     PetscInt dim;
335     PetscCall(DMGetDimension(dm_deg, &dim));
336     PetscCall(SetupDMByDegree(dm_deg, rp->degree, rp->q_extra, rp->num_comp_u, dim, bp_options[rp->bp_choice].enforce_bc));
337     for (PetscInt r = 0; r < num_resources; r++) {
338       PetscCall(RunWithDM(rp, dm_deg, ceed_resources[r]));
339     }
340     PetscCall(DMDestroy(&dm_deg));
341     rp->q_extra = q_extra;
342   }
343 
344   PetscCall(DMDestroy(&dm));
345   PetscFunctionReturn(0);
346 }
347 
348 int main(int argc, char **argv) {
349   PetscInt  comm_size;
350   RunParams rp;
351   MPI_Comm  comm;
352   char      filename[PETSC_MAX_PATH_LEN];
353   char     *ceed_resources[30];
354   PetscInt  num_ceed_resources = 30;
355   char      hostname[PETSC_MAX_PATH_LEN];
356 
357   PetscInt    dim = 3, mesh_elem[3] = {3, 3, 3};
358   PetscInt    num_degrees = 30, degree[30] = {}, num_local_nodes = 2, local_nodes[2] = {};
359   PetscMPIInt ranks_per_node;
360   PetscBool   degree_set;
361   BPType      bp_choices[10];
362   PetscInt    num_bp_choices = 10;
363 
364   // Initialize PETSc
365   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
366   comm = PETSC_COMM_WORLD;
367   PetscCall(MPI_Comm_size(comm, &comm_size));
368 #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
369   {
370     MPI_Comm splitcomm;
371     PetscCall(MPI_Comm_split_type(comm, MPI_COMM_TYPE_SHARED, 0, MPI_INFO_NULL, &splitcomm));
372     PetscCall(MPI_Comm_size(splitcomm, &ranks_per_node));
373     PetscCall(MPI_Comm_free(&splitcomm));
374   }
375 #else
376   ranks_per_node = -1;  // Unknown
377 #endif
378 
379   // Setup all parameters needed in Run()
380   PetscCall(PetscMalloc1(1, &rp));
381   rp->comm = comm;
382 
383   // Read command line options
384   PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
385   {
386     PetscBool set;
387     PetscCall(PetscOptionsEnumArray("-problem", "CEED benchmark problem to solve", NULL, bp_types, (PetscEnum *)bp_choices, &num_bp_choices, &set));
388     if (!set) {
389       bp_choices[0]  = CEED_BP1;
390       num_bp_choices = 1;
391     }
392   }
393   rp->test_mode = PETSC_FALSE;
394   PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, rp->test_mode, &rp->test_mode, NULL));
395   rp->write_solution = PETSC_FALSE;
396   PetscCall(PetscOptionsBool("-write_solution", "Write solution for visualization", NULL, rp->write_solution, &rp->write_solution, NULL));
397   rp->simplex = PETSC_FALSE;
398   PetscCall(PetscOptionsBool("-simplex", "Element topology (default:hex)", NULL, rp->simplex, &rp->simplex, NULL));
399   if ((bp_choices[0] == CEED_BP5 || bp_choices[0] == CEED_BP6) && (rp->simplex)) {
400     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "BP5/6 is not supported with simplex");
401   }
402   degree[0] = rp->test_mode ? 3 : 2;
403   PetscCall(PetscOptionsIntArray("-degree", "Polynomial degree of tensor product basis", NULL, degree, &num_degrees, &degree_set));
404   if (!degree_set) num_degrees = 1;
405   rp->q_extra = PETSC_DECIDE;
406   PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points (-1 for auto)", NULL, rp->q_extra, &rp->q_extra, NULL));
407   {
408     PetscBool set;
409     PetscCall(PetscOptionsStringArray("-ceed", "CEED resource specifier (comma-separated list)", NULL, ceed_resources, &num_ceed_resources, &set));
410     if (!set) {
411       PetscCall(PetscStrallocpy("/cpu/self", &ceed_resources[0]));
412       num_ceed_resources = 1;
413     }
414   }
415   PetscCall(PetscGetHostName(hostname, sizeof hostname));
416   PetscCall(PetscOptionsString("-hostname", "Hostname for output", NULL, hostname, hostname, sizeof(hostname), NULL));
417   rp->read_mesh = PETSC_FALSE;
418   PetscCall(PetscOptionsString("-mesh", "Read mesh from file", NULL, filename, filename, sizeof(filename), &rp->read_mesh));
419   rp->filename = filename;
420   if (!rp->read_mesh) {
421     PetscInt tmp = dim;
422     PetscCall(PetscOptionsIntArray("-cells", "Number of cells per dimension", NULL, mesh_elem, &tmp, NULL));
423   }
424   local_nodes[0] = 1000;
425   PetscCall(PetscOptionsIntArray("-local_nodes",
426                                  "Target number of locally owned nodes per "
427                                  "process (single value or min,max)",
428                                  NULL, local_nodes, &num_local_nodes, &rp->user_l_nodes));
429   if (num_local_nodes < 2) local_nodes[1] = 2 * local_nodes[0];
430   {
431     PetscInt two           = 2;
432     rp->ksp_max_it_clip[0] = 5;
433     rp->ksp_max_it_clip[1] = 20;
434     PetscCall(PetscOptionsIntArray("-ksp_max_it_clip", "Min and max number of iterations to use during benchmarking", NULL, rp->ksp_max_it_clip, &two,
435                                    NULL));
436   }
437   if (!degree_set) {
438     PetscInt max_degree = 8;
439     PetscCall(PetscOptionsInt("-max_degree", "Range of degrees [1, max_degree] to run with", NULL, max_degree, &max_degree, NULL));
440     for (PetscInt i = 0; i < max_degree; i++) degree[i] = i + 1;
441     num_degrees = max_degree;
442   }
443   {
444     PetscBool flg;
445     PetscInt  p = ranks_per_node;
446     PetscCall(PetscOptionsInt("-p", "Number of MPI ranks per node", NULL, p, &p, &flg));
447     if (flg) ranks_per_node = p;
448   }
449 
450   PetscOptionsEnd();
451 
452   // Register PETSc logging stage
453   PetscCall(PetscLogStageRegister("Solve Stage", &rp->solve_stage));
454 
455   rp->hostname       = hostname;
456   rp->dim            = dim;
457   rp->mesh_elem      = mesh_elem;
458   rp->ranks_per_node = ranks_per_node;
459 
460   for (PetscInt d = 0; d < num_degrees; d++) {
461     PetscInt deg = degree[d];
462     for (PetscInt n = local_nodes[0]; n < local_nodes[1]; n *= 2) {
463       rp->degree      = deg;
464       rp->local_nodes = n;
465       PetscCall(Run(rp, num_ceed_resources, ceed_resources, num_bp_choices, bp_choices));
466     }
467   }
468   // Clear memory
469   PetscCall(PetscFree(rp));
470   for (PetscInt i = 0; i < num_ceed_resources; i++) {
471     PetscCall(PetscFree(ceed_resources[i]));
472   }
473   return PetscFinalize();
474 }
475