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