xref: /libCEED/examples/petsc/bpsraw.c (revision 36d2312c6a9f418dd1dc278db7b470a7527039c0)
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 is intentionally "raw", using only low-level communication
23 // primitives.
24 //
25 // Build with:
26 //
27 //     make bpsraw [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
28 //
29 // Sample runs:
30 //
31 //     ./bpsraw -problem bp1
32 //     ./bpsraw -problem bp2 -ceed /cpu/self
33 //     ./bpsraw -problem bp3 -ceed /gpu/occa
34 //     ./bpsraw -problem bp4 -ceed /cpu/occa
35 //     ./bpsraw -problem bp5 -ceed /omp/occa
36 //     ./bpsraw -problem bp6 -ceed /ocl/occa
37 //
38 //TESTARGS -ceed {ceed_resource} -test -problem bp2 -degree 5 -qextra 5
39 
40 /// @file
41 /// CEED BPs example using PETSc
42 /// See bps.c for an implementation using DMPlex unstructured grids.
43 const char help[] = "Solve CEED BPs using PETSc\n";
44 
45 #include <stdbool.h>
46 #include <string.h>
47 #include <petscksp.h>
48 #include <ceed.h>
49 #include "qfunctions/bps/common.h"
50 #include "qfunctions/bps/bp1.h"
51 #include "qfunctions/bps/bp2.h"
52 #include "qfunctions/bps/bp3.h"
53 #include "qfunctions/bps/bp4.h"
54 
55 #include <petscsys.h>
56 #if PETSC_VERSION_LT(3,12,0)
57 #ifdef PETSC_HAVE_CUDA
58 #include <petsccuda.h>
59 // Note: With PETSc prior to version 3.12.0, providing the source path to
60 //       include 'cublas_v2.h' will be needed to use 'petsccuda.h'.
61 #endif
62 #endif
63 
64 static void Split3(PetscInt size, PetscInt m[3], bool reverse) {
65   for (PetscInt d=0,sizeleft=size; d<3; d++) {
66     PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(sizeleft, 1./(3 - d)));
67     while (try * (sizeleft / try) != sizeleft) try++;
68     m[reverse ? 2-d : d] = try;
69     sizeleft /= try;
70   }
71 }
72 
73 static PetscInt Max3(const PetscInt a[3]) {
74   return PetscMax(a[0], PetscMax(a[1], a[2]));
75 }
76 static PetscInt Min3(const PetscInt a[3]) {
77   return PetscMin(a[0], PetscMin(a[1], a[2]));
78 }
79 static void GlobalNodes(const PetscInt p[3], const PetscInt irank[3],
80                         PetscInt degree, const PetscInt melem[3],
81                         PetscInt mnodes[3]) {
82   for (int d=0; d<3; d++)
83     mnodes[d] = degree*melem[d] + (irank[d] == p[d]-1);
84 }
85 static PetscInt GlobalStart(const PetscInt p[3], const PetscInt irank[3],
86                             PetscInt degree, const PetscInt melem[3]) {
87   PetscInt start = 0;
88   // Dumb brute-force is easier to read
89   for (PetscInt i=0; i<p[0]; i++) {
90     for (PetscInt j=0; j<p[1]; j++) {
91       for (PetscInt k=0; k<p[2]; k++) {
92         PetscInt mnodes[3], ijkrank[] = {i,j,k};
93         if (i == irank[0] && j == irank[1] && k == irank[2]) return start;
94         GlobalNodes(p, ijkrank, degree, melem, mnodes);
95         start += mnodes[0] * mnodes[1] * mnodes[2];
96       }
97     }
98   }
99   return -1;
100 }
101 static int CreateRestriction(Ceed ceed, CeedInterlaceMode imode,
102                              const CeedInt melem[3], CeedInt P, CeedInt ncomp,
103                              CeedElemRestriction *Erestrict) {
104   const PetscInt nelem = melem[0]*melem[1]*melem[2];
105   PetscInt mnodes[3], *idx, *idxp;
106 
107   // Get indicies
108   for (int d=0; d<3; d++) mnodes[d] = melem[d]*(P-1) + 1;
109   idxp = idx = malloc(nelem*P*P*P*sizeof idx[0]);
110   for (CeedInt i=0; i<melem[0]; i++)
111     for (CeedInt j=0; j<melem[1]; j++)
112       for (CeedInt k=0; k<melem[2]; k++,idxp += P*P*P)
113         for (CeedInt ii=0; ii<P; ii++)
114           for (CeedInt jj=0; jj<P; jj++)
115             for (CeedInt kk=0; kk<P; kk++) {
116               if (0) { // This is the C-style (i,j,k) ordering that I prefer
117                 idxp[(ii*P+jj)*P+kk] = (((i*(P-1)+ii)*mnodes[1]
118                                          + (j*(P-1)+jj))*mnodes[2]
119                                         + (k*(P-1)+kk));
120               } else { // (k,j,i) ordering for consistency with MFEM example
121                 idxp[ii+P*(jj+P*kk)] = (((i*(P-1)+ii)*mnodes[1]
122                                          + (j*(P-1)+jj))*mnodes[2]
123                                         + (k*(P-1)+kk));
124               }
125             }
126 
127   // Setup CEED restriction
128   CeedElemRestrictionCreate(ceed, imode, nelem, P*P*P,
129                             mnodes[0]*mnodes[1]*mnodes[2], ncomp,
130                             CEED_MEM_HOST, CEED_OWN_POINTER, idx, Erestrict);
131 
132   PetscFunctionReturn(0);
133 }
134 
135 // Data for PETSc
136 typedef struct User_ *User;
137 struct User_ {
138   MPI_Comm comm;
139   VecScatter ltog;              // Scatter for all entries
140   VecScatter ltog0;             // Skip Dirichlet values
141   VecScatter gtogD;             // global-to-global; only Dirichlet values
142   Vec Xloc, Yloc;
143   CeedVector xceed, yceed;
144   CeedOperator op;
145   CeedVector qdata;
146   Ceed ceed;
147   CeedMemType memtype;
148   int (*VecGetArray)(Vec, PetscScalar **);
149   int (*VecGetArrayRead)(Vec, const PetscScalar **);
150   int (*VecRestoreArray)(Vec, PetscScalar **);
151   int (*VecRestoreArrayRead)(Vec, const PetscScalar **);
152 };
153 
154 // MemType Options
155 static const char *const memTypes[] = {"host","device","memType",
156                                        "CEED_MEM_",0
157                                       };
158 
159 // BP Options
160 typedef enum {
161   CEED_BP1 = 0, CEED_BP2 = 1, CEED_BP3 = 2,
162   CEED_BP4 = 3, CEED_BP5 = 4, CEED_BP6 = 5
163 } bpType;
164 static const char *const bpTypes[] = {"bp1","bp2","bp3","bp4","bp5","bp6",
165                                       "bpType","CEED_BP",0
166                                      };
167 
168 // BP specific data
169 typedef struct {
170   CeedInt ncompu, qdatasize, qextra;
171   CeedQFunctionUser setupgeo, setuprhs, apply, error;
172   const char *setupgeofname, *setuprhsfname, *applyfname, *errorfname;
173   CeedEvalMode inmode, outmode;
174   CeedQuadMode qmode;
175 } bpData;
176 
177 bpData bpOptions[6] = {
178   [CEED_BP1] = {
179     .ncompu = 1,
180     .qdatasize = 1,
181     .qextra = 1,
182     .setupgeo = SetupMassGeo,
183     .setuprhs = SetupMassRhs,
184     .apply = Mass,
185     .error = Error,
186     .setupgeofname = SetupMassGeo_loc,
187     .setuprhsfname = SetupMassRhs_loc,
188     .applyfname = Mass_loc,
189     .errorfname = Error_loc,
190     .inmode = CEED_EVAL_INTERP,
191     .outmode = CEED_EVAL_INTERP,
192     .qmode = CEED_GAUSS
193   },
194   [CEED_BP2] = {
195     .ncompu = 3,
196     .qdatasize = 1,
197     .qextra = 1,
198     .setupgeo = SetupMassGeo,
199     .setuprhs = SetupMassRhs3,
200     .apply = Mass3,
201     .error = Error3,
202     .setupgeofname = SetupMassGeo_loc,
203     .setuprhsfname = SetupMassRhs3_loc,
204     .applyfname = Mass3_loc,
205     .errorfname = Error3_loc,
206     .inmode = CEED_EVAL_INTERP,
207     .outmode = CEED_EVAL_INTERP,
208     .qmode = CEED_GAUSS
209   },
210   [CEED_BP3] = {
211     .ncompu = 1,
212     .qdatasize = 6,
213     .qextra = 1,
214     .setupgeo = SetupDiffGeo,
215     .setuprhs = SetupDiffRhs,
216     .apply = Diff,
217     .error = Error,
218     .setupgeofname = SetupDiffGeo_loc,
219     .setuprhsfname = SetupDiffRhs_loc,
220     .applyfname = Diff_loc,
221     .errorfname = Error_loc,
222     .inmode = CEED_EVAL_GRAD,
223     .outmode = CEED_EVAL_GRAD,
224     .qmode = CEED_GAUSS
225   },
226   [CEED_BP4] = {
227     .ncompu = 3,
228     .qdatasize = 6,
229     .qextra = 1,
230     .setupgeo = SetupDiffGeo,
231     .setuprhs = SetupDiffRhs3,
232     .apply = Diff3,
233     .error = Error3,
234     .setupgeofname = SetupDiffGeo_loc,
235     .setuprhsfname = SetupDiffRhs3_loc,
236     .applyfname = Diff_loc,
237     .errorfname = Error3_loc,
238     .inmode = CEED_EVAL_GRAD,
239     .outmode = CEED_EVAL_GRAD,
240     .qmode = CEED_GAUSS
241   },
242   [CEED_BP5] = {
243     .ncompu = 1,
244     .qdatasize = 6,
245     .qextra = 0,
246     .setupgeo = SetupDiffGeo,
247     .setuprhs = SetupDiffRhs,
248     .apply = Diff,
249     .error = Error,
250     .setupgeofname = SetupDiffGeo_loc,
251     .setuprhsfname = SetupDiffRhs_loc,
252     .applyfname = Diff_loc,
253     .errorfname = Error_loc,
254     .inmode = CEED_EVAL_GRAD,
255     .outmode = CEED_EVAL_GRAD,
256     .qmode = CEED_GAUSS_LOBATTO
257   },
258   [CEED_BP6] = {
259     .ncompu = 3,
260     .qdatasize = 6,
261     .qextra = 0,
262     .setupgeo = SetupDiffGeo,
263     .setuprhs = SetupDiffRhs3,
264     .apply = Diff3,
265     .error = Error3,
266     .setupgeofname = SetupDiffGeo_loc,
267     .setuprhsfname = SetupDiffRhs3_loc,
268     .applyfname = Diff_loc,
269     .errorfname = Error3_loc,
270     .inmode = CEED_EVAL_GRAD,
271     .outmode = CEED_EVAL_GRAD,
272     .qmode = CEED_GAUSS_LOBATTO
273   }
274 };
275 
276 // This function uses libCEED to compute the action of the mass matrix
277 static PetscErrorCode MatMult_Mass(Mat A, Vec X, Vec Y) {
278   PetscErrorCode ierr;
279   User user;
280   PetscScalar *x, *y;
281 
282   PetscFunctionBeginUser;
283 
284   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
285 
286   // Global-to-local
287   ierr = VecScatterBegin(user->ltog, X, user->Xloc, INSERT_VALUES,
288                          SCATTER_REVERSE); CHKERRQ(ierr);
289   ierr = VecScatterEnd(user->ltog, X, user->Xloc, INSERT_VALUES,
290                        SCATTER_REVERSE); CHKERRQ(ierr);
291 
292   // Setup libCEED vectors
293   ierr = user->VecGetArrayRead(user->Xloc, (const PetscScalar **)&x);
294   CHKERRQ(ierr);
295   ierr = user->VecGetArray(user->Yloc, &y); CHKERRQ(ierr);
296   CeedVectorSetArray(user->xceed, user->memtype, CEED_USE_POINTER, x);
297   CeedVectorSetArray(user->yceed, user->memtype, CEED_USE_POINTER, y);
298 
299   // Apply libCEED operator
300   CeedOperatorApply(user->op, user->xceed, user->yceed,
301                     CEED_REQUEST_IMMEDIATE);
302   CeedVectorSyncArray(user->yceed, user->memtype);
303 
304   // Restore PETSc vectors
305   ierr = user->VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x);
306   CHKERRQ(ierr);
307   ierr = user->VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr);
308 
309   // Local-to-global
310   if (Y) {
311     ierr = VecZeroEntries(Y); CHKERRQ(ierr);
312     ierr = VecScatterBegin(user->ltog, user->Yloc, Y, ADD_VALUES,
313                            SCATTER_FORWARD); CHKERRQ(ierr);
314     ierr = VecScatterEnd(user->ltog, user->Yloc, Y, ADD_VALUES,
315                          SCATTER_FORWARD); CHKERRQ(ierr);
316   }
317   PetscFunctionReturn(0);
318 }
319 
320 // This function uses libCEED to compute the action of the Laplacian with
321 // Dirichlet boundary conditions
322 static PetscErrorCode MatMult_Diff(Mat A, Vec X, Vec Y) {
323   PetscErrorCode ierr;
324   User user;
325   PetscScalar *x, *y;
326 
327   PetscFunctionBeginUser;
328 
329   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
330 
331   // Global-to-local
332   ierr = VecScatterBegin(user->ltog0, X, user->Xloc, INSERT_VALUES,
333                          SCATTER_REVERSE); CHKERRQ(ierr);
334   ierr = VecScatterEnd(user->ltog0, X, user->Xloc, INSERT_VALUES,
335                        SCATTER_REVERSE);
336   CHKERRQ(ierr);
337 
338   // Setup libCEED vectors
339   ierr = user->VecGetArrayRead(user->Xloc, (const PetscScalar **)&x);
340   CHKERRQ(ierr);
341   ierr = user->VecGetArray(user->Yloc, &y); CHKERRQ(ierr);
342   CeedVectorSetArray(user->xceed, user->memtype, CEED_USE_POINTER, x);
343   CeedVectorSetArray(user->yceed, user->memtype, CEED_USE_POINTER, y);
344 
345   // Apply libCEED operator
346   CeedOperatorApply(user->op, user->xceed, user->yceed,
347                     CEED_REQUEST_IMMEDIATE);
348   CeedVectorSyncArray(user->yceed, user->memtype);
349 
350   // Restore PETSc vectors
351   ierr = user->VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x);
352   CHKERRQ(ierr);
353   ierr = user->VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr);
354 
355   // Local-to-global
356   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
357   ierr = VecScatterBegin(user->gtogD, X, Y, INSERT_VALUES, SCATTER_FORWARD);
358   CHKERRQ(ierr);
359   ierr = VecScatterEnd(user->gtogD, X, Y, INSERT_VALUES, SCATTER_FORWARD);
360   CHKERRQ(ierr);
361   ierr = VecScatterBegin(user->ltog0, user->Yloc, Y, ADD_VALUES,
362                          SCATTER_FORWARD); CHKERRQ(ierr);
363   ierr = VecScatterEnd(user->ltog0, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD);
364   CHKERRQ(ierr);
365 
366   PetscFunctionReturn(0);
367 }
368 
369 // This function calculates the error in the final solution
370 static PetscErrorCode ComputeErrorMax(User user, CeedOperator operror, Vec X,
371                                       CeedVector target, PetscReal *maxerror) {
372   PetscErrorCode ierr;
373   PetscScalar *x;
374   CeedVector collocated_error;
375   CeedInt length;
376 
377   PetscFunctionBeginUser;
378 
379   CeedVectorGetLength(target, &length);
380   CeedVectorCreate(user->ceed, length, &collocated_error);
381 
382   // Global-to-local
383   ierr = VecScatterBegin(user->ltog, X, user->Xloc, INSERT_VALUES,
384                          SCATTER_REVERSE); CHKERRQ(ierr);
385   ierr = VecScatterEnd(user->ltog, X, user->Xloc, INSERT_VALUES,
386                        SCATTER_REVERSE); CHKERRQ(ierr);
387 
388   // Setup libCEED vector
389   ierr = user->VecGetArrayRead(user->Xloc, (const PetscScalar **)&x);
390   CHKERRQ(ierr);
391   CeedVectorSetArray(user->xceed, user->memtype, CEED_USE_POINTER, x);
392 
393   // Apply libCEED operator
394   CeedOperatorApply(operror, user->xceed, collocated_error,
395                     CEED_REQUEST_IMMEDIATE);
396 
397   // Restore PETSc vector
398   ierr = user->VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x);
399   CHKERRQ(ierr);
400 
401   // Reduce max error
402   *maxerror = 0;
403   const CeedScalar *e;
404   CeedVectorGetArrayRead(collocated_error, CEED_MEM_HOST, &e);
405   for (CeedInt i=0; i<length; i++) {
406     *maxerror = PetscMax(*maxerror, PetscAbsScalar(e[i]));
407   }
408   CeedVectorRestoreArrayRead(collocated_error, &e);
409   ierr = MPI_Allreduce(MPI_IN_PLACE, maxerror, 1, MPIU_REAL, MPIU_MAX,
410                        user->comm); CHKERRQ(ierr);
411 
412   // Cleanup
413   CeedVectorDestroy(&collocated_error);
414 
415   PetscFunctionReturn(0);
416 }
417 
418 int main(int argc, char **argv) {
419   PetscInt ierr;
420   MPI_Comm comm;
421   char ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self";
422   double my_rt_start, my_rt, rt_min, rt_max;
423   PetscInt degree, qextra, localnodes, localelem, melem[3], mnodes[3], p[3],
424            irank[3], lnodes[3], lsize, ncompu = 1;
425   PetscScalar *r;
426   PetscBool test_mode, benchmark_mode, write_solution;
427   PetscMPIInt size, rank;
428   PetscLogStage solvestage;
429   Vec X, Xloc, rhs, rhsloc;
430   Mat mat;
431   KSP ksp;
432   VecScatter ltog, ltog0, gtogD;
433   User user;
434   Ceed ceed;
435   CeedBasis basisx, basisu;
436   CeedElemRestriction Erestrictx, Erestrictu, Erestrictui, Erestrictqdi;
437   CeedQFunction qfsetupgeo, qfsetuprhs, qfapply, qferror;
438   CeedOperator opsetupgeo, opsetuprhs, opapply, operror;
439   CeedVector xcoord, qdata, rhsceed, target;
440   CeedMemType memtyperequested;
441   CeedInt P, Q;
442   const CeedInt dim = 3, ncompx = 3;
443   bpType bpchoice;
444 
445   // Check for PETSc CUDA availability
446   PetscBool petschavecuda, setmemtyperequest = PETSC_FALSE;
447   // *INDENT-OFF*
448   #ifdef PETSC_HAVE_CUDA
449   petschavecuda = PETSC_TRUE;
450   #else
451   petschavecuda = PETSC_FALSE;
452   #endif
453   // *INDENT-ON*
454 
455   ierr = PetscInitialize(&argc, &argv, NULL, help);
456   if (ierr) return ierr;
457   comm = PETSC_COMM_WORLD;
458 
459   // Read command line options
460   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr);
461   bpchoice = CEED_BP1;
462   ierr = PetscOptionsEnum("-problem",
463                           "CEED benchmark problem to solve", NULL,
464                           bpTypes, (PetscEnum)bpchoice, (PetscEnum *)&bpchoice,
465                           NULL); CHKERRQ(ierr);
466   ncompu = bpOptions[bpchoice].ncompu;
467   test_mode = PETSC_FALSE;
468   ierr = PetscOptionsBool("-test",
469                           "Testing mode (do not print unless error is large)",
470                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
471   benchmark_mode = PETSC_FALSE;
472   ierr = PetscOptionsBool("-benchmark",
473                           "Benchmarking mode (prints benchmark statistics)",
474                           NULL, benchmark_mode, &benchmark_mode, NULL);
475   CHKERRQ(ierr);
476   write_solution = PETSC_FALSE;
477   ierr = PetscOptionsBool("-write_solution",
478                           "Write solution for visualization",
479                           NULL, write_solution, &write_solution, NULL);
480   CHKERRQ(ierr);
481   degree = test_mode ? 3 : 1;
482   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
483                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
484   qextra = bpOptions[bpchoice].qextra;
485   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
486                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
487   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
488                             NULL, ceedresource, ceedresource,
489                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
490   localnodes = 1000;
491   ierr = PetscOptionsInt("-local",
492                          "Target number of locally owned nodes per process",
493                          NULL, localnodes, &localnodes, NULL); CHKERRQ(ierr);
494   memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST;
495   ierr = PetscOptionsEnum("-memtype",
496                           "CEED MemType requested", NULL,
497                           memTypes, (PetscEnum)memtyperequested,
498                           (PetscEnum *)&memtyperequested, &setmemtyperequest);
499   CHKERRQ(ierr);
500   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
501   P = degree + 1;
502   Q = P + qextra;
503 
504   // Set up libCEED
505   CeedInit(ceedresource, &ceed);
506   CeedMemType memtypebackend;
507   CeedGetPreferredMemType(ceed, &memtypebackend);
508 
509   // Check memtype compatibility
510   if (!setmemtyperequest)
511     memtyperequested = memtypebackend;
512   else if (!petschavecuda && memtyperequested == CEED_MEM_DEVICE)
513     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS,
514              "PETSc was not built with CUDA. "
515              "Requested MemType CEED_MEM_DEVICE is not supported.", NULL);
516 
517   // Determine size of process grid
518   ierr = MPI_Comm_size(comm, &size); CHKERRQ(ierr);
519   Split3(size, p, false);
520 
521   // Find a nicely composite number of elements no less than localnodes
522   for (localelem = PetscMax(1, localnodes / (degree*degree*degree)); ;
523        localelem++) {
524     Split3(localelem, melem, true);
525     if (Max3(melem) / Min3(melem) <= 2) break;
526   }
527 
528   // Find my location in the process grid
529   ierr = MPI_Comm_rank(comm, &rank); CHKERRQ(ierr);
530   for (int d=0,rankleft=rank; d<dim; d++) {
531     const int pstride[3] = {p[1] *p[2], p[2], 1};
532     irank[d] = rankleft / pstride[d];
533     rankleft -= irank[d] * pstride[d];
534   }
535 
536   GlobalNodes(p, irank, degree, melem, mnodes);
537 
538   // Setup global vector
539   ierr = VecCreate(comm, &X); CHKERRQ(ierr);
540   if (memtyperequested == CEED_MEM_DEVICE) {
541     ierr = VecSetType(X, VECCUDA); CHKERRQ(ierr);
542   }
543   ierr = VecSetSizes(X, mnodes[0]*mnodes[1]*mnodes[2]*ncompu, PETSC_DECIDE);
544   CHKERRQ(ierr);
545   ierr = VecSetUp(X); CHKERRQ(ierr);
546 
547   // Set up libCEED
548   CeedInit(ceedresource, &ceed);
549 
550   // Print summary
551   CeedInt gsize;
552   ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
553   if (!test_mode) {
554     const char *usedresource;
555     CeedGetResource(ceed, &usedresource);
556 
557     VecType vectype;
558     ierr = VecGetType(X, &vectype); CHKERRQ(ierr);
559 
560     ierr = PetscPrintf(comm,
561                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n"
562                        "  PETSc:\n"
563                        "    PETSc Vec Type                     : %s\n"
564                        "  libCEED:\n"
565                        "    libCEED Backend                    : %s\n"
566                        "    libCEED Backend MemType            : %s\n"
567                        "    libCEED User Requested MemType     : %s\n"
568                        "  Mesh:\n"
569                        "    Number of 1D Basis Nodes (p)       : %d\n"
570                        "    Number of 1D Quadrature Points (q) : %d\n"
571                        "    Global nodes                       : %D\n"
572                        "    Process Decomposition              : %D %D %D\n"
573                        "    Local Elements                     : %D = %D %D %D\n"
574                        "    Owned nodes                        : %D = %D %D %D\n"
575                        "    DoF per node                       : %D\n",
576                        bpchoice+1, vectype, usedresource,
577                        CeedMemTypes[memtypebackend],
578                        (setmemtyperequest) ?
579                        CeedMemTypes[memtyperequested] : "none",
580                        P, Q,  gsize/ncompu, p[0], p[1], p[2], localelem,
581                        melem[0], melem[1], melem[2],
582                        mnodes[0]*mnodes[1]*mnodes[2], mnodes[0], mnodes[1],
583                        mnodes[2], ncompu); CHKERRQ(ierr);
584   }
585 
586   {
587     lsize = 1;
588     for (int d=0; d<dim; d++) {
589       lnodes[d] = melem[d]*degree + 1;
590       lsize *= lnodes[d];
591     }
592     ierr = VecCreate(PETSC_COMM_SELF, &Xloc); CHKERRQ(ierr);
593     if (memtyperequested == CEED_MEM_DEVICE) {
594       ierr = VecSetType(Xloc, VECCUDA); CHKERRQ(ierr);
595     }
596     ierr = VecSetSizes(Xloc, lsize*ncompu, PETSC_DECIDE); CHKERRQ(ierr);
597     ierr = VecSetUp(Xloc); CHKERRQ(ierr);
598 
599     // Create local-to-global scatter
600     PetscInt *ltogind, *ltogind0, *locind, l0count;
601     IS ltogis, ltogis0, locis;
602     PetscInt gstart[2][2][2], gmnodes[2][2][2][dim];
603 
604     for (int i=0; i<2; i++) {
605       for (int j=0; j<2; j++) {
606         for (int k=0; k<2; k++) {
607           PetscInt ijkrank[3] = {irank[0]+i, irank[1]+j, irank[2]+k};
608           gstart[i][j][k] = GlobalStart(p, ijkrank, degree, melem);
609           GlobalNodes(p, ijkrank, degree, melem, gmnodes[i][j][k]);
610         }
611       }
612     }
613 
614     ierr = PetscMalloc1(lsize, &ltogind); CHKERRQ(ierr);
615     ierr = PetscMalloc1(lsize, &ltogind0); CHKERRQ(ierr);
616     ierr = PetscMalloc1(lsize, &locind); CHKERRQ(ierr);
617     l0count = 0;
618     for (PetscInt i=0,ir,ii; ir=i>=mnodes[0], ii=i-ir*mnodes[0], i<lnodes[0]; i++)
619       for (PetscInt j=0,jr,jj; jr=j>=mnodes[1], jj=j-jr*mnodes[1], j<lnodes[1]; j++)
620         for (PetscInt k=0,kr,kk; kr=k>=mnodes[2], kk=k-kr*mnodes[2], k<lnodes[2]; k++) {
621           PetscInt here = (i*lnodes[1]+j)*lnodes[2]+k;
622           ltogind[here] =
623             gstart[ir][jr][kr] + (ii*gmnodes[ir][jr][kr][1]+jj)*gmnodes[ir][jr][kr][2]+kk;
624           if ((irank[0] == 0 && i == 0)
625               || (irank[1] == 0 && j == 0)
626               || (irank[2] == 0 && k == 0)
627               || (irank[0]+1 == p[0] && i+1 == lnodes[0])
628               || (irank[1]+1 == p[1] && j+1 == lnodes[1])
629               || (irank[2]+1 == p[2] && k+1 == lnodes[2]))
630             continue;
631           ltogind0[l0count] = ltogind[here];
632           locind[l0count++] = here;
633         }
634     ierr = ISCreateBlock(comm, ncompu, lsize, ltogind, PETSC_OWN_POINTER,
635                          &ltogis); CHKERRQ(ierr);
636     ierr = VecScatterCreate(Xloc, NULL, X, ltogis, &ltog); CHKERRQ(ierr);
637     CHKERRQ(ierr);
638     ierr = ISCreateBlock(comm, ncompu, l0count, ltogind0, PETSC_OWN_POINTER,
639                          &ltogis0); CHKERRQ(ierr);
640     ierr = ISCreateBlock(comm, ncompu, l0count, locind, PETSC_OWN_POINTER,
641                          &locis); CHKERRQ(ierr);
642     ierr = VecScatterCreate(Xloc, locis, X, ltogis0, &ltog0); CHKERRQ(ierr);
643     {
644       // Create global-to-global scatter for Dirichlet values (everything not in
645       // ltogis0, which is the range of ltog0)
646       PetscInt xstart, xend, *indD, countD = 0;
647       IS isD;
648       const PetscScalar *x;
649       ierr = VecZeroEntries(Xloc); CHKERRQ(ierr);
650       ierr = VecSet(X, 1.0); CHKERRQ(ierr);
651       ierr = VecScatterBegin(ltog0, Xloc, X, INSERT_VALUES, SCATTER_FORWARD);
652       CHKERRQ(ierr);
653       ierr = VecScatterEnd(ltog0, Xloc, X, INSERT_VALUES, SCATTER_FORWARD);
654       CHKERRQ(ierr);
655       ierr = VecGetOwnershipRange(X, &xstart, &xend); CHKERRQ(ierr);
656       ierr = PetscMalloc1(xend-xstart, &indD); CHKERRQ(ierr);
657       ierr = VecGetArrayRead(X, &x); CHKERRQ(ierr);
658       for (PetscInt i=0; i<xend-xstart; i++) {
659         if (x[i] == 1.) indD[countD++] = xstart + i;
660       }
661       ierr = VecRestoreArrayRead(X, &x); CHKERRQ(ierr);
662       ierr = ISCreateGeneral(comm, countD, indD, PETSC_COPY_VALUES, &isD);
663       CHKERRQ(ierr);
664       ierr = PetscFree(indD); CHKERRQ(ierr);
665       ierr = VecScatterCreate(X, isD, X, isD, &gtogD); CHKERRQ(ierr);
666       ierr = ISDestroy(&isD); CHKERRQ(ierr);
667     }
668     ierr = ISDestroy(&ltogis); CHKERRQ(ierr);
669     ierr = ISDestroy(&ltogis0); CHKERRQ(ierr);
670     ierr = ISDestroy(&locis); CHKERRQ(ierr);
671   }
672 
673   // CEED bases
674   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompu, P, Q,
675                                   bpOptions[bpchoice].qmode, &basisu);
676   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, Q,
677                                   bpOptions[bpchoice].qmode, &basisx);
678 
679   // CEED restrictions
680   CreateRestriction(ceed, CEED_INTERLACED, melem, P, ncompu, &Erestrictu);
681   CreateRestriction(ceed, CEED_NONINTERLACED, melem, 2, dim, &Erestrictx);
682   CeedInt nelem = melem[0]*melem[1]*melem[2];
683   CeedElemRestrictionCreateStrided(ceed, nelem, Q*Q*Q, nelem*Q*Q*Q, ncompu,
684                                    CEED_STRIDES_BACKEND, &Erestrictui);
685   CeedElemRestrictionCreateStrided(ceed, nelem, Q*Q*Q, nelem*Q*Q*Q,
686                                    bpOptions[bpchoice].qdatasize,
687                                    CEED_STRIDES_BACKEND, &Erestrictqdi);
688   {
689     CeedScalar *xloc;
690     CeedInt shape[3] = {melem[0]+1, melem[1]+1, melem[2]+1}, len =
691                          shape[0]*shape[1]*shape[2];
692     xloc = malloc(len*ncompx*sizeof xloc[0]);
693     for (CeedInt i=0; i<shape[0]; i++) {
694       for (CeedInt j=0; j<shape[1]; j++) {
695         for (CeedInt k=0; k<shape[2]; k++) {
696           xloc[((i*shape[1]+j)*shape[2]+k) + 0*len] = 1.*(irank[0]*melem[0]+i) /
697               (p[0]*melem[0]);
698           xloc[((i*shape[1]+j)*shape[2]+k) + 1*len] = 1.*(irank[1]*melem[1]+j) /
699               (p[1]*melem[1]);
700           xloc[((i*shape[1]+j)*shape[2]+k) + 2*len] = 1.*(irank[2]*melem[2]+k) /
701               (p[2]*melem[2]);
702         }
703       }
704     }
705     CeedVectorCreate(ceed, len*ncompx, &xcoord);
706     CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_OWN_POINTER, xloc);
707   }
708 
709   // Create the Qfunction that builds the operator quadrature data
710   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].setupgeo,
711                               bpOptions[bpchoice].setupgeofname, &qfsetupgeo);
712   CeedQFunctionAddInput(qfsetupgeo, "dx", ncompx*dim, CEED_EVAL_GRAD);
713   CeedQFunctionAddInput(qfsetupgeo, "weight", 1, CEED_EVAL_WEIGHT);
714   CeedQFunctionAddOutput(qfsetupgeo, "qdata", bpOptions[bpchoice].qdatasize,
715                          CEED_EVAL_NONE);
716 
717   // Create the Qfunction that sets up the RHS and true solution
718   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].setuprhs,
719                               bpOptions[bpchoice].setuprhsfname, &qfsetuprhs);
720   CeedQFunctionAddInput(qfsetuprhs, "x", ncompx, CEED_EVAL_INTERP);
721   CeedQFunctionAddInput(qfsetuprhs, "dx", ncompx*dim, CEED_EVAL_GRAD);
722   CeedQFunctionAddInput(qfsetuprhs, "weight", 1, CEED_EVAL_WEIGHT);
723   CeedQFunctionAddOutput(qfsetuprhs, "true_soln", ncompu, CEED_EVAL_NONE);
724   CeedQFunctionAddOutput(qfsetuprhs, "rhs", ncompu, CEED_EVAL_INTERP);
725 
726   // Set up PDE operator
727   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].apply,
728                               bpOptions[bpchoice].applyfname, &qfapply);
729   // Add inputs and outputs
730   CeedInt inscale = bpOptions[bpchoice].inmode==CEED_EVAL_GRAD ? 3 : 1;
731   CeedInt outscale = bpOptions[bpchoice].outmode==CEED_EVAL_GRAD ? 3 : 1;
732   CeedQFunctionAddInput(qfapply, "u", ncompu*inscale,
733                         bpOptions[bpchoice].inmode);
734   CeedQFunctionAddInput(qfapply, "qdata", bpOptions[bpchoice].qdatasize,
735                         CEED_EVAL_NONE);
736   CeedQFunctionAddOutput(qfapply, "v", ncompu*outscale,
737                          bpOptions[bpchoice].outmode);
738 
739   // Create the error qfunction
740   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].error,
741                               bpOptions[bpchoice].errorfname, &qferror);
742   CeedQFunctionAddInput(qferror, "u", ncompu, CEED_EVAL_INTERP);
743   CeedQFunctionAddInput(qferror, "true_soln", ncompu, CEED_EVAL_NONE);
744   CeedQFunctionAddOutput(qferror, "error", ncompu, CEED_EVAL_NONE);
745 
746   // Create the persistent vectors that will be needed in setup
747   CeedInt nqpts;
748   CeedBasisGetNumQuadraturePoints(basisu, &nqpts);
749   CeedVectorCreate(ceed, bpOptions[bpchoice].qdatasize*nelem*nqpts, &qdata);
750   CeedVectorCreate(ceed, nelem*nqpts*ncompu, &target);
751   CeedVectorCreate(ceed, lsize*ncompu, &rhsceed);
752 
753   // Create the operator that builds the quadrature data for the ceed operator
754   CeedOperatorCreate(ceed, qfsetupgeo, CEED_QFUNCTION_NONE,
755                      CEED_QFUNCTION_NONE, &opsetupgeo);
756   CeedOperatorSetField(opsetupgeo, "dx", Erestrictx, basisx,
757                        CEED_VECTOR_ACTIVE);
758   CeedOperatorSetField(opsetupgeo, "weight", CEED_ELEMRESTRICTION_NONE, basisx,
759                        CEED_VECTOR_NONE);
760   CeedOperatorSetField(opsetupgeo, "qdata", Erestrictqdi,
761                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
762 
763   // Create the operator that builds the RHS and true solution
764   CeedOperatorCreate(ceed, qfsetuprhs, CEED_QFUNCTION_NONE,
765                      CEED_QFUNCTION_NONE, &opsetuprhs);
766   CeedOperatorSetField(opsetuprhs, "x", Erestrictx, basisx,
767                        CEED_VECTOR_ACTIVE);
768   CeedOperatorSetField(opsetuprhs, "dx", Erestrictx, basisx,
769                        CEED_VECTOR_ACTIVE);
770   CeedOperatorSetField(opsetuprhs, "weight", CEED_ELEMRESTRICTION_NONE, basisx,
771                        CEED_VECTOR_NONE);
772   CeedOperatorSetField(opsetuprhs, "true_soln", Erestrictui,
773                        CEED_BASIS_COLLOCATED, target);
774   CeedOperatorSetField(opsetuprhs, "rhs", Erestrictu, basisu,
775                        CEED_VECTOR_ACTIVE);
776 
777   // Create the mass or diff operator
778   CeedOperatorCreate(ceed, qfapply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
779                      &opapply);
780   CeedOperatorSetField(opapply, "u", Erestrictu, basisu, CEED_VECTOR_ACTIVE);
781   CeedOperatorSetField(opapply, "qdata", Erestrictqdi, CEED_BASIS_COLLOCATED,
782                        qdata);
783   CeedOperatorSetField(opapply, "v", Erestrictu, basisu, CEED_VECTOR_ACTIVE);
784 
785   // Create the error operator
786   CeedOperatorCreate(ceed, qferror, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
787                      &operror);
788   CeedOperatorSetField(operror, "u", Erestrictu, basisu, CEED_VECTOR_ACTIVE);
789   CeedOperatorSetField(operror, "true_soln", Erestrictui,
790                        CEED_BASIS_COLLOCATED, target);
791   CeedOperatorSetField(operror, "error", Erestrictui, CEED_BASIS_COLLOCATED,
792                        CEED_VECTOR_ACTIVE);
793 
794   // Set up Mat
795   ierr = PetscMalloc1(1, &user); CHKERRQ(ierr);
796   user->comm = comm;
797   user->ltog = ltog;
798   if (bpchoice != CEED_BP1 && bpchoice != CEED_BP2) {
799     user->ltog0 = ltog0;
800     user->gtogD = gtogD;
801   }
802   user->Xloc = Xloc;
803   ierr = VecDuplicate(Xloc, &user->Yloc); CHKERRQ(ierr);
804   CeedVectorCreate(ceed, lsize*ncompu, &user->xceed);
805   CeedVectorCreate(ceed, lsize*ncompu, &user->yceed);
806   user->op = opapply;
807   user->qdata = qdata;
808   user->ceed = ceed;
809   user->memtype = memtyperequested;
810   if (memtyperequested == CEED_MEM_HOST) {
811     user->VecGetArray = VecGetArray;
812     user->VecGetArrayRead = VecGetArrayRead;
813     user->VecRestoreArray = VecRestoreArray;
814     user->VecRestoreArrayRead = VecRestoreArrayRead;
815   } else {
816     user->VecGetArray = VecCUDAGetArray;
817     user->VecGetArrayRead = VecCUDAGetArrayRead;
818     user->VecRestoreArray = VecCUDARestoreArray;
819     user->VecRestoreArrayRead = VecCUDARestoreArrayRead;
820   }
821 
822   ierr = MatCreateShell(comm, mnodes[0]*mnodes[1]*mnodes[2]*ncompu,
823                         mnodes[0]*mnodes[1]*mnodes[2]*ncompu,
824                         PETSC_DECIDE, PETSC_DECIDE, user, &mat); CHKERRQ(ierr);
825   if (bpchoice == CEED_BP1 || bpchoice == CEED_BP2) {
826     ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Mass);
827     CHKERRQ(ierr);
828   } else {
829     ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Diff);
830     CHKERRQ(ierr);
831   }
832   if (user->memtype == CEED_MEM_DEVICE) {
833     ierr = MatShellSetVecType(mat, VECCUDA); CHKERRQ(ierr);
834   }
835 
836   // Get RHS vector
837   ierr = VecDuplicate(X, &rhs); CHKERRQ(ierr);
838   ierr = VecDuplicate(Xloc, &rhsloc); CHKERRQ(ierr);
839   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
840   ierr = user->VecGetArray(rhsloc, &r); CHKERRQ(ierr);
841   CeedVectorSetArray(rhsceed, user->memtype, CEED_USE_POINTER, r);
842 
843   // Setup qdata, rhs, and target
844   CeedOperatorApply(opsetupgeo, xcoord, qdata, CEED_REQUEST_IMMEDIATE);
845   CeedOperatorApply(opsetuprhs, xcoord, rhsceed, CEED_REQUEST_IMMEDIATE);
846   ierr = CeedVectorSyncArray(rhsceed, user->memtype); CHKERRQ(ierr);
847   CeedVectorDestroy(&xcoord);
848 
849   // Gather RHS
850   ierr = user->VecRestoreArray(rhsloc, &r); CHKERRQ(ierr);
851   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
852   ierr = VecScatterBegin(ltog, rhsloc, rhs, ADD_VALUES, SCATTER_FORWARD);
853   CHKERRQ(ierr);
854   ierr = VecScatterEnd(ltog, rhsloc, rhs, ADD_VALUES, SCATTER_FORWARD);
855   CHKERRQ(ierr);
856   CeedVectorDestroy(&rhsceed);
857 
858   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
859   {
860     PC pc;
861     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
862     if (bpchoice == CEED_BP1 || bpchoice == CEED_BP2) {
863       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
864       ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr);
865     } else {
866       ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr);
867     }
868     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
869     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
870     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
871                             PETSC_DEFAULT); CHKERRQ(ierr);
872   }
873   ierr = KSPSetOperators(ksp, mat, mat); CHKERRQ(ierr);
874   // First run, if benchmarking
875   if (benchmark_mode) {
876     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
877     CHKERRQ(ierr);
878     my_rt_start = MPI_Wtime();
879     ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
880     my_rt = MPI_Wtime() - my_rt_start;
881     ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
882     CHKERRQ(ierr);
883     // Set maxits based on first iteration timing
884     if (my_rt > 0.02) {
885       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
886       CHKERRQ(ierr);
887     } else {
888       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
889       CHKERRQ(ierr);
890     }
891   }
892   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
893 
894   // Timed solve
895   ierr = VecZeroEntries(X); CHKERRQ(ierr);
896   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
897 
898   // -- Performance logging
899   ierr = PetscLogStageRegister("Solve Stage", &solvestage); CHKERRQ(ierr);
900   ierr = PetscLogStagePush(solvestage); CHKERRQ(ierr);
901 
902   // -- Solve
903   my_rt_start = MPI_Wtime();
904   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
905   my_rt = MPI_Wtime() - my_rt_start;
906 
907   // -- Performance logging
908   ierr = PetscLogStagePop();
909 
910   // Output results
911   {
912     KSPType ksptype;
913     KSPConvergedReason reason;
914     PetscReal rnorm;
915     PetscInt its;
916     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
917     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
918     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
919     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
920     if (!test_mode || reason < 0 || rnorm > 1e-8) {
921       ierr = PetscPrintf(comm,
922                          "  KSP:\n"
923                          "    KSP Type                           : %s\n"
924                          "    KSP Convergence                    : %s\n"
925                          "    Total KSP Iterations               : %D\n"
926                          "    Final rnorm                        : %e\n",
927                          ksptype, KSPConvergedReasons[reason], its,
928                          (double)rnorm); CHKERRQ(ierr);
929     }
930     if (!test_mode) {
931       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
932     }
933     {
934       PetscReal maxerror;
935       ierr = ComputeErrorMax(user, operror, X, target, &maxerror);
936       CHKERRQ(ierr);
937       PetscReal tol = 5e-2;
938       if (!test_mode || maxerror > tol) {
939         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
940         CHKERRQ(ierr);
941         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
942         CHKERRQ(ierr);
943         ierr = PetscPrintf(comm,
944                            "    Pointwise Error (max)              : %e\n"
945                            "    CG Solve Time                      : %g (%g) sec\n",
946                            (double)maxerror, rt_max, rt_min); CHKERRQ(ierr);
947       }
948     }
949     if (benchmark_mode && (!test_mode)) {
950       ierr = PetscPrintf(comm,
951                          "    DoFs/Sec in CG                     : %g (%g) million\n",
952                          1e-6*gsize*its/rt_max,
953                          1e-6*gsize*its/rt_min); CHKERRQ(ierr);
954     }
955   }
956 
957   if (write_solution) {
958     PetscViewer vtkviewersoln;
959 
960     ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr);
961     ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr);
962     ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr);
963     ierr = VecView(X, vtkviewersoln); CHKERRQ(ierr);
964     ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr);
965   }
966 
967   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
968   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
969   ierr = VecDestroy(&X); CHKERRQ(ierr);
970   ierr = VecDestroy(&user->Xloc); CHKERRQ(ierr);
971   ierr = VecDestroy(&user->Yloc); CHKERRQ(ierr);
972   ierr = VecScatterDestroy(&ltog); CHKERRQ(ierr);
973   ierr = VecScatterDestroy(&ltog0); CHKERRQ(ierr);
974   ierr = VecScatterDestroy(&gtogD); CHKERRQ(ierr);
975   ierr = MatDestroy(&mat); CHKERRQ(ierr);
976   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
977 
978   CeedVectorDestroy(&user->xceed);
979   CeedVectorDestroy(&user->yceed);
980   CeedVectorDestroy(&user->qdata);
981   CeedVectorDestroy(&target);
982   CeedOperatorDestroy(&opsetupgeo);
983   CeedOperatorDestroy(&opsetuprhs);
984   CeedOperatorDestroy(&opapply);
985   CeedOperatorDestroy(&operror);
986   CeedElemRestrictionDestroy(&Erestrictu);
987   CeedElemRestrictionDestroy(&Erestrictx);
988   CeedElemRestrictionDestroy(&Erestrictui);
989   CeedElemRestrictionDestroy(&Erestrictqdi);
990   CeedQFunctionDestroy(&qfsetupgeo);
991   CeedQFunctionDestroy(&qfsetuprhs);
992   CeedQFunctionDestroy(&qfapply);
993   CeedQFunctionDestroy(&qferror);
994   CeedBasisDestroy(&basisu);
995   CeedBasisDestroy(&basisx);
996   CeedDestroy(&ceed);
997   ierr = PetscFree(user); CHKERRQ(ierr);
998   return PetscFinalize();
999 }
1000