xref: /libCEED/examples/petsc/bps.c (revision 7bfe0f0e497883534890a072c2a8b865352898b0)
1 // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3 // reserved. See files LICENSE and NOTICE for details.
4 //
5 // This file is part of CEED, a collection of benchmarks, miniapps, software
6 // libraries and APIs for efficient high-order finite element and spectral
7 // element discretizations for exascale applications. For more information and
8 // source code availability see http://github.com/ceed.
9 //
10 // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11 // a collaborative effort of two U.S. Department of Energy organizations (Office
12 // of Science and the National Nuclear Security Administration) responsible for
13 // the planning and preparation of a capable exascale ecosystem, including
14 // software, applications, hardware, advanced system engineering and early
15 // testbed platforms, in support of the nation's exascale computing imperative.
16 
17 //                        libCEED + PETSc Example: CEED BPs
18 //
19 // This example demonstrates a simple usage of libCEED with PETSc to solve the
20 // CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps.
21 //
22 // The code uses higher level communication protocols in DMPlex.
23 //
24 // Build with:
25 //
26 //     make bps [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
27 //
28 // Sample runs:
29 //
30 //     ./bps -problem bp1 -degree 3
31 //     ./bps -problem bp2 -degree 3
32 //     ./bps -problem bp3 -degree 3
33 //     ./bps -problem bp4 -degree 3
34 //     ./bps -problem bp5 -degree 3 -ceed /cpu/self
35 //     ./bps -problem bp6 -degree 3 -ceed /gpu/cuda
36 //
37 //TESTARGS -ceed {ceed_resource} -test -problem bp5 -degree 3 -ksp_max_it_clip 15,15
38 
39 /// @file
40 /// CEED BPs example using PETSc with DMPlex
41 /// See bpsraw.c for a "raw" implementation using a structured grid.
42 const char help[] = "Solve CEED BPs using PETSc with DMPlex\n";
43 
44 #include <stdbool.h>
45 #include <string.h>
46 #include <ceed.h>
47 #include <petsc.h>
48 #include <petscdmplex.h>
49 #include <petscksp.h>
50 #include <petscsys.h>
51 
52 #include "bps.h"
53 #include "include/bpsproblemdata.h"
54 #include "include/petscutils.h"
55 #include "include/petscversion.h"
56 #include "include/matops.h"
57 #include "include/structs.h"
58 #include "include/libceedsetup.h"
59 
60 #if PETSC_VERSION_LT(3,12,0)
61 #ifdef PETSC_HAVE_CUDA
62 #include <petsccuda.h>
63 // Note: With PETSc prior to version 3.12.0, providing the source path to
64 //       include 'cublas_v2.h' will be needed to use 'petsccuda.h'.
65 #endif
66 #endif
67 
68 // -----------------------------------------------------------------------------
69 // Utilities
70 // -----------------------------------------------------------------------------
71 
72 // Utility function, compute three factors of an integer
73 static void Split3(PetscInt size, PetscInt m[3], bool reverse) {
74   for (PetscInt d=0, size_left=size; d<3; d++) {
75     PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(size_left, 1./(3 - d)));
76     while (try * (size_left / try) != size_left) try++;
77     m[reverse ? 2-d : d] = try;
78     size_left /= try;
79   }
80 }
81 
82 static int Max3(const PetscInt a[3]) {
83   return PetscMax(a[0], PetscMax(a[1], a[2]));
84 }
85 
86 static int Min3(const PetscInt a[3]) {
87   return PetscMin(a[0], PetscMin(a[1], a[2]));
88 }
89 
90 // -----------------------------------------------------------------------------
91 // Parameter structure for running problems
92 // -----------------------------------------------------------------------------
93 typedef struct RunParams_ *RunParams;
94 struct RunParams_ {
95   MPI_Comm comm;
96   PetscBool test_mode, read_mesh, user_l_nodes, write_solution;
97   char *filename, *hostname;
98   PetscInt local_nodes, degree, q_extra, dim, num_comp_u, *mesh_elem;
99   PetscInt ksp_max_it_clip[2];
100   PetscMPIInt ranks_per_node;
101   BPType bp_choice;
102   PetscLogStage solve_stage;
103 };
104 
105 // -----------------------------------------------------------------------------
106 // Main body of program, called in a loop for performance benchmarking purposes
107 // -----------------------------------------------------------------------------
108 static PetscErrorCode RunWithDM(RunParams rp, DM dm,
109                                 const char *ceed_resource) {
110   PetscErrorCode ierr;
111   double my_rt_start, my_rt, rt_min, rt_max;
112   PetscInt xl_size, l_size, g_size;
113   PetscScalar *r;
114   Vec X, X_loc, rhs, rhs_loc;
115   Mat mat_O;
116   KSP ksp;
117   UserO user_O;
118   Ceed ceed;
119   CeedData ceed_data;
120   CeedQFunction qf_error;
121   CeedOperator op_error;
122   CeedVector rhs_ceed, target;
123   VecType vec_type;
124   PetscMemType mem_type;
125 
126   PetscFunctionBeginUser;
127   // Set up libCEED
128   CeedInit(ceed_resource, &ceed);
129   CeedMemType mem_type_backend;
130   CeedGetPreferredMemType(ceed, &mem_type_backend);
131 
132   ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr);
133   if (!vec_type) { // Not yet set by user -dm_vec_type
134     switch (mem_type_backend) {
135     case CEED_MEM_HOST: vec_type = VECSTANDARD; break;
136     case CEED_MEM_DEVICE: {
137       const char *resolved;
138       CeedGetResource(ceed, &resolved);
139       if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA;
140       else if (strstr(resolved, "/gpu/hip/occa"))
141         vec_type = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678
142       else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP;
143       else vec_type = VECSTANDARD;
144     }
145     }
146     ierr = DMSetVecType(dm, vec_type); CHKERRQ(ierr);
147   }
148 
149   // Create global and local solution vectors
150   ierr = DMCreateGlobalVector(dm, &X); CHKERRQ(ierr);
151   ierr = VecGetLocalSize(X, &l_size); CHKERRQ(ierr);
152   ierr = VecGetSize(X, &g_size); CHKERRQ(ierr);
153   ierr = DMCreateLocalVector(dm, &X_loc); CHKERRQ(ierr);
154   ierr = VecGetSize(X_loc, &xl_size); CHKERRQ(ierr);
155   ierr = VecDuplicate(X, &rhs); CHKERRQ(ierr);
156 
157   // Operator
158   ierr = PetscMalloc1(1, &user_O); CHKERRQ(ierr);
159   ierr = MatCreateShell(rp->comm, l_size, l_size, g_size, g_size,
160                         user_O, &mat_O); CHKERRQ(ierr);
161   ierr = MatShellSetOperation(mat_O, MATOP_MULT,
162                               (void(*)(void))MatMult_Ceed); CHKERRQ(ierr);
163   ierr = MatShellSetOperation(mat_O, MATOP_GET_DIAGONAL,
164                               (void(*)(void))MatGetDiag); CHKERRQ(ierr);
165   ierr = MatShellSetVecType(mat_O, vec_type); CHKERRQ(ierr);
166 
167   // Print summary
168   if (!rp->test_mode) {
169     PetscInt P = rp->degree + 1, Q = P + rp->q_extra;
170 
171     const char *used_resource;
172     CeedGetResource(ceed, &used_resource);
173 
174     VecType vec_type;
175     ierr = VecGetType(X, &vec_type); CHKERRQ(ierr);
176 
177     PetscInt c_start, c_end;
178     ierr = DMPlexGetHeightStratum(dm, 0, &c_start, &c_end); CHKERRQ(ierr);
179     PetscMPIInt comm_size;
180     ierr = MPI_Comm_size(rp->comm, &comm_size); CHKERRQ(ierr);
181     ierr = PetscPrintf(rp->comm,
182                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n"
183                        "  MPI:\n"
184                        "    Hostname                           : %s\n"
185                        "    Total ranks                        : %d\n"
186                        "    Ranks per compute node             : %d\n"
187                        "  PETSc:\n"
188                        "    PETSc Vec Type                     : %s\n"
189                        "  libCEED:\n"
190                        "    libCEED Backend                    : %s\n"
191                        "    libCEED Backend MemType            : %s\n"
192                        "  Mesh:\n"
193                        "    Number of 1D Basis Nodes (P)       : %d\n"
194                        "    Number of 1D Quadrature Points (Q) : %d\n"
195                        "    Global nodes                       : %D\n"
196                        "    Local Elements                     : %D\n"
197                        "    Owned nodes                        : %D\n"
198                        "    DoF per node                       : %D\n",
199                        rp->bp_choice+1, rp->hostname, comm_size,
200                        rp->ranks_per_node, vec_type, used_resource,
201                        CeedMemTypes[mem_type_backend],
202                        P, Q, g_size/rp->num_comp_u, c_end - c_start, l_size/rp->num_comp_u,
203                        rp->num_comp_u);
204     CHKERRQ(ierr);
205   }
206 
207   // Create RHS vector
208   ierr = VecDuplicate(X_loc, &rhs_loc); CHKERRQ(ierr);
209   ierr = VecZeroEntries(rhs_loc); CHKERRQ(ierr);
210   ierr = VecGetArrayAndMemType(rhs_loc, &r, &mem_type); CHKERRQ(ierr);
211   CeedVectorCreate(ceed, xl_size, &rhs_ceed);
212   CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r);
213 
214   ierr = PetscMalloc1(1, &ceed_data); CHKERRQ(ierr);
215   ierr = SetupLibceedByDegree(dm, ceed, rp->degree, rp->dim, rp->q_extra,
216                               rp->dim, rp->num_comp_u, g_size, xl_size, bp_options[rp->bp_choice],
217                               ceed_data, true, rhs_ceed, &target); CHKERRQ(ierr);
218 
219   // Gather RHS
220   CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL);
221   ierr = VecRestoreArrayAndMemType(rhs_loc, &r); CHKERRQ(ierr);
222   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
223   ierr = DMLocalToGlobal(dm, rhs_loc, ADD_VALUES, rhs); CHKERRQ(ierr);
224   CeedVectorDestroy(&rhs_ceed);
225 
226   // Create the error QFunction
227   CeedQFunctionCreateInterior(ceed, 1, bp_options[rp->bp_choice].error,
228                               bp_options[rp->bp_choice].error_loc, &qf_error);
229   CeedQFunctionAddInput(qf_error, "u", rp->num_comp_u, CEED_EVAL_INTERP);
230   CeedQFunctionAddInput(qf_error, "true_soln", rp->num_comp_u, CEED_EVAL_NONE);
231   CeedQFunctionAddOutput(qf_error, "error", rp->num_comp_u, CEED_EVAL_NONE);
232 
233   // Create the error operator
234   CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
235                      &op_error);
236   CeedOperatorSetField(op_error, "u", ceed_data->elem_restr_u,
237                        ceed_data->basis_u, CEED_VECTOR_ACTIVE);
238   CeedOperatorSetField(op_error, "true_soln", ceed_data->elem_restr_u_i,
239                        CEED_BASIS_COLLOCATED, target);
240   CeedOperatorSetField(op_error, "error", ceed_data->elem_restr_u_i,
241                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
242 
243   // Set up Mat
244   user_O->comm = rp->comm;
245   user_O->dm = dm;
246   user_O->X_loc = X_loc;
247   ierr = VecDuplicate(X_loc, &user_O->Y_loc); CHKERRQ(ierr);
248   user_O->x_ceed = ceed_data->x_ceed;
249   user_O->y_ceed = ceed_data->y_ceed;
250   user_O->op = ceed_data->op_apply;
251   user_O->ceed = ceed;
252 
253   ierr = KSPCreate(rp->comm, &ksp); CHKERRQ(ierr);
254   {
255     PC pc;
256     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
257     if (rp->bp_choice == CEED_BP1 || rp->bp_choice == CEED_BP2) {
258       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
259       ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr);
260     } else {
261       ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr);
262     }
263     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
264     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
265     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
266                             PETSC_DEFAULT); CHKERRQ(ierr);
267   }
268   ierr = KSPSetOperators(ksp, mat_O, mat_O); CHKERRQ(ierr);
269 
270   // First run's performance log is not considered for benchmarking purposes
271   ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
272   CHKERRQ(ierr);
273   my_rt_start = MPI_Wtime();
274   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
275   my_rt = MPI_Wtime() - my_rt_start;
276   ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, rp->comm);
277   CHKERRQ(ierr);
278   // Set maxits based on first iteration timing
279   if (my_rt > 0.02) {
280     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
281                             rp->ksp_max_it_clip[0]);
282     CHKERRQ(ierr);
283   } else {
284     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
285                             rp->ksp_max_it_clip[1]);
286     CHKERRQ(ierr);
287   }
288   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
289 
290   // Timed solve
291   ierr = VecZeroEntries(X); CHKERRQ(ierr);
292   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
293 
294   // -- Performance logging
295   ierr = PetscLogStagePush(rp->solve_stage); CHKERRQ(ierr);
296 
297   // -- Solve
298   my_rt_start = MPI_Wtime();
299   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
300   my_rt = MPI_Wtime() - my_rt_start;
301 
302   // -- Performance logging
303   ierr = PetscLogStagePop();
304 
305   // Output results
306   {
307     KSPType ksp_type;
308     KSPConvergedReason reason;
309     PetscReal rnorm;
310     PetscInt its;
311     ierr = KSPGetType(ksp, &ksp_type); CHKERRQ(ierr);
312     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
313     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
314     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
315     if (!rp->test_mode || reason < 0 || rnorm > 1e-8) {
316       ierr = PetscPrintf(rp->comm,
317                          "  KSP:\n"
318                          "    KSP Type                           : %s\n"
319                          "    KSP Convergence                    : %s\n"
320                          "    Total KSP Iterations               : %D\n"
321                          "    Final rnorm                        : %e\n",
322                          ksp_type, KSPConvergedReasons[reason], its,
323                          (double)rnorm); CHKERRQ(ierr);
324     }
325     if (!rp->test_mode) {
326       ierr = PetscPrintf(rp->comm,"  Performance:\n"); CHKERRQ(ierr);
327     }
328     {
329       PetscReal max_error;
330       ierr = ComputeErrorMax(user_O, op_error, X, target, &max_error);
331       CHKERRQ(ierr);
332       PetscReal tol = 5e-2;
333       if (!rp->test_mode || max_error > tol) {
334         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, rp->comm);
335         CHKERRQ(ierr);
336         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, rp->comm);
337         CHKERRQ(ierr);
338         ierr = PetscPrintf(rp->comm,
339                            "    Pointwise Error (max)              : %e\n"
340                            "    CG Solve Time                      : %g (%g) sec\n",
341                            (double)max_error, rt_max, rt_min); CHKERRQ(ierr);
342       }
343     }
344     if (!rp->test_mode) {
345       ierr = PetscPrintf(rp->comm,
346                          "    DoFs/Sec in CG                     : %g (%g) million\n",
347                          1e-6*g_size*its/rt_max,
348                          1e-6*g_size*its/rt_min); CHKERRQ(ierr);
349     }
350   }
351 
352   if (rp->write_solution) {
353     PetscViewer vtk_viewer_soln;
354 
355     ierr = PetscViewerCreate(rp->comm, &vtk_viewer_soln); CHKERRQ(ierr);
356     ierr = PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK); CHKERRQ(ierr);
357     ierr = PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"); CHKERRQ(ierr);
358     ierr = VecView(X, vtk_viewer_soln); CHKERRQ(ierr);
359     ierr = PetscViewerDestroy(&vtk_viewer_soln); CHKERRQ(ierr);
360   }
361 
362   // Cleanup
363   ierr = VecDestroy(&X); CHKERRQ(ierr);
364   ierr = VecDestroy(&X_loc); CHKERRQ(ierr);
365   ierr = VecDestroy(&user_O->Y_loc); CHKERRQ(ierr);
366   ierr = MatDestroy(&mat_O); CHKERRQ(ierr);
367   ierr = PetscFree(user_O); CHKERRQ(ierr);
368   ierr = CeedDataDestroy(0, ceed_data); CHKERRQ(ierr);
369 
370   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
371   ierr = VecDestroy(&rhs_loc); CHKERRQ(ierr);
372   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
373   CeedVectorDestroy(&target);
374   CeedQFunctionDestroy(&qf_error);
375   CeedOperatorDestroy(&op_error);
376   CeedDestroy(&ceed);
377   PetscFunctionReturn(0);
378 }
379 
380 static PetscErrorCode Run(RunParams rp, PetscInt num_resources,
381                           char *const *ceed_resources, PetscInt num_bp_choices,
382                           const BPType *bp_choices) {
383   PetscInt ierr;
384   DM dm;
385 
386   PetscFunctionBeginUser;
387   // Setup DM
388   if (rp->read_mesh) {
389     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, rp->filename, NULL, PETSC_TRUE,
390                                 &dm);
391     CHKERRQ(ierr);
392   } else {
393     if (rp->user_l_nodes) {
394       // Find a nicely composite number of elements no less than global nodes
395       PetscMPIInt size;
396       ierr = MPI_Comm_size(rp->comm, &size); CHKERRQ(ierr);
397       for (PetscInt g_elem =
398              PetscMax(1, size * rp->local_nodes / PetscPowInt(rp->degree, rp->dim));
399            ;
400            g_elem++) {
401         Split3(g_elem, rp->mesh_elem, true);
402         if (Max3(rp->mesh_elem) / Min3(rp->mesh_elem) <= 2) break;
403       }
404     }
405     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, rp->dim, PETSC_FALSE,
406                                rp->mesh_elem,
407                                NULL, NULL, NULL, PETSC_TRUE, &dm); CHKERRQ(ierr);
408   }
409 
410   {
411     DM dm_dist = NULL;
412     PetscPartitioner part;
413 
414     ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
415     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
416     ierr = DMPlexDistribute(dm, 0, NULL, &dm_dist); CHKERRQ(ierr);
417     if (dm_dist) {
418       ierr = DMDestroy(&dm); CHKERRQ(ierr);
419       dm  = dm_dist;
420     }
421   }
422   // Disable default VECSTANDARD *after* distribution (which creates a Vec)
423   ierr = DMSetVecType(dm, NULL); CHKERRQ(ierr);
424 
425   for (PetscInt b = 0; b < num_bp_choices; b++) {
426     DM dm_deg;
427     VecType vec_type;
428     PetscInt q_extra = rp->q_extra;
429     rp->bp_choice = bp_choices[b];
430     rp->num_comp_u = bp_options[rp->bp_choice].num_comp_u;
431     rp->q_extra = q_extra < 0 ? bp_options[rp->bp_choice].q_extra : q_extra;
432     ierr = DMClone(dm, &dm_deg); CHKERRQ(ierr);
433     ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr);
434     ierr = DMSetVecType(dm_deg, vec_type); CHKERRQ(ierr);
435     // Create DM
436     PetscInt dim;
437     ierr = DMGetDimension(dm_deg, &dim); CHKERRQ(ierr);
438     ierr = SetupDMByDegree(dm_deg, rp->degree, rp->num_comp_u, dim,
439                            bp_options[rp->bp_choice].enforce_bc,
440                            bp_options[rp->bp_choice].bc_func); CHKERRQ(ierr);
441     for (PetscInt r = 0; r < num_resources; r++) {
442       ierr = RunWithDM(rp, dm_deg, ceed_resources[r]); CHKERRQ(ierr);
443     }
444     ierr = DMDestroy(&dm_deg); CHKERRQ(ierr);
445     rp->q_extra = q_extra;
446   }
447 
448   ierr = DMDestroy(&dm); CHKERRQ(ierr);
449   PetscFunctionReturn(0);
450 }
451 
452 int main(int argc, char **argv) {
453   PetscInt ierr, comm_size;
454   RunParams rp;
455   MPI_Comm comm;
456   char filename[PETSC_MAX_PATH_LEN];
457   char *ceed_resources[30];
458   PetscInt num_ceed_resources = 30;
459   char hostname[PETSC_MAX_PATH_LEN];
460 
461   PetscInt dim = 3, mesh_elem[3] = {3, 3, 3};
462   PetscInt num_degrees = 30, degree[30] = {}, num_local_nodes = 2,
463                                           local_nodes[2] = {};
464   PetscMPIInt ranks_per_node;
465   PetscBool degree_set;
466   BPType bp_choices[10];
467   PetscInt num_bp_choices = 10;
468 
469   // Initialize PETSc
470   ierr = PetscInitialize(&argc, &argv, NULL, help);
471   if (ierr) return ierr;
472   comm = PETSC_COMM_WORLD;
473   ierr = MPI_Comm_size(comm, &comm_size);
474   if (ierr != MPI_SUCCESS) return ierr;
475   #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
476   {
477     MPI_Comm splitcomm;
478     ierr = MPI_Comm_split_type(comm, MPI_COMM_TYPE_SHARED, 0, MPI_INFO_NULL,
479                                &splitcomm);
480     CHKERRQ(ierr);
481     ierr = MPI_Comm_size(splitcomm, &ranks_per_node); CHKERRQ(ierr);
482     ierr = MPI_Comm_free(&splitcomm); CHKERRQ(ierr);
483   }
484   #else
485   ranks_per_node = -1; // Unknown
486   #endif
487 
488   // Setup all parameters needed in Run()
489   ierr = PetscMalloc1(1, &rp); CHKERRQ(ierr);
490   rp->comm = comm;
491 
492   // Read command line options
493   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
494   CHKERRQ(ierr);
495   {
496     PetscBool set;
497     ierr = PetscOptionsEnumArray("-problem", "CEED benchmark problem to solve",
498                                  NULL,
499                                  bp_types, (PetscEnum *)bp_choices, &num_bp_choices, &set);
500     CHKERRQ(ierr);
501     if (!set) {
502       bp_choices[0] = CEED_BP1;
503       num_bp_choices = 1;
504     }
505   }
506   rp->test_mode = PETSC_FALSE;
507   ierr = PetscOptionsBool("-test",
508                           "Testing mode (do not print unless error is large)",
509                           NULL, rp->test_mode, &rp->test_mode, NULL); CHKERRQ(ierr);
510   rp->write_solution = PETSC_FALSE;
511   ierr = PetscOptionsBool("-write_solution", "Write solution for visualization",
512                           NULL, rp->write_solution, &rp->write_solution, NULL);
513   CHKERRQ(ierr);
514   degree[0] = rp->test_mode ? 3 : 2;
515   ierr = PetscOptionsIntArray("-degree",
516                               "Polynomial degree of tensor product basis", NULL,
517                               degree, &num_degrees, &degree_set); CHKERRQ(ierr);
518   if (!degree_set)
519     num_degrees = 1;
520   rp->q_extra = PETSC_DECIDE;
521   ierr = PetscOptionsInt("-q_extra",
522                          "Number of extra quadrature points (-1 for auto)", NULL,
523                          rp->q_extra, &rp->q_extra, NULL); CHKERRQ(ierr);
524   {
525     PetscBool set;
526     ierr = PetscOptionsStringArray("-ceed",
527                                    "CEED resource specifier (comma-separated list)", NULL,
528                                    ceed_resources, &num_ceed_resources, &set); CHKERRQ(ierr);
529     if (!set) {
530       ierr = PetscStrallocpy( "/cpu/self", &ceed_resources[0]); CHKERRQ(ierr);
531       num_ceed_resources = 1;
532     }
533   }
534   ierr = PetscGetHostName(hostname, sizeof hostname); CHKERRQ(ierr);
535   ierr = PetscOptionsString("-hostname", "Hostname for output", NULL, hostname,
536                             hostname, sizeof(hostname), NULL); CHKERRQ(ierr);
537   rp->read_mesh = PETSC_FALSE;
538   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL, filename,
539                             filename, sizeof(filename), &rp->read_mesh);
540   CHKERRQ(ierr);
541   rp->filename = filename;
542   if (!rp->read_mesh) {
543     PetscInt tmp = dim;
544     ierr = PetscOptionsIntArray("-cells", "Number of cells per dimension", NULL,
545                                 mesh_elem, &tmp, NULL); CHKERRQ(ierr);
546   }
547   local_nodes[0] = 1000;
548   ierr = PetscOptionsIntArray("-local_nodes",
549                               "Target number of locally owned nodes per "
550                               "process (single value or min,max)",
551                               NULL, local_nodes, &num_local_nodes, &rp->user_l_nodes);
552   CHKERRQ(ierr);
553   if (num_local_nodes < 2)
554     local_nodes[1] = 2 * local_nodes[0];
555   {
556     PetscInt two = 2;
557     rp->ksp_max_it_clip[0] = 5;
558     rp->ksp_max_it_clip[1] = 20;
559     ierr = PetscOptionsIntArray("-ksp_max_it_clip",
560                                 "Min and max number of iterations to use during benchmarking",
561                                 NULL, rp->ksp_max_it_clip, &two, NULL); CHKERRQ(ierr);
562   }
563   if (!degree_set) {
564     PetscInt max_degree = 8;
565     ierr = PetscOptionsInt("-max_degree",
566                            "Range of degrees [1, max_degree] to run with",
567                            NULL, max_degree, &max_degree, NULL);
568     CHKERRQ(ierr);
569     for (PetscInt i = 0; i < max_degree; i++)
570       degree[i] = i + 1;
571     num_degrees = max_degree;
572   }
573   {
574     PetscBool flg;
575     PetscInt p = ranks_per_node;
576     ierr = PetscOptionsInt("-p", "Number of MPI ranks per node", NULL,
577                            p, &p, &flg);
578     CHKERRQ(ierr);
579     if (flg) ranks_per_node = p;
580   }
581 
582   ierr = PetscOptionsEnd();
583   CHKERRQ(ierr);
584 
585   // Register PETSc logging stage
586   ierr = PetscLogStageRegister("Solve Stage", &rp->solve_stage);
587   CHKERRQ(ierr);
588 
589   rp->hostname = hostname;
590   rp->dim = dim;
591   rp->mesh_elem = mesh_elem;
592   rp->ranks_per_node = ranks_per_node;
593 
594   for (PetscInt d = 0; d < num_degrees; d++) {
595     PetscInt deg = degree[d];
596     for (PetscInt n = local_nodes[0]; n < local_nodes[1]; n *= 2) {
597       rp->degree = deg;
598       rp->local_nodes = n;
599       ierr = Run(rp, num_ceed_resources, ceed_resources,
600                  num_bp_choices, bp_choices); CHKERRQ(ierr);
601     }
602   }
603   // Clear memory
604   ierr = PetscFree(rp); CHKERRQ(ierr);
605   for (PetscInt i=0; i<num_ceed_resources; i++) {
606     ierr = PetscFree(ceed_resources[i]); CHKERRQ(ierr);
607   }
608   return PetscFinalize();
609 }
610