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