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