xref: /libCEED/examples/petsc/bps.c (revision e83e87a50f1b2e8810b376a735430565127e4d25)
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 <ceed.h>
45 #include <petscdmplex.h>
46 #include <petscksp.h>
47 #include <stdbool.h>
48 #include <string.h>
49 #include "bps.h"
50 
51 // -----------------------------------------------------------------------------
52 // Utilities
53 // -----------------------------------------------------------------------------
54 
55 // Utility function, compute three factors of an integer
56 static void Split3(PetscInt size, PetscInt m[3], bool reverse) {
57   for (PetscInt d=0,sizeleft=size; d<3; d++) {
58     PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(sizeleft, 1./(3 - d)));
59     while (try * (sizeleft / try) != sizeleft) try++;
60     m[reverse ? 2-d : d] = try;
61     sizeleft /= try;
62   }
63 }
64 
65 static int Max3(const PetscInt a[3]) {
66   return PetscMax(a[0], PetscMax(a[1], a[2]));
67 }
68 
69 static int Min3(const PetscInt a[3]) {
70   return PetscMin(a[0], PetscMin(a[1], a[2]));
71 }
72 
73 // -----------------------------------------------------------------------------
74 // Parameter structure for running problems
75 // -----------------------------------------------------------------------------
76 typedef struct RunParams_ *RunParams;
77 struct RunParams_ {
78   MPI_Comm comm;
79   PetscBool test_mode, read_mesh, userlnodes, write_solution;
80   char *filename, *hostname;
81   PetscInt localnodes, degree, qextra, dim, ncompu, *melem;
82   PetscInt ksp_max_it_clip[2];
83   PetscMPIInt rankspernode;
84   bpType bpchoice;
85   PetscLogStage solvestage;
86 };
87 
88 // -----------------------------------------------------------------------------
89 // Main body of program, called in a loop for performance benchmarking purposes
90 // -----------------------------------------------------------------------------
91 static PetscErrorCode RunWithDM(RunParams rp, DM dm,
92                                 const char *ceedresource) {
93   PetscErrorCode ierr;
94   double my_rt_start, my_rt, rt_min, rt_max;
95   PetscInt xlsize, lsize, gsize;
96   PetscScalar *r;
97   Vec X, Xloc, rhs, rhsloc;
98   Mat matO;
99   KSP ksp;
100   UserO userO;
101   Ceed ceed;
102   CeedData ceeddata;
103   CeedQFunction qfError;
104   CeedOperator opError;
105   CeedVector rhsceed, target;
106   VecType vectype;
107   PetscMemType memtype;
108 
109   PetscFunctionBeginUser;
110   // Set up libCEED
111   CeedInit(ceedresource, &ceed);
112   CeedMemType memtypebackend;
113   CeedGetPreferredMemType(ceed, &memtypebackend);
114 
115   ierr = DMGetVecType(dm, &vectype); CHKERRQ(ierr);
116   if (!vectype) { // Not yet set by user -dm_vec_type
117     switch (memtypebackend) {
118     case CEED_MEM_HOST: vectype = VECSTANDARD; break;
119     case CEED_MEM_DEVICE: {
120       const char *resolved;
121       CeedGetResource(ceed, &resolved);
122       if (strstr(resolved, "/gpu/cuda")) vectype = VECCUDA;
123       else if (strstr(resolved, "/gpu/hip/occa"))
124         vectype = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678
125       else if (strstr(resolved, "/gpu/hip")) vectype = VECHIP;
126       else vectype = VECSTANDARD;
127     }
128     }
129     ierr = DMSetVecType(dm, vectype); CHKERRQ(ierr);
130   }
131 
132   // Create global and local solution vectors
133   ierr = DMCreateGlobalVector(dm, &X); CHKERRQ(ierr);
134   ierr = VecGetLocalSize(X, &lsize); CHKERRQ(ierr);
135   ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
136   ierr = DMCreateLocalVector(dm, &Xloc); CHKERRQ(ierr);
137   ierr = VecGetSize(Xloc, &xlsize); CHKERRQ(ierr);
138   ierr = VecDuplicate(X, &rhs); CHKERRQ(ierr);
139 
140   // Operator
141   ierr = PetscMalloc1(1, &userO); CHKERRQ(ierr);
142   ierr = MatCreateShell(rp->comm, lsize, lsize, gsize, gsize,
143                         userO, &matO); CHKERRQ(ierr);
144   ierr = MatShellSetOperation(matO, MATOP_MULT,
145                               (void(*)(void))MatMult_Ceed); CHKERRQ(ierr);
146   ierr = MatShellSetOperation(matO, MATOP_GET_DIAGONAL,
147                               (void(*)(void))MatGetDiag); CHKERRQ(ierr);
148   ierr = MatShellSetVecType(matO, vectype); CHKERRQ(ierr);
149 
150   // Print summary
151   if (!rp->test_mode) {
152     PetscInt P = rp->degree + 1, Q = P + rp->qextra;
153 
154     const char *usedresource;
155     CeedGetResource(ceed, &usedresource);
156 
157     VecType vectype;
158     ierr = VecGetType(X, &vectype); CHKERRQ(ierr);
159 
160     PetscInt cStart, cEnd;
161     ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); CHKERRQ(ierr);
162     PetscMPIInt comm_size;
163     ierr = MPI_Comm_size(rp->comm, &comm_size); CHKERRQ(ierr);
164     ierr = PetscPrintf(rp->comm,
165                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n"
166                        "  MPI:\n"
167                        "    Hostname                           : %s\n"
168                        "    Total ranks                        : %d\n"
169                        "    Ranks per compute node             : %d\n"
170                        "  PETSc:\n"
171                        "    PETSc Vec Type                     : %s\n"
172                        "  libCEED:\n"
173                        "    libCEED Backend                    : %s\n"
174                        "    libCEED Backend MemType            : %s\n"
175                        "  Mesh:\n"
176                        "    Number of 1D Basis Nodes (P)       : %d\n"
177                        "    Number of 1D Quadrature Points (Q) : %d\n"
178                        "    Global nodes                       : %D\n"
179                        "    Local Elements                     : %D\n"
180                        "    Owned nodes                        : %D\n"
181                        "    DoF per node                       : %D\n",
182                        rp->bpchoice+1,
183                        rp->hostname,
184                        comm_size, rp->rankspernode,
185                        vectype, usedresource,
186                        CeedMemTypes[memtypebackend],
187                        P, Q, gsize/rp->ncompu, cEnd - cStart, lsize/rp->ncompu,
188                        rp->ncompu);
189     CHKERRQ(ierr);
190   }
191 
192   // Create RHS vector
193   ierr = VecDuplicate(Xloc, &rhsloc); CHKERRQ(ierr);
194   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
195   ierr = VecGetArrayAndMemType(rhsloc, &r, &memtype); CHKERRQ(ierr);
196   CeedVectorCreate(ceed, xlsize, &rhsceed);
197   CeedVectorSetArray(rhsceed, MemTypeP2C(memtype), CEED_USE_POINTER, r);
198 
199   ierr = PetscMalloc1(1, &ceeddata); CHKERRQ(ierr);
200   ierr = SetupLibceedByDegree(dm, ceed, rp->degree, rp->dim, rp->qextra,
201                               rp->dim, rp->ncompu, gsize, xlsize, bpOptions[rp->bpchoice],
202                               ceeddata, true, rhsceed, &target); CHKERRQ(ierr);
203 
204   // Gather RHS
205   CeedVectorTakeArray(rhsceed, MemTypeP2C(memtype), NULL);
206   ierr = VecRestoreArrayAndMemType(rhsloc, &r); CHKERRQ(ierr);
207   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
208   ierr = DMLocalToGlobal(dm, rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr);
209   CeedVectorDestroy(&rhsceed);
210 
211   // Create the error QFunction
212   CeedQFunctionCreateInterior(ceed, 1, bpOptions[rp->bpchoice].error,
213                               bpOptions[rp->bpchoice].errorfname, &qfError);
214   CeedQFunctionAddInput(qfError, "u", rp->ncompu, CEED_EVAL_INTERP);
215   CeedQFunctionAddInput(qfError, "true_soln", rp->ncompu, CEED_EVAL_NONE);
216   CeedQFunctionAddOutput(qfError, "error", rp->ncompu, CEED_EVAL_NONE);
217 
218   // Create the error operator
219   CeedOperatorCreate(ceed, qfError, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
220                      &opError);
221   CeedOperatorSetField(opError, "u", ceeddata->Erestrictu,
222                        ceeddata->basisu, CEED_VECTOR_ACTIVE);
223   CeedOperatorSetField(opError, "true_soln", ceeddata->Erestrictui,
224                        CEED_BASIS_COLLOCATED, target);
225   CeedOperatorSetField(opError, "error", ceeddata->Erestrictui,
226                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
227 
228   // Set up Mat
229   userO->comm = rp->comm;
230   userO->dm = dm;
231   userO->Xloc = Xloc;
232   ierr = VecDuplicate(Xloc, &userO->Yloc); CHKERRQ(ierr);
233   userO->Xceed = ceeddata->Xceed;
234   userO->Yceed = ceeddata->Yceed;
235   userO->op = ceeddata->opApply;
236   userO->ceed = ceed;
237 
238   ierr = KSPCreate(rp->comm, &ksp); CHKERRQ(ierr);
239   {
240     PC pc;
241     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
242     if (rp->bpchoice == CEED_BP1 || rp->bpchoice == CEED_BP2) {
243       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
244       ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr);
245     } else {
246       ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr);
247     }
248     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
249     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
250     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
251                             PETSC_DEFAULT); CHKERRQ(ierr);
252   }
253   ierr = KSPSetOperators(ksp, matO, matO); CHKERRQ(ierr);
254 
255   // First run's performance log is not considered for benchmarking purposes
256   ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
257   CHKERRQ(ierr);
258   my_rt_start = MPI_Wtime();
259   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
260   my_rt = MPI_Wtime() - my_rt_start;
261   ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, rp->comm);
262   CHKERRQ(ierr);
263   // Set maxits based on first iteration timing
264   if (my_rt > 0.02) {
265     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
266                             rp->ksp_max_it_clip[0]);
267     CHKERRQ(ierr);
268   } else {
269     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
270                             rp->ksp_max_it_clip[1]);
271     CHKERRQ(ierr);
272   }
273   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
274 
275   // Timed solve
276   ierr = VecZeroEntries(X); CHKERRQ(ierr);
277   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
278 
279   // -- Performance logging
280   ierr = PetscLogStagePush(rp->solvestage); CHKERRQ(ierr);
281 
282   // -- Solve
283   my_rt_start = MPI_Wtime();
284   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
285   my_rt = MPI_Wtime() - my_rt_start;
286 
287   // -- Performance logging
288   ierr = PetscLogStagePop();
289 
290   // Output results
291   {
292     KSPType ksptype;
293     KSPConvergedReason reason;
294     PetscReal rnorm;
295     PetscInt its;
296     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
297     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
298     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
299     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
300     if (!rp->test_mode || reason < 0 || rnorm > 1e-8) {
301       ierr = PetscPrintf(rp->comm,
302                          "  KSP:\n"
303                          "    KSP Type                           : %s\n"
304                          "    KSP Convergence                    : %s\n"
305                          "    Total KSP Iterations               : %D\n"
306                          "    Final rnorm                        : %e\n",
307                          ksptype, KSPConvergedReasons[reason], its,
308                          (double)rnorm); CHKERRQ(ierr);
309     }
310     if (!rp->test_mode) {
311       ierr = PetscPrintf(rp->comm,"  Performance:\n"); CHKERRQ(ierr);
312     }
313     {
314       PetscReal maxerror;
315       ierr = ComputeErrorMax(userO, opError, X, target, &maxerror);
316       CHKERRQ(ierr);
317       PetscReal tol = 5e-2;
318       if (!rp->test_mode || maxerror > tol) {
319         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, rp->comm);
320         CHKERRQ(ierr);
321         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, rp->comm);
322         CHKERRQ(ierr);
323         ierr = PetscPrintf(rp->comm,
324                            "    Pointwise Error (max)              : %e\n"
325                            "    CG Solve Time                      : %g (%g) sec\n",
326                            (double)maxerror, rt_max, rt_min); CHKERRQ(ierr);
327       }
328     }
329     if (!rp->test_mode) {
330       ierr = PetscPrintf(rp->comm,
331                          "    DoFs/Sec in CG                     : %g (%g) million\n",
332                          1e-6*gsize*its/rt_max,
333                          1e-6*gsize*its/rt_min); CHKERRQ(ierr);
334     }
335   }
336 
337   if (rp->write_solution) {
338     PetscViewer vtkviewersoln;
339 
340     ierr = PetscViewerCreate(rp->comm, &vtkviewersoln); CHKERRQ(ierr);
341     ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr);
342     ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtu"); CHKERRQ(ierr);
343     ierr = VecView(X, vtkviewersoln); CHKERRQ(ierr);
344     ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr);
345   }
346 
347   // Cleanup
348   ierr = VecDestroy(&X); CHKERRQ(ierr);
349   ierr = VecDestroy(&Xloc); CHKERRQ(ierr);
350   ierr = VecDestroy(&userO->Yloc); CHKERRQ(ierr);
351   ierr = MatDestroy(&matO); CHKERRQ(ierr);
352   ierr = PetscFree(userO); CHKERRQ(ierr);
353   ierr = CeedDataDestroy(0, ceeddata); CHKERRQ(ierr);
354 
355   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
356   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
357   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
358   CeedVectorDestroy(&target);
359   CeedQFunctionDestroy(&qfError);
360   CeedOperatorDestroy(&opError);
361   CeedDestroy(&ceed);
362   PetscFunctionReturn(0);
363 }
364 
365 static PetscErrorCode Run(RunParams rp, PetscInt num_resources,
366                           char *const *ceedresources, PetscInt num_bpchoices,
367                           const bpType *bpchoices) {
368   PetscInt ierr;
369   DM  dm;
370 
371   PetscFunctionBeginUser;
372   // Setup DM
373   if (rp->read_mesh) {
374     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, rp->filename, PETSC_TRUE, &dm);
375     CHKERRQ(ierr);
376   } else {
377     if (rp->userlnodes) {
378       // Find a nicely composite number of elements no less than global nodes
379       PetscMPIInt size;
380       ierr = MPI_Comm_size(rp->comm, &size); CHKERRQ(ierr);
381       for (PetscInt gelem =
382              PetscMax(1, size * rp->localnodes / PetscPowInt(rp->degree, rp->dim));
383            ;
384            gelem++) {
385         Split3(gelem, rp->melem, true);
386         if (Max3(rp->melem) / Min3(rp->melem) <= 2) break;
387       }
388     }
389     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, rp->dim, PETSC_FALSE, rp->melem,
390                                NULL, NULL, NULL, PETSC_TRUE, &dm); CHKERRQ(ierr);
391   }
392 
393   {
394     DM dmDist = NULL;
395     PetscPartitioner part;
396 
397     ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
398     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
399     ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
400     if (dmDist) {
401       ierr = DMDestroy(&dm); CHKERRQ(ierr);
402       dm  = dmDist;
403     }
404   }
405   // Disable default VECSTANDARD *after* distribution (which creates a Vec)
406   ierr = DMSetVecType(dm, NULL); CHKERRQ(ierr);
407 
408   for (PetscInt b = 0; b < num_bpchoices; b++) {
409     DM dm_deg;
410     VecType vectype;
411     PetscInt qextra = rp->qextra;
412     rp->bpchoice = bpchoices[b];
413     rp->ncompu = bpOptions[rp->bpchoice].ncompu;
414     rp->qextra = qextra < 0 ? bpOptions[rp->bpchoice].qextra : qextra;
415     ierr = DMClone(dm, &dm_deg); CHKERRQ(ierr);
416     ierr = DMGetVecType(dm, &vectype); CHKERRQ(ierr);
417     ierr = DMSetVecType(dm_deg, vectype); CHKERRQ(ierr);
418     // Create DM
419     PetscInt dim;
420     ierr = DMGetDimension(dm_deg, &dim); CHKERRQ(ierr);
421     ierr = SetupDMByDegree(dm_deg, rp->degree, rp->ncompu, dim,
422                            bpOptions[rp->bpchoice].enforcebc,
423                            bpOptions[rp->bpchoice].bcsfunc); CHKERRQ(ierr);
424     for (PetscInt r = 0; r < num_resources; r++) {
425       ierr = RunWithDM(rp, dm_deg, ceedresources[r]); CHKERRQ(ierr);
426     }
427     ierr = DMDestroy(&dm_deg); CHKERRQ(ierr);
428     rp->qextra = qextra;
429   }
430 
431   ierr = DMDestroy(&dm); CHKERRQ(ierr);
432   PetscFunctionReturn(0);
433 }
434 
435 int main(int argc, char **argv) {
436   PetscInt ierr, commsize;
437   RunParams rp;
438   MPI_Comm comm;
439   char filename[PETSC_MAX_PATH_LEN];
440   char *ceedresources[30];
441   PetscInt num_ceedresources = 30;
442   char hostname[PETSC_MAX_PATH_LEN];
443 
444   PetscInt dim = 3, melem[3] = {3, 3, 3};
445   PetscInt num_degrees = 30, degree[30] = {}, num_localnodes = 2, localnodes[2] = {};
446   PetscMPIInt rankspernode;
447   PetscBool degree_set;
448   bpType bpchoices[10];
449   PetscInt num_bpchoices = 10;
450 
451   // Initialize PETSc
452   ierr = PetscInitialize(&argc, &argv, NULL, help);
453   if (ierr) return ierr;
454   comm = PETSC_COMM_WORLD;
455   ierr = MPI_Comm_size(comm, &commsize);
456   if (ierr != MPI_SUCCESS) return ierr;
457   #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
458   {
459     MPI_Comm splitcomm;
460     ierr = MPI_Comm_split_type(comm, MPI_COMM_TYPE_SHARED, 0, MPI_INFO_NULL,
461                                &splitcomm);
462     CHKERRQ(ierr);
463     ierr = MPI_Comm_size(splitcomm, &rankspernode); CHKERRQ(ierr);
464     ierr = MPI_Comm_free(&splitcomm); CHKERRQ(ierr);
465   }
466   #else
467   rankspernode = -1; // Unknown
468   #endif
469 
470   // Setup all parameters needed in Run()
471   ierr = PetscMalloc1(1, &rp); CHKERRQ(ierr);
472   rp->comm = comm;
473 
474   // Read command line options
475   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
476   CHKERRQ(ierr);
477   {
478     PetscBool set;
479     ierr = PetscOptionsEnumArray("-problem", "CEED benchmark problem to solve",
480                                  NULL,
481                                  bpTypes, (PetscEnum *)bpchoices, &num_bpchoices, &set);
482     CHKERRQ(ierr);
483     if (!set) {
484       bpchoices[0] = CEED_BP1;
485       num_bpchoices = 1;
486     }
487   }
488   rp->test_mode = PETSC_FALSE;
489   ierr = PetscOptionsBool("-test",
490                           "Testing mode (do not print unless error is large)",
491                           NULL, rp->test_mode, &rp->test_mode, NULL); CHKERRQ(ierr);
492   rp->write_solution = PETSC_FALSE;
493   ierr = PetscOptionsBool("-write_solution", "Write solution for visualization",
494                           NULL, rp->write_solution, &rp->write_solution, NULL);
495   CHKERRQ(ierr);
496   degree[0] = rp->test_mode ? 3 : 2;
497   ierr = PetscOptionsIntArray("-degree",
498                               "Polynomial degree of tensor product basis", NULL,
499                               degree, &num_degrees, &degree_set); CHKERRQ(ierr);
500   if (!degree_set)
501     num_degrees = 1;
502   rp->qextra = PETSC_DECIDE;
503   ierr = PetscOptionsInt("-qextra",
504                          "Number of extra quadrature points (-1 for auto)", NULL,
505                          rp->qextra, &rp->qextra, NULL); CHKERRQ(ierr);
506   {
507     PetscBool set;
508     ierr = PetscOptionsStringArray("-ceed",
509                                    "CEED resource specifier (comma-separated list)", NULL,
510                                    ceedresources, &num_ceedresources, &set); CHKERRQ(ierr);
511     if (!set) {
512       ierr = PetscStrallocpy( "/cpu/self", &ceedresources[0]); CHKERRQ(ierr);
513       num_ceedresources = 1;
514     }
515   }
516   ierr = PetscGetHostName(hostname, sizeof hostname); CHKERRQ(ierr);
517   ierr = PetscOptionsString("-hostname", "Hostname for output", NULL, hostname,
518                             hostname, sizeof(hostname), NULL); CHKERRQ(ierr);
519   rp->read_mesh = PETSC_FALSE;
520   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL, filename,
521                             filename, sizeof(filename), &rp->read_mesh);
522   CHKERRQ(ierr);
523   rp->filename = filename;
524   if (!rp->read_mesh) {
525     PetscInt tmp = dim;
526     ierr = PetscOptionsIntArray("-cells", "Number of cells per dimension", NULL,
527                                 melem, &tmp, NULL); CHKERRQ(ierr);
528   }
529   localnodes[0] = 1000;
530   ierr = PetscOptionsIntArray("-local_nodes",
531                               "Target number of locally owned nodes per "
532                               "process (single value or min,max)",
533                               NULL, localnodes, &num_localnodes, &rp->userlnodes);
534   CHKERRQ(ierr);
535   if (num_localnodes < 2)
536     localnodes[1] = 2 * localnodes[0];
537   {
538     PetscInt two = 2;
539     rp->ksp_max_it_clip[0] = 5;
540     rp->ksp_max_it_clip[1] = 20;
541     ierr = PetscOptionsIntArray("-ksp_max_it_clip",
542                                 "Min and max number of iterations to use during benchmarking",
543                                 NULL, rp->ksp_max_it_clip, &two, NULL); CHKERRQ(ierr);
544   }
545   if (!degree_set) {
546     PetscInt maxdegree = 8;
547     ierr = PetscOptionsInt("-max_degree",
548                            "Range of degrees [1, maxdegree] to run with",
549                            NULL, maxdegree, &maxdegree, NULL);
550     CHKERRQ(ierr);
551     for (PetscInt i = 0; i < maxdegree; i++)
552       degree[i] = i + 1;
553     num_degrees = maxdegree;
554   }
555   {
556     PetscBool flg;
557     PetscInt p = rankspernode;
558     ierr = PetscOptionsInt("-p", "Number of MPI ranks per node", NULL,
559                            p, &p, &flg);
560     CHKERRQ(ierr);
561     if (flg) rankspernode = p;
562   }
563 
564   ierr = PetscOptionsEnd();
565   CHKERRQ(ierr);
566 
567   // Register PETSc logging stage
568   ierr = PetscLogStageRegister("Solve Stage", &rp->solvestage);
569   CHKERRQ(ierr);
570 
571   rp->hostname = hostname;
572   rp->dim = dim;
573   rp->melem = melem;
574   rp->rankspernode = rankspernode;
575 
576   for (PetscInt d = 0; d < num_degrees; d++) {
577     PetscInt deg = degree[d];
578     for (PetscInt n = localnodes[0]; n < localnodes[1]; n *= 2) {
579       rp->degree = deg;
580       rp->localnodes = n;
581       ierr = Run(rp, num_ceedresources, ceedresources,
582                  num_bpchoices, bpchoices); CHKERRQ(ierr);
583     }
584   }
585   // Clear memory
586   ierr = PetscFree(rp); CHKERRQ(ierr);
587   for (PetscInt i=0; i<num_ceedresources; i++) {
588     ierr = PetscFree(ceedresources[i]); CHKERRQ(ierr);
589   }
590   return PetscFinalize();
591 }
592