xref: /libCEED/examples/petsc/bps.c (revision 6bd9afcaf313bb7ac3bbee04731c2cb873b0c958)
1819eb1b3Sjeremylt // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2819eb1b3Sjeremylt // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3819eb1b3Sjeremylt // reserved. See files LICENSE and NOTICE for details.
4819eb1b3Sjeremylt //
5819eb1b3Sjeremylt // This file is part of CEED, a collection of benchmarks, miniapps, software
6819eb1b3Sjeremylt // libraries and APIs for efficient high-order finite element and spectral
7819eb1b3Sjeremylt // element discretizations for exascale applications. For more information and
8819eb1b3Sjeremylt // source code availability see http://github.com/ceed.
9819eb1b3Sjeremylt //
10819eb1b3Sjeremylt // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11819eb1b3Sjeremylt // a collaborative effort of two U.S. Department of Energy organizations (Office
12819eb1b3Sjeremylt // of Science and the National Nuclear Security Administration) responsible for
13819eb1b3Sjeremylt // the planning and preparation of a capable exascale ecosystem, including
14819eb1b3Sjeremylt // software, applications, hardware, advanced system engineering and early
15819eb1b3Sjeremylt // testbed platforms, in support of the nation's exascale computing imperative.
16819eb1b3Sjeremylt 
170c59ef15SJeremy L Thompson //                        libCEED + PETSc Example: CEED BPs
180c59ef15SJeremy L Thompson //
190c59ef15SJeremy L Thompson // This example demonstrates a simple usage of libCEED with PETSc to solve the
200c59ef15SJeremy L Thompson // CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps.
210c59ef15SJeremy L Thompson //
220c59ef15SJeremy L Thompson // The code is intentionally "raw", using only low-level communication
230c59ef15SJeremy L Thompson // primitives.
240c59ef15SJeremy L Thompson //
250c59ef15SJeremy L Thompson // Build with:
260c59ef15SJeremy L Thompson //
270c59ef15SJeremy L Thompson //     make bps [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
280c59ef15SJeremy L Thompson //
290c59ef15SJeremy L Thompson // Sample runs:
300c59ef15SJeremy L Thompson //
310c59ef15SJeremy L Thompson //     bps -problem bp1
320c59ef15SJeremy L Thompson //     bps -problem bp2 -ceed /cpu/self
330c59ef15SJeremy L Thompson //     bps -problem bp3 -ceed /gpu/occa
340c59ef15SJeremy L Thompson //     bps -problem bp4 -ceed /cpu/occa
350c59ef15SJeremy L Thompson //     bps -problem bp5 -ceed /omp/occa
360c59ef15SJeremy L Thompson //     bps -problem bp6 -ceed /ocl/occa
370c59ef15SJeremy L Thompson //
380c59ef15SJeremy L Thompson //TESTARGS -ceed {ceed_resource} -test -problem bp3
390c59ef15SJeremy L Thompson 
400c59ef15SJeremy L Thompson /// @file
410c59ef15SJeremy L Thompson /// CEED BPs example using PETSc
42819eb1b3Sjeremylt /// See bpsdmplex.c for an implementation using DMPlex unstructured grids.
430c59ef15SJeremy L Thompson const char help[] = "Solve CEED BPs using PETSc\n";
440c59ef15SJeremy L Thompson 
450c59ef15SJeremy L Thompson #include <stdbool.h>
460c59ef15SJeremy L Thompson #include <string.h>
470c59ef15SJeremy L Thompson #include "common.h"
480c59ef15SJeremy L Thompson #include "bp1.h"
490c59ef15SJeremy L Thompson #include "bp2.h"
500c59ef15SJeremy L Thompson #include "bp3.h"
510c59ef15SJeremy L Thompson #include "bp4.h"
520c59ef15SJeremy L Thompson 
530c59ef15SJeremy L Thompson #define PATH(BASE) __DIR__ #BASE
540c59ef15SJeremy L Thompson 
550c59ef15SJeremy L Thompson static void Split3(PetscInt size, PetscInt m[3], bool reverse) {
560c59ef15SJeremy L Thompson   for (PetscInt d=0,sizeleft=size; d<3; d++) {
570c59ef15SJeremy L Thompson     PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(sizeleft, 1./(3 - d)));
580c59ef15SJeremy L Thompson     while (try * (sizeleft / try) != sizeleft) try++;
590c59ef15SJeremy L Thompson     m[reverse ? 2-d : d] = try;
600c59ef15SJeremy L Thompson     sizeleft /= try;
610c59ef15SJeremy L Thompson   }
620c59ef15SJeremy L Thompson }
630c59ef15SJeremy L Thompson 
640c59ef15SJeremy L Thompson static PetscInt Max3(const PetscInt a[3]) {
650c59ef15SJeremy L Thompson   return PetscMax(a[0], PetscMax(a[1], a[2]));
660c59ef15SJeremy L Thompson }
670c59ef15SJeremy L Thompson static PetscInt Min3(const PetscInt a[3]) {
680c59ef15SJeremy L Thompson   return PetscMin(a[0], PetscMin(a[1], a[2]));
690c59ef15SJeremy L Thompson }
700c59ef15SJeremy L Thompson static void GlobalDof(const PetscInt p[3], const PetscInt irank[3],
710c59ef15SJeremy L Thompson                       PetscInt degree, const PetscInt melem[3],
720c59ef15SJeremy L Thompson                       PetscInt mdof[3]) {
730c59ef15SJeremy L Thompson   for (int d=0; d<3; d++)
740c59ef15SJeremy L Thompson     mdof[d] = degree*melem[d] + (irank[d] == p[d]-1);
750c59ef15SJeremy L Thompson }
760c59ef15SJeremy L Thompson static PetscInt GlobalStart(const PetscInt p[3], const PetscInt irank[3],
770c59ef15SJeremy L Thompson                             PetscInt degree, const PetscInt melem[3]) {
780c59ef15SJeremy L Thompson   PetscInt start = 0;
790c59ef15SJeremy L Thompson   // Dumb brute-force is easier to read
800c59ef15SJeremy L Thompson   for (PetscInt i=0; i<p[0]; i++) {
810c59ef15SJeremy L Thompson     for (PetscInt j=0; j<p[1]; j++) {
820c59ef15SJeremy L Thompson       for (PetscInt k=0; k<p[2]; k++) {
830c59ef15SJeremy L Thompson         PetscInt mdof[3], ijkrank[] = {i,j,k};
840c59ef15SJeremy L Thompson         if (i == irank[0] && j == irank[1] && k == irank[2]) return start;
850c59ef15SJeremy L Thompson         GlobalDof(p, ijkrank, degree, melem, mdof);
860c59ef15SJeremy L Thompson         start += mdof[0] * mdof[1] * mdof[2];
870c59ef15SJeremy L Thompson       }
880c59ef15SJeremy L Thompson     }
890c59ef15SJeremy L Thompson   }
900c59ef15SJeremy L Thompson   return -1;
910c59ef15SJeremy L Thompson }
920c59ef15SJeremy L Thompson static int CreateRestriction(Ceed ceed, const CeedInt melem[3],
930c59ef15SJeremy L Thompson                              CeedInt P, CeedInt ncomp,
940c59ef15SJeremy L Thompson                              CeedElemRestriction *Erestrict) {
95819eb1b3Sjeremylt   const PetscInt nelem = melem[0]*melem[1]*melem[2];
960c59ef15SJeremy L Thompson   PetscInt mdof[3], *idx, *idxp;
970c59ef15SJeremy L Thompson 
98819eb1b3Sjeremylt   // Get indicies
990c59ef15SJeremy L Thompson   for (int d=0; d<3; d++) mdof[d] = melem[d]*(P-1) + 1;
100819eb1b3Sjeremylt   idxp = idx = malloc(nelem*P*P*P*sizeof idx[0]);
1010c59ef15SJeremy L Thompson   for (CeedInt i=0; i<melem[0]; i++) {
1020c59ef15SJeremy L Thompson     for (CeedInt j=0; j<melem[1]; j++) {
1030c59ef15SJeremy L Thompson       for (CeedInt k=0; k<melem[2]; k++,idxp += P*P*P) {
1040c59ef15SJeremy L Thompson         for (CeedInt ii=0; ii<P; ii++) {
1050c59ef15SJeremy L Thompson           for (CeedInt jj=0; jj<P; jj++) {
1060c59ef15SJeremy L Thompson             for (CeedInt kk=0; kk<P; kk++) {
1070c59ef15SJeremy L Thompson               if (0) { // This is the C-style (i,j,k) ordering that I prefer
1080c59ef15SJeremy L Thompson                 idxp[(ii*P+jj)*P+kk] = (((i*(P-1)+ii)*mdof[1]
1090c59ef15SJeremy L Thompson                                          + (j*(P-1)+jj))*mdof[2]
1100c59ef15SJeremy L Thompson                                         + (k*(P-1)+kk));
1110c59ef15SJeremy L Thompson               } else { // (k,j,i) ordering for consistency with MFEM example
1120c59ef15SJeremy L Thompson                 idxp[ii+P*(jj+P*kk)] = (((i*(P-1)+ii)*mdof[1]
1130c59ef15SJeremy L Thompson                                          + (j*(P-1)+jj))*mdof[2]
1140c59ef15SJeremy L Thompson                                         + (k*(P-1)+kk));
1150c59ef15SJeremy L Thompson               }
1160c59ef15SJeremy L Thompson             }
1170c59ef15SJeremy L Thompson           }
1180c59ef15SJeremy L Thompson         }
1190c59ef15SJeremy L Thompson       }
1200c59ef15SJeremy L Thompson     }
1210c59ef15SJeremy L Thompson   }
122819eb1b3Sjeremylt 
123819eb1b3Sjeremylt   // Setup CEED restriction
124819eb1b3Sjeremylt   CeedElemRestrictionCreate(ceed, nelem, P*P*P, mdof[0]*mdof[1]*mdof[2], ncomp,
1250c59ef15SJeremy L Thompson                             CEED_MEM_HOST, CEED_OWN_POINTER, idx, Erestrict);
126819eb1b3Sjeremylt 
1270c59ef15SJeremy L Thompson   PetscFunctionReturn(0);
1280c59ef15SJeremy L Thompson }
1290c59ef15SJeremy L Thompson 
1300c59ef15SJeremy L Thompson // Data for PETSc
1310c59ef15SJeremy L Thompson typedef struct User_ *User;
1320c59ef15SJeremy L Thompson struct User_ {
1330c59ef15SJeremy L Thompson   MPI_Comm comm;
1340c59ef15SJeremy L Thompson   VecScatter ltog;              // Scatter for all entries
1350c59ef15SJeremy L Thompson   VecScatter ltog0;             // Skip Dirichlet values
1360c59ef15SJeremy L Thompson   VecScatter gtogD;             // global-to-global; only Dirichlet values
1370c59ef15SJeremy L Thompson   Vec Xloc, Yloc;
1380c59ef15SJeremy L Thompson   CeedVector xceed, yceed;
1390c59ef15SJeremy L Thompson   CeedOperator op;
1400c59ef15SJeremy L Thompson   CeedVector rho;
1410c59ef15SJeremy L Thompson   Ceed ceed;
1420c59ef15SJeremy L Thompson };
1430c59ef15SJeremy L Thompson 
1440c59ef15SJeremy L Thompson // BP Options
1450c59ef15SJeremy L Thompson typedef enum {
1460c59ef15SJeremy L Thompson   CEED_BP1 = 0, CEED_BP2 = 1, CEED_BP3 = 2,
1470c59ef15SJeremy L Thompson   CEED_BP4 = 3, CEED_BP5 = 4, CEED_BP6 = 5
1480c59ef15SJeremy L Thompson } bpType;
1490c59ef15SJeremy L Thompson static const char *const bpTypes[] = {"bp1","bp2","bp3","bp4","bp5","bp6",
1500c59ef15SJeremy L Thompson                                       "bpType","CEED_BP",0};
1510c59ef15SJeremy L Thompson 
1520c59ef15SJeremy L Thompson // BP specific data
1530c59ef15SJeremy L Thompson typedef struct {
1540c59ef15SJeremy L Thompson   CeedInt vscale, qdatasize, qextra;
1550c59ef15SJeremy L Thompson   CeedQFunctionUser setup, apply, error;
1560c59ef15SJeremy L Thompson   const char setupfname[PETSC_MAX_PATH_LEN], applyfname[PETSC_MAX_PATH_LEN],
1570c59ef15SJeremy L Thompson              errorfname[PETSC_MAX_PATH_LEN];
1580c59ef15SJeremy L Thompson   CeedEvalMode inmode, outmode;
159e07e9ddcSjeremylt   CeedQuadMode qmode;
1600c59ef15SJeremy L Thompson } bpData;
1610c59ef15SJeremy L Thompson 
1620c59ef15SJeremy L Thompson bpData bpOptions[6] = {
1630c59ef15SJeremy L Thompson   [CEED_BP1] = {
1640c59ef15SJeremy L Thompson     .vscale = 1,
1650c59ef15SJeremy L Thompson     .qdatasize = 1,
166194c25f7Sjeremylt     .qextra = 1,
1670c59ef15SJeremy L Thompson     .setup = SetupMass,
1680c59ef15SJeremy L Thompson     .apply = Mass,
1690c59ef15SJeremy L Thompson     .error = Error,
1700c59ef15SJeremy L Thompson     .setupfname = PATH(bp1.h:SetupMass),
1710c59ef15SJeremy L Thompson     .applyfname = PATH(bp1.h:Mass),
1720c59ef15SJeremy L Thompson     .errorfname = PATH(common.h:Error),
1730c59ef15SJeremy L Thompson     .inmode = CEED_EVAL_INTERP,
174e07e9ddcSjeremylt     .outmode = CEED_EVAL_INTERP,
175e07e9ddcSjeremylt     .qmode = CEED_GAUSS},
1760c59ef15SJeremy L Thompson   [CEED_BP2] = {
1770c59ef15SJeremy L Thompson     .vscale = 3,
1780c59ef15SJeremy L Thompson     .qdatasize = 1,
179194c25f7Sjeremylt     .qextra = 1,
1800c59ef15SJeremy L Thompson     .setup = SetupMass3,
1810c59ef15SJeremy L Thompson     .apply = Mass3,
1820c59ef15SJeremy L Thompson     .error = Error3,
1830c59ef15SJeremy L Thompson     .setupfname = PATH(bp2.h:SetupMass3),
1840c59ef15SJeremy L Thompson     .applyfname = PATH(bp2.h:Mass3),
1850c59ef15SJeremy L Thompson     .errorfname = PATH(common.h:Error3),
1860c59ef15SJeremy L Thompson     .inmode = CEED_EVAL_INTERP,
187e07e9ddcSjeremylt     .outmode = CEED_EVAL_INTERP,
188e07e9ddcSjeremylt     .qmode = CEED_GAUSS},
1890c59ef15SJeremy L Thompson   [CEED_BP3] = {
1900c59ef15SJeremy L Thompson     .vscale = 1,
1910c59ef15SJeremy L Thompson     .qdatasize = 6,
192194c25f7Sjeremylt     .qextra = 1,
1930c59ef15SJeremy L Thompson     .setup = SetupDiff,
1940c59ef15SJeremy L Thompson     .apply = Diff,
1950c59ef15SJeremy L Thompson     .error = Error,
1960c59ef15SJeremy L Thompson     .setupfname = PATH(bp3.h:SetupDiff),
1970c59ef15SJeremy L Thompson     .applyfname = PATH(bp3.h:Diff),
1980c59ef15SJeremy L Thompson     .errorfname = PATH(common.h:Error),
1990c59ef15SJeremy L Thompson     .inmode = CEED_EVAL_GRAD,
200e07e9ddcSjeremylt     .outmode = CEED_EVAL_GRAD,
201e07e9ddcSjeremylt     .qmode = CEED_GAUSS},
2020c59ef15SJeremy L Thompson   [CEED_BP4] = {
2030c59ef15SJeremy L Thompson     .vscale = 3,
2040c59ef15SJeremy L Thompson     .qdatasize = 6,
205194c25f7Sjeremylt     .qextra = 1,
2060c59ef15SJeremy L Thompson     .setup = SetupDiff3,
2070c59ef15SJeremy L Thompson     .apply = Diff3,
2084b5b4ec1Sjeremylt     .error = Error3,
2090c59ef15SJeremy L Thompson     .setupfname = PATH(bp4.h:SetupDiff3),
210241a4b83SYohann     .applyfname = PATH(bp4.h:Diff3),
2110c59ef15SJeremy L Thompson     .errorfname = PATH(common.h:Error3),
2120c59ef15SJeremy L Thompson     .inmode = CEED_EVAL_GRAD,
213e07e9ddcSjeremylt     .outmode = CEED_EVAL_GRAD,
214e07e9ddcSjeremylt     .qmode = CEED_GAUSS},
2150c59ef15SJeremy L Thompson   [CEED_BP5] = {
2160c59ef15SJeremy L Thompson     .vscale = 1,
2170c59ef15SJeremy L Thompson     .qdatasize = 6,
218194c25f7Sjeremylt     .qextra = 0,
2190c59ef15SJeremy L Thompson     .setup = SetupDiff,
2200c59ef15SJeremy L Thompson     .apply = Diff,
2210c59ef15SJeremy L Thompson     .error = Error,
2220c59ef15SJeremy L Thompson     .setupfname = PATH(bp3.h:SetupDiff),
2230c59ef15SJeremy L Thompson     .applyfname = PATH(bp3.h:Diff),
2240c59ef15SJeremy L Thompson     .errorfname = PATH(common.h:Error),
2250c59ef15SJeremy L Thompson     .inmode = CEED_EVAL_GRAD,
226e07e9ddcSjeremylt     .outmode = CEED_EVAL_GRAD,
227e07e9ddcSjeremylt     .qmode = CEED_GAUSS_LOBATTO},
2280c59ef15SJeremy L Thompson   [CEED_BP6] = {
2290c59ef15SJeremy L Thompson     .vscale = 3,
2300c59ef15SJeremy L Thompson     .qdatasize = 6,
231194c25f7Sjeremylt     .qextra = 0,
2320c59ef15SJeremy L Thompson     .setup = SetupDiff3,
2330c59ef15SJeremy L Thompson     .apply = Diff3,
2344b5b4ec1Sjeremylt     .error = Error3,
2350c59ef15SJeremy L Thompson     .setupfname = PATH(bp4.h:SetupDiff3),
236241a4b83SYohann     .applyfname = PATH(bp4.h:Diff3),
2370c59ef15SJeremy L Thompson     .errorfname = PATH(common.h:Error3),
2380c59ef15SJeremy L Thompson     .inmode = CEED_EVAL_GRAD,
239e07e9ddcSjeremylt     .outmode = CEED_EVAL_GRAD,
240e07e9ddcSjeremylt     .qmode = CEED_GAUSS_LOBATTO}
2410c59ef15SJeremy L Thompson };
2420c59ef15SJeremy L Thompson 
2430c59ef15SJeremy L Thompson // This function uses libCEED to compute the action of the mass matrix
2440c59ef15SJeremy L Thompson static PetscErrorCode MatMult_Mass(Mat A, Vec X, Vec Y) {
2450c59ef15SJeremy L Thompson   PetscErrorCode ierr;
2460c59ef15SJeremy L Thompson   User user;
2470c59ef15SJeremy L Thompson   PetscScalar *x, *y;
2480c59ef15SJeremy L Thompson 
2490c59ef15SJeremy L Thompson   PetscFunctionBeginUser;
2500c59ef15SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
2510c59ef15SJeremy L Thompson   ierr = VecScatterBegin(user->ltog, X, user->Xloc, INSERT_VALUES,
2520c59ef15SJeremy L Thompson                          SCATTER_REVERSE); CHKERRQ(ierr);
2530c59ef15SJeremy L Thompson   ierr = VecScatterEnd(user->ltog, X, user->Xloc, INSERT_VALUES, SCATTER_REVERSE);
2540c59ef15SJeremy L Thompson   CHKERRQ(ierr);
2550c59ef15SJeremy L Thompson   ierr = VecZeroEntries(user->Yloc); CHKERRQ(ierr);
2560c59ef15SJeremy L Thompson 
2570c59ef15SJeremy L Thompson   ierr = VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr);
2580c59ef15SJeremy L Thompson   ierr = VecGetArray(user->Yloc, &y); CHKERRQ(ierr);
2590c59ef15SJeremy L Thompson   CeedVectorSetArray(user->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
2600c59ef15SJeremy L Thompson   CeedVectorSetArray(user->yceed, CEED_MEM_HOST, CEED_USE_POINTER, y);
2610c59ef15SJeremy L Thompson 
2620c59ef15SJeremy L Thompson   CeedOperatorApply(user->op, user->xceed, user->yceed,
2630c59ef15SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
2640c59ef15SJeremy L Thompson   ierr = CeedVectorSyncArray(user->yceed, CEED_MEM_HOST); CHKERRQ(ierr);
2650c59ef15SJeremy L Thompson 
2660c59ef15SJeremy L Thompson   ierr = VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr);
2670c59ef15SJeremy L Thompson   ierr = VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr);
2680c59ef15SJeremy L Thompson 
2690c59ef15SJeremy L Thompson   if (Y) {
2700c59ef15SJeremy L Thompson     ierr = VecZeroEntries(Y); CHKERRQ(ierr);
2710c59ef15SJeremy L Thompson     ierr = VecScatterBegin(user->ltog, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD);
2720c59ef15SJeremy L Thompson     CHKERRQ(ierr);
2730c59ef15SJeremy L Thompson     ierr = VecScatterEnd(user->ltog, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD);
2740c59ef15SJeremy L Thompson     CHKERRQ(ierr);
2750c59ef15SJeremy L Thompson   }
2760c59ef15SJeremy L Thompson   PetscFunctionReturn(0);
2770c59ef15SJeremy L Thompson }
2780c59ef15SJeremy L Thompson 
2790c59ef15SJeremy L Thompson // This function uses libCEED to compute the action of the Laplacian with
2800c59ef15SJeremy L Thompson // Dirichlet boundary conditions
2810c59ef15SJeremy L Thompson static PetscErrorCode MatMult_Diff(Mat A, Vec X, Vec Y) {
2820c59ef15SJeremy L Thompson   PetscErrorCode ierr;
2830c59ef15SJeremy L Thompson   User user;
2840c59ef15SJeremy L Thompson   PetscScalar *x, *y;
2850c59ef15SJeremy L Thompson 
2860c59ef15SJeremy L Thompson   PetscFunctionBeginUser;
2870c59ef15SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
288819eb1b3Sjeremylt 
289819eb1b3Sjeremylt   // Global-to-local
2900c59ef15SJeremy L Thompson   ierr = VecScatterBegin(user->ltog0, X, user->Xloc, INSERT_VALUES,
2910c59ef15SJeremy L Thompson                          SCATTER_REVERSE); CHKERRQ(ierr);
2920c59ef15SJeremy L Thompson   ierr = VecScatterEnd(user->ltog0, X, user->Xloc, INSERT_VALUES,
2930c59ef15SJeremy L Thompson                        SCATTER_REVERSE);
2940c59ef15SJeremy L Thompson   CHKERRQ(ierr);
2950c59ef15SJeremy L Thompson   ierr = VecZeroEntries(user->Yloc); CHKERRQ(ierr);
2960c59ef15SJeremy L Thompson 
297819eb1b3Sjeremylt   // Setup CEED vectors
2980c59ef15SJeremy L Thompson   ierr = VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr);
2990c59ef15SJeremy L Thompson   ierr = VecGetArray(user->Yloc, &y); CHKERRQ(ierr);
3000c59ef15SJeremy L Thompson   CeedVectorSetArray(user->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
3010c59ef15SJeremy L Thompson   CeedVectorSetArray(user->yceed, CEED_MEM_HOST, CEED_USE_POINTER, y);
3020c59ef15SJeremy L Thompson 
303819eb1b3Sjeremylt   // Apply CEED operator
3040c59ef15SJeremy L Thompson   CeedOperatorApply(user->op, user->xceed, user->yceed,
3050c59ef15SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
3060c59ef15SJeremy L Thompson   ierr = CeedVectorSyncArray(user->yceed, CEED_MEM_HOST); CHKERRQ(ierr);
3070c59ef15SJeremy L Thompson 
308819eb1b3Sjeremylt   // Restore PETSc vectors
3090c59ef15SJeremy L Thompson   ierr = VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr);
3100c59ef15SJeremy L Thompson   ierr = VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr);
3110c59ef15SJeremy L Thompson 
312819eb1b3Sjeremylt   // Local-to-global
3130c59ef15SJeremy L Thompson   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
3140c59ef15SJeremy L Thompson   ierr = VecScatterBegin(user->gtogD, X, Y, INSERT_VALUES, SCATTER_FORWARD);
3150c59ef15SJeremy L Thompson   CHKERRQ(ierr);
3160c59ef15SJeremy L Thompson   ierr = VecScatterEnd(user->gtogD, X, Y, INSERT_VALUES, SCATTER_FORWARD);
3170c59ef15SJeremy L Thompson   CHKERRQ(ierr);
3180c59ef15SJeremy L Thompson   ierr = VecScatterBegin(user->ltog0, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD);
3190c59ef15SJeremy L Thompson   CHKERRQ(ierr);
3200c59ef15SJeremy L Thompson   ierr = VecScatterEnd(user->ltog0, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD);
3210c59ef15SJeremy L Thompson   CHKERRQ(ierr);
322819eb1b3Sjeremylt 
3230c59ef15SJeremy L Thompson   PetscFunctionReturn(0);
3240c59ef15SJeremy L Thompson }
3250c59ef15SJeremy L Thompson 
3260c59ef15SJeremy L Thompson // This function calculates the error in the final solution
3270c59ef15SJeremy L Thompson static PetscErrorCode ComputeErrorMax(User user, CeedOperator op_error, Vec X,
3280c59ef15SJeremy L Thompson                                       CeedVector target, PetscReal *maxerror) {
3290c59ef15SJeremy L Thompson   PetscErrorCode ierr;
3300c59ef15SJeremy L Thompson   PetscScalar *x;
3310c59ef15SJeremy L Thompson   CeedVector collocated_error;
3320c59ef15SJeremy L Thompson   CeedInt length;
3330c59ef15SJeremy L Thompson 
3340c59ef15SJeremy L Thompson   PetscFunctionBeginUser;
3350c59ef15SJeremy L Thompson   CeedVectorGetLength(target, &length);
3360c59ef15SJeremy L Thompson   CeedVectorCreate(user->ceed, length, &collocated_error);
337819eb1b3Sjeremylt 
338819eb1b3Sjeremylt   // Global-to-local
3390c59ef15SJeremy L Thompson   ierr = VecScatterBegin(user->ltog, X, user->Xloc, INSERT_VALUES,
3400c59ef15SJeremy L Thompson                          SCATTER_REVERSE); CHKERRQ(ierr);
3410c59ef15SJeremy L Thompson   ierr = VecScatterEnd(user->ltog, X, user->Xloc, INSERT_VALUES, SCATTER_REVERSE);
3420c59ef15SJeremy L Thompson   CHKERRQ(ierr);
343819eb1b3Sjeremylt 
344819eb1b3Sjeremylt   // Setup CEED vector
3450c59ef15SJeremy L Thompson   ierr = VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr);
3460c59ef15SJeremy L Thompson   CeedVectorSetArray(user->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
347819eb1b3Sjeremylt 
348819eb1b3Sjeremylt   // Apply CEED operator
3490c59ef15SJeremy L Thompson   CeedOperatorApply(op_error, user->xceed, collocated_error,
3500c59ef15SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
351819eb1b3Sjeremylt 
352819eb1b3Sjeremylt   // Restore PETSc vector
3530c59ef15SJeremy L Thompson   VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr);
3540c59ef15SJeremy L Thompson 
355819eb1b3Sjeremylt   // Reduce max error
3560c59ef15SJeremy L Thompson   *maxerror = 0;
3570c59ef15SJeremy L Thompson   const CeedScalar *e;
3580c59ef15SJeremy L Thompson   CeedVectorGetArrayRead(collocated_error, CEED_MEM_HOST, &e);
3590c59ef15SJeremy L Thompson   for (CeedInt i=0; i<length; i++) {
3600c59ef15SJeremy L Thompson     *maxerror = PetscMax(*maxerror, PetscAbsScalar(e[i]));
3610c59ef15SJeremy L Thompson   }
3620c59ef15SJeremy L Thompson   CeedVectorRestoreArrayRead(collocated_error, &e);
36300742fa8SJed Brown   ierr = MPI_Allreduce(MPI_IN_PLACE, maxerror,
36400742fa8SJed Brown                        1, MPIU_REAL, MPIU_MAX, user->comm); CHKERRQ(ierr);
365819eb1b3Sjeremylt 
366819eb1b3Sjeremylt   // Cleanup
3670c59ef15SJeremy L Thompson   CeedVectorDestroy(&collocated_error);
368819eb1b3Sjeremylt 
3690c59ef15SJeremy L Thompson   PetscFunctionReturn(0);
3700c59ef15SJeremy L Thompson }
3710c59ef15SJeremy L Thompson 
3720c59ef15SJeremy L Thompson int main(int argc, char **argv) {
3730c59ef15SJeremy L Thompson   PetscInt ierr;
3740c59ef15SJeremy L Thompson   MPI_Comm comm;
3750c59ef15SJeremy L Thompson   char ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self";
376819eb1b3Sjeremylt   double my_rt_start, my_rt, rt_min, rt_max;
3770c59ef15SJeremy L Thompson   PetscInt degree, qextra, localdof, localelem, melem[3], mdof[3], p[3],
3780c59ef15SJeremy L Thompson            irank[3], ldof[3], lsize, vscale = 1;
3790c59ef15SJeremy L Thompson   PetscScalar *r;
380*6bd9afcaSjeremylt   PetscBool test_mode, benchmark_mode, write_solution;
3810c59ef15SJeremy L Thompson   PetscMPIInt size, rank;
382819eb1b3Sjeremylt   Vec X, Xloc, rhs, rhsloc;
383819eb1b3Sjeremylt   Mat mat;
384819eb1b3Sjeremylt   KSP ksp;
3854b5b4ec1Sjeremylt   VecScatter ltog, ltog0, gtogD;
386819eb1b3Sjeremylt   User user;
3870c59ef15SJeremy L Thompson   Ceed ceed;
3880c59ef15SJeremy L Thompson   CeedBasis basisx, basisu;
3890c59ef15SJeremy L Thompson   CeedElemRestriction Erestrictx, Erestrictu, Erestrictxi, Erestrictui,
3900c59ef15SJeremy L Thompson                       Erestrictqdi;
3910c59ef15SJeremy L Thompson   CeedQFunction qf_setup, qf_apply, qf_error;
3920c59ef15SJeremy L Thompson   CeedOperator op_setup, op_apply, op_error;
3930c59ef15SJeremy L Thompson   CeedVector xcoord, rho, rhsceed, target;
3940c59ef15SJeremy L Thompson   CeedInt P, Q;
3950c59ef15SJeremy L Thompson   bpType bpChoice;
3960c59ef15SJeremy L Thompson 
3970c59ef15SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
3980c59ef15SJeremy L Thompson   if (ierr) return ierr;
3990c59ef15SJeremy L Thompson   comm = PETSC_COMM_WORLD;
4000c59ef15SJeremy L Thompson   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr);
4010c59ef15SJeremy L Thompson   bpChoice = CEED_BP1;
4020c59ef15SJeremy L Thompson   ierr = PetscOptionsEnum("-problem",
4030c59ef15SJeremy L Thompson                           "CEED benchmark problem to solve", NULL,
4040c59ef15SJeremy L Thompson                           bpTypes, (PetscEnum)bpChoice, (PetscEnum*)&bpChoice,
4050c59ef15SJeremy L Thompson                           NULL); CHKERRQ(ierr);
4060c59ef15SJeremy L Thompson   vscale = bpOptions[bpChoice].vscale;
4070c59ef15SJeremy L Thompson   test_mode = PETSC_FALSE;
4080c59ef15SJeremy L Thompson   ierr = PetscOptionsBool("-test",
4090c59ef15SJeremy L Thompson                           "Testing mode (do not print unless error is large)",
4100c59ef15SJeremy L Thompson                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
4110c59ef15SJeremy L Thompson   benchmark_mode = PETSC_FALSE;
4120c59ef15SJeremy L Thompson   ierr = PetscOptionsBool("-benchmark",
4130c59ef15SJeremy L Thompson                           "Benchmarking mode (prints benchmark statistics)",
4140c59ef15SJeremy L Thompson                           NULL, benchmark_mode, &benchmark_mode, NULL);
4150c59ef15SJeremy L Thompson   CHKERRQ(ierr);
416*6bd9afcaSjeremylt   write_solution = PETSC_FALSE;
417*6bd9afcaSjeremylt   ierr = PetscOptionsBool("-write_solution",
418*6bd9afcaSjeremylt                           "Write solution for visualization",
419*6bd9afcaSjeremylt                           NULL, write_solution, &write_solution, NULL);
420*6bd9afcaSjeremylt   CHKERRQ(ierr);
4210c59ef15SJeremy L Thompson   degree = test_mode ? 3 : 1;
4220c59ef15SJeremy L Thompson   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
4230c59ef15SJeremy L Thompson                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
4240c59ef15SJeremy L Thompson   qextra = bpOptions[bpChoice].qextra;
4250c59ef15SJeremy L Thompson   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
4260c59ef15SJeremy L Thompson                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
4270c59ef15SJeremy L Thompson   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
4280c59ef15SJeremy L Thompson                             NULL, ceedresource, ceedresource,
4290c59ef15SJeremy L Thompson                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
4300c59ef15SJeremy L Thompson   localdof = 1000;
4310c59ef15SJeremy L Thompson   ierr = PetscOptionsInt("-local",
4320c59ef15SJeremy L Thompson                          "Target number of locally owned degrees of freedom per process",
4330c59ef15SJeremy L Thompson                          NULL, localdof, &localdof, NULL); CHKERRQ(ierr);
4340c59ef15SJeremy L Thompson   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
4353ce86546SJeremy L Thompson   P = degree + 1;
4363ce86546SJeremy L Thompson   Q = P + qextra;
4370c59ef15SJeremy L Thompson 
4380c59ef15SJeremy L Thompson   // Determine size of process grid
4390c59ef15SJeremy L Thompson   ierr = MPI_Comm_size(comm, &size); CHKERRQ(ierr);
4400c59ef15SJeremy L Thompson   Split3(size, p, false);
4410c59ef15SJeremy L Thompson 
4420c59ef15SJeremy L Thompson   // Find a nicely composite number of elements no less than localdof
4430c59ef15SJeremy L Thompson   for (localelem = PetscMax(1, localdof / (degree*degree*degree)); ;
4440c59ef15SJeremy L Thompson        localelem++) {
4450c59ef15SJeremy L Thompson     Split3(localelem, melem, true);
4460c59ef15SJeremy L Thompson     if (Max3(melem) / Min3(melem) <= 2) break;
4470c59ef15SJeremy L Thompson   }
4480c59ef15SJeremy L Thompson 
4490c59ef15SJeremy L Thompson   // Find my location in the process grid
4500c59ef15SJeremy L Thompson   ierr = MPI_Comm_rank(comm, &rank); CHKERRQ(ierr);
4510c59ef15SJeremy L Thompson   for (int d=0,rankleft=rank; d<3; d++) {
4520c59ef15SJeremy L Thompson     const int pstride[3] = {p[1] *p[2], p[2], 1};
4530c59ef15SJeremy L Thompson     irank[d] = rankleft / pstride[d];
4540c59ef15SJeremy L Thompson     rankleft -= irank[d] * pstride[d];
4550c59ef15SJeremy L Thompson   }
4560c59ef15SJeremy L Thompson 
4570c59ef15SJeremy L Thompson   GlobalDof(p, irank, degree, melem, mdof);
4580c59ef15SJeremy L Thompson 
459819eb1b3Sjeremylt   // Setup global vector
4600c59ef15SJeremy L Thompson   ierr = VecCreate(comm, &X); CHKERRQ(ierr);
4610c59ef15SJeremy L Thompson   ierr = VecSetSizes(X, mdof[0]*mdof[1]*mdof[2]*vscale, PETSC_DECIDE);
4620c59ef15SJeremy L Thompson   CHKERRQ(ierr);
4630c59ef15SJeremy L Thompson   ierr = VecSetUp(X); CHKERRQ(ierr);
4640c59ef15SJeremy L Thompson 
465819eb1b3Sjeremylt   // Print summary
4660c59ef15SJeremy L Thompson   if (!test_mode) {
4670c59ef15SJeremy L Thompson     CeedInt gsize;
4680c59ef15SJeremy L Thompson     ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
4693ce86546SJeremy L Thompson     ierr = PetscPrintf(comm,
4703ce86546SJeremy L Thompson                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n"
4713ce86546SJeremy L Thompson                        "  libCEED:\n"
4723ce86546SJeremy L Thompson                        "    libCEED Backend                    : %s\n"
4733ce86546SJeremy L Thompson                        "  Mesh:\n"
4743ce86546SJeremy L Thompson                        "    Number of 1D Basis Nodes (p)       : %d\n"
4753ce86546SJeremy L Thompson                        "    Number of 1D Quadrature Points (q) : %d\n"
4763ce86546SJeremy L Thompson                        "    Global DOFs                        : %D\n"
4773ce86546SJeremy L Thompson                        "    Process Decomposition              : %D %D %D\n"
4783ce86546SJeremy L Thompson                        "    Local Elements                     : %D = %D %D %D\n"
4793ce86546SJeremy L Thompson                        "    Owned DOFs                         : %D = %D %D %D\n",
4803ce86546SJeremy L Thompson                        bpChoice+1, ceedresource, P, Q,  gsize/vscale, p[0],
4813ce86546SJeremy L Thompson                        p[1], p[2], localelem, melem[0], melem[1], melem[2],
4823ce86546SJeremy L Thompson                        mdof[0]*mdof[1]*mdof[2], mdof[0], mdof[1], mdof[2]);
4833ce86546SJeremy L Thompson     CHKERRQ(ierr);
4840c59ef15SJeremy L Thompson   }
4850c59ef15SJeremy L Thompson 
4860c59ef15SJeremy L Thompson   {
4870c59ef15SJeremy L Thompson     lsize = 1;
4880c59ef15SJeremy L Thompson     for (int d=0; d<3; d++) {
4890c59ef15SJeremy L Thompson       ldof[d] = melem[d]*degree + 1;
4900c59ef15SJeremy L Thompson       lsize *= ldof[d];
4910c59ef15SJeremy L Thompson     }
4920c59ef15SJeremy L Thompson     ierr = VecCreate(PETSC_COMM_SELF, &Xloc); CHKERRQ(ierr);
4930c59ef15SJeremy L Thompson     ierr = VecSetSizes(Xloc, lsize*vscale, PETSC_DECIDE); CHKERRQ(ierr);
4940c59ef15SJeremy L Thompson     ierr = VecSetUp(Xloc); CHKERRQ(ierr);
4950c59ef15SJeremy L Thompson 
4960c59ef15SJeremy L Thompson     // Create local-to-global scatter
4970c59ef15SJeremy L Thompson     PetscInt *ltogind, *ltogind0, *locind, l0count;
4980c59ef15SJeremy L Thompson     IS ltogis, ltogis0, locis;
4990c59ef15SJeremy L Thompson     PetscInt gstart[2][2][2], gmdof[2][2][2][3];
5000c59ef15SJeremy L Thompson 
5010c59ef15SJeremy L Thompson     for (int i=0; i<2; i++) {
5020c59ef15SJeremy L Thompson       for (int j=0; j<2; j++) {
5030c59ef15SJeremy L Thompson         for (int k=0; k<2; k++) {
5040c59ef15SJeremy L Thompson           PetscInt ijkrank[3] = {irank[0]+i, irank[1]+j, irank[2]+k};
5050c59ef15SJeremy L Thompson           gstart[i][j][k] = GlobalStart(p, ijkrank, degree, melem);
5060c59ef15SJeremy L Thompson           GlobalDof(p, ijkrank, degree, melem, gmdof[i][j][k]);
5070c59ef15SJeremy L Thompson         }
5080c59ef15SJeremy L Thompson       }
5090c59ef15SJeremy L Thompson     }
5100c59ef15SJeremy L Thompson 
5110c59ef15SJeremy L Thompson     ierr = PetscMalloc1(lsize, &ltogind); CHKERRQ(ierr);
5120c59ef15SJeremy L Thompson     ierr = PetscMalloc1(lsize, &ltogind0); CHKERRQ(ierr);
5130c59ef15SJeremy L Thompson     ierr = PetscMalloc1(lsize, &locind); CHKERRQ(ierr);
5140c59ef15SJeremy L Thompson     l0count = 0;
5150c59ef15SJeremy L Thompson     for (PetscInt i=0,ir,ii; ir=i>=mdof[0], ii=i-ir*mdof[0], i<ldof[0]; i++) {
5160c59ef15SJeremy L Thompson       for (PetscInt j=0,jr,jj; jr=j>=mdof[1], jj=j-jr*mdof[1], j<ldof[1]; j++) {
5170c59ef15SJeremy L Thompson         for (PetscInt k=0,kr,kk; kr=k>=mdof[2], kk=k-kr*mdof[2], k<ldof[2]; k++) {
5180c59ef15SJeremy L Thompson           PetscInt here = (i*ldof[1]+j)*ldof[2]+k;
5190c59ef15SJeremy L Thompson           ltogind[here] =
5200c59ef15SJeremy L Thompson             gstart[ir][jr][kr] + (ii*gmdof[ir][jr][kr][1]+jj)*gmdof[ir][jr][kr][2]+kk;
5210c59ef15SJeremy L Thompson           if ((irank[0] == 0 && i == 0)
5220c59ef15SJeremy L Thompson               || (irank[1] == 0 && j == 0)
5230c59ef15SJeremy L Thompson               || (irank[2] == 0 && k == 0)
5240c59ef15SJeremy L Thompson               || (irank[0]+1 == p[0] && i+1 == ldof[0])
5250c59ef15SJeremy L Thompson               || (irank[1]+1 == p[1] && j+1 == ldof[1])
5260c59ef15SJeremy L Thompson               || (irank[2]+1 == p[2] && k+1 == ldof[2]))
5270c59ef15SJeremy L Thompson             continue;
5280c59ef15SJeremy L Thompson           ltogind0[l0count] = ltogind[here];
5290c59ef15SJeremy L Thompson           locind[l0count++] = here;
5300c59ef15SJeremy L Thompson         }
5310c59ef15SJeremy L Thompson       }
5320c59ef15SJeremy L Thompson     }
5330c59ef15SJeremy L Thompson     ierr = ISCreateBlock(comm, vscale, lsize, ltogind, PETSC_OWN_POINTER,
5340c59ef15SJeremy L Thompson                          &ltogis); CHKERRQ(ierr);
5350c59ef15SJeremy L Thompson     ierr = VecScatterCreate(Xloc, NULL, X, ltogis, &ltog); CHKERRQ(ierr);
5360c59ef15SJeremy L Thompson     CHKERRQ(ierr);
5370c59ef15SJeremy L Thompson     ierr = ISCreateBlock(comm, vscale, l0count, ltogind0, PETSC_OWN_POINTER,
5380c59ef15SJeremy L Thompson                          &ltogis0); CHKERRQ(ierr);
5390c59ef15SJeremy L Thompson     ierr = ISCreateBlock(comm, vscale, l0count, locind, PETSC_OWN_POINTER,
5400c59ef15SJeremy L Thompson                          &locis); CHKERRQ(ierr);
5410c59ef15SJeremy L Thompson     ierr = VecScatterCreate(Xloc, locis, X, ltogis0, &ltog0); CHKERRQ(ierr);
5420c59ef15SJeremy L Thompson     {
5430c59ef15SJeremy L Thompson       // Create global-to-global scatter for Dirichlet values (everything not in
5440c59ef15SJeremy L Thompson       // ltogis0, which is the range of ltog0)
5450c59ef15SJeremy L Thompson       PetscInt xstart, xend, *indD, countD = 0;
5460c59ef15SJeremy L Thompson       IS isD;
5470c59ef15SJeremy L Thompson       const PetscScalar *x;
5480c59ef15SJeremy L Thompson       ierr = VecZeroEntries(Xloc); CHKERRQ(ierr);
5490c59ef15SJeremy L Thompson       ierr = VecSet(X, 1.0); CHKERRQ(ierr);
5500c59ef15SJeremy L Thompson       ierr = VecScatterBegin(ltog0, Xloc, X, INSERT_VALUES, SCATTER_FORWARD);
5510c59ef15SJeremy L Thompson       CHKERRQ(ierr);
5520c59ef15SJeremy L Thompson       ierr = VecScatterEnd(ltog0, Xloc, X, INSERT_VALUES, SCATTER_FORWARD);
5530c59ef15SJeremy L Thompson       CHKERRQ(ierr);
5540c59ef15SJeremy L Thompson       ierr = VecGetOwnershipRange(X, &xstart, &xend); CHKERRQ(ierr);
5550c59ef15SJeremy L Thompson       ierr = PetscMalloc1(xend-xstart, &indD); CHKERRQ(ierr);
5560c59ef15SJeremy L Thompson       ierr = VecGetArrayRead(X, &x); CHKERRQ(ierr);
5570c59ef15SJeremy L Thompson       for (PetscInt i=0; i<xend-xstart; i++) {
5580c59ef15SJeremy L Thompson         if (x[i] == 1.) indD[countD++] = xstart + i;
5590c59ef15SJeremy L Thompson       }
5600c59ef15SJeremy L Thompson       ierr = VecRestoreArrayRead(X, &x); CHKERRQ(ierr);
5610c59ef15SJeremy L Thompson       ierr = ISCreateGeneral(comm, countD, indD, PETSC_COPY_VALUES, &isD);
5620c59ef15SJeremy L Thompson       CHKERRQ(ierr);
5630c59ef15SJeremy L Thompson       ierr = PetscFree(indD); CHKERRQ(ierr);
5640c59ef15SJeremy L Thompson       ierr = VecScatterCreate(X, isD, X, isD, &gtogD); CHKERRQ(ierr);
5650c59ef15SJeremy L Thompson       ierr = ISDestroy(&isD); CHKERRQ(ierr);
5660c59ef15SJeremy L Thompson     }
5670c59ef15SJeremy L Thompson     ierr = ISDestroy(&ltogis); CHKERRQ(ierr);
5680c59ef15SJeremy L Thompson     ierr = ISDestroy(&ltogis0); CHKERRQ(ierr);
5690c59ef15SJeremy L Thompson     ierr = ISDestroy(&locis); CHKERRQ(ierr);
5700c59ef15SJeremy L Thompson   }
5710c59ef15SJeremy L Thompson 
5720c59ef15SJeremy L Thompson   // Set up libCEED
5730c59ef15SJeremy L Thompson   CeedInit(ceedresource, &ceed);
574e07e9ddcSjeremylt   CeedBasisCreateTensorH1Lagrange(ceed, 3, vscale, P, Q,
575e07e9ddcSjeremylt                                   bpOptions[bpChoice].qmode, &basisu);
576e07e9ddcSjeremylt   CeedBasisCreateTensorH1Lagrange(ceed, 3, 3, 2, Q,
577e07e9ddcSjeremylt                                   bpOptions[bpChoice].qmode, &basisx);
5780c59ef15SJeremy L Thompson 
5790c59ef15SJeremy L Thompson   CreateRestriction(ceed, melem, P, vscale, &Erestrictu);
5800c59ef15SJeremy L Thompson   CreateRestriction(ceed, melem, 2, 3, &Erestrictx);
5810c59ef15SJeremy L Thompson   CeedInt nelem = melem[0]*melem[1]*melem[2];
5820c59ef15SJeremy L Thompson   CeedElemRestrictionCreateIdentity(ceed, nelem, Q*Q*Q, nelem*Q*Q*Q, vscale,
5830c59ef15SJeremy L Thompson                                     &Erestrictui);
5840c59ef15SJeremy L Thompson   CeedElemRestrictionCreateIdentity(ceed, nelem,
5850c59ef15SJeremy L Thompson                                     bpOptions[bpChoice].qdatasize*Q*Q*Q,
5860c59ef15SJeremy L Thompson                                     bpOptions[bpChoice].qdatasize*nelem*Q*Q*Q,
5870c59ef15SJeremy L Thompson                                     1, &Erestrictqdi);
5880c59ef15SJeremy L Thompson   CeedElemRestrictionCreateIdentity(ceed, nelem, Q*Q*Q, nelem*Q*Q*Q, 1,
5890c59ef15SJeremy L Thompson                                     &Erestrictxi);
5900c59ef15SJeremy L Thompson   {
5910c59ef15SJeremy L Thompson     CeedScalar *xloc;
5920c59ef15SJeremy L Thompson     CeedInt shape[3] = {melem[0]+1, melem[1]+1, melem[2]+1}, len =
5930c59ef15SJeremy L Thompson                          shape[0]*shape[1]*shape[2];
5940c59ef15SJeremy L Thompson     xloc = malloc(len*3*sizeof xloc[0]);
5950c59ef15SJeremy L Thompson     for (CeedInt i=0; i<shape[0]; i++) {
5960c59ef15SJeremy L Thompson       for (CeedInt j=0; j<shape[1]; j++) {
5970c59ef15SJeremy L Thompson         for (CeedInt k=0; k<shape[2]; k++) {
5980c59ef15SJeremy L Thompson           xloc[((i*shape[1]+j)*shape[2]+k) + 0*len] = 1.*(irank[0]*melem[0]+i) /
5990c59ef15SJeremy L Thompson               (p[0]*melem[0]);
6000c59ef15SJeremy L Thompson           xloc[((i*shape[1]+j)*shape[2]+k) + 1*len] = 1.*(irank[1]*melem[1]+j) /
6010c59ef15SJeremy L Thompson               (p[1]*melem[1]);
6020c59ef15SJeremy L Thompson           xloc[((i*shape[1]+j)*shape[2]+k) + 2*len] = 1.*(irank[2]*melem[2]+k) /
6030c59ef15SJeremy L Thompson               (p[2]*melem[2]);
6040c59ef15SJeremy L Thompson         }
6050c59ef15SJeremy L Thompson       }
6060c59ef15SJeremy L Thompson     }
6070c59ef15SJeremy L Thompson     CeedVectorCreate(ceed, len*3, &xcoord);
6080c59ef15SJeremy L Thompson     CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_OWN_POINTER, xloc);
6090c59ef15SJeremy L Thompson   }
6100c59ef15SJeremy L Thompson 
6110c59ef15SJeremy L Thompson   // Create the Q-function that builds the operator (i.e. computes its
6124b5b4ec1Sjeremylt   // quadrature data) and set its context data
6130c59ef15SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].setup,
6140c59ef15SJeremy L Thompson                               bpOptions[bpChoice].setupfname, &qf_setup);
6150c59ef15SJeremy L Thompson   CeedQFunctionAddInput(qf_setup, "x", 3, CEED_EVAL_INTERP);
6160c59ef15SJeremy L Thompson   CeedQFunctionAddInput(qf_setup, "dx", 3, CEED_EVAL_GRAD);
6170c59ef15SJeremy L Thompson   CeedQFunctionAddInput(qf_setup, "weight", 1, CEED_EVAL_WEIGHT);
6180c59ef15SJeremy L Thompson   CeedQFunctionAddOutput(qf_setup, "rho", bpOptions[bpChoice].qdatasize,
6190c59ef15SJeremy L Thompson                          CEED_EVAL_NONE);
6200c59ef15SJeremy L Thompson   CeedQFunctionAddOutput(qf_setup, "true_soln", vscale, CEED_EVAL_NONE);
6210c59ef15SJeremy L Thompson   CeedQFunctionAddOutput(qf_setup, "rhs", vscale, CEED_EVAL_INTERP);
6220c59ef15SJeremy L Thompson 
6230c59ef15SJeremy L Thompson   // Set up PDE operator
6240c59ef15SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].apply,
6250c59ef15SJeremy L Thompson                               bpOptions[bpChoice].applyfname, &qf_apply);
6260c59ef15SJeremy L Thompson   // Add inputs and outputs
6270c59ef15SJeremy L Thompson   CeedQFunctionAddInput(qf_apply, "u", vscale, bpOptions[bpChoice].inmode);
6280c59ef15SJeremy L Thompson   CeedQFunctionAddInput(qf_apply, "rho", bpOptions[bpChoice].qdatasize,
6290c59ef15SJeremy L Thompson                         CEED_EVAL_NONE);
6300c59ef15SJeremy L Thompson   CeedQFunctionAddOutput(qf_apply, "v", vscale, bpOptions[bpChoice].outmode);
6310c59ef15SJeremy L Thompson 
6320c59ef15SJeremy L Thompson   // Create the error qfunction
6330c59ef15SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].error,
6340c59ef15SJeremy L Thompson                               bpOptions[bpChoice].errorfname, &qf_error);
6350c59ef15SJeremy L Thompson   CeedQFunctionAddInput(qf_error, "u", vscale, CEED_EVAL_INTERP);
6360c59ef15SJeremy L Thompson   CeedQFunctionAddInput(qf_error, "true_soln", vscale, CEED_EVAL_NONE);
6370c59ef15SJeremy L Thompson   CeedQFunctionAddOutput(qf_error, "error", vscale, CEED_EVAL_NONE);
6380c59ef15SJeremy L Thompson 
6390c59ef15SJeremy L Thompson   // Create the persistent vectors that will be needed in setup
640819eb1b3Sjeremylt   CeedInt nqpts;
641819eb1b3Sjeremylt   CeedBasisGetNumQuadraturePoints(basisu, &nqpts);
642819eb1b3Sjeremylt   CeedVectorCreate(ceed, bpOptions[bpChoice].qdatasize*nelem*nqpts, &rho);
643819eb1b3Sjeremylt   CeedVectorCreate(ceed, nelem*nqpts*vscale, &target);
6440c59ef15SJeremy L Thompson   CeedVectorCreate(ceed, lsize*vscale, &rhsceed);
6450c59ef15SJeremy L Thompson 
6464b5b4ec1Sjeremylt   // Create the operator that builds the quadrature data for the ceed operator
6470c59ef15SJeremy L Thompson   CeedOperatorCreate(ceed, qf_setup, NULL, NULL, &op_setup);
6480c59ef15SJeremy L Thompson   CeedOperatorSetField(op_setup, "x", Erestrictx, CEED_NOTRANSPOSE,
6490c59ef15SJeremy L Thompson                        basisx, CEED_VECTOR_ACTIVE);
6500c59ef15SJeremy L Thompson   CeedOperatorSetField(op_setup, "dx", Erestrictx, CEED_NOTRANSPOSE,
6510c59ef15SJeremy L Thompson                        basisx, CEED_VECTOR_ACTIVE);
6520c59ef15SJeremy L Thompson   CeedOperatorSetField(op_setup, "weight", Erestrictxi, CEED_NOTRANSPOSE,
6530c59ef15SJeremy L Thompson                        basisx, CEED_VECTOR_NONE);
6540c59ef15SJeremy L Thompson   CeedOperatorSetField(op_setup, "rho", Erestrictqdi, CEED_NOTRANSPOSE,
6550c59ef15SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
6560c59ef15SJeremy L Thompson   CeedOperatorSetField(op_setup, "true_soln", Erestrictui, CEED_NOTRANSPOSE,
6570c59ef15SJeremy L Thompson                        CEED_BASIS_COLLOCATED, target);
6580c59ef15SJeremy L Thompson   CeedOperatorSetField(op_setup, "rhs", Erestrictu, CEED_TRANSPOSE,
6590c59ef15SJeremy L Thompson                        basisu, rhsceed);
6600c59ef15SJeremy L Thompson 
6614b5b4ec1Sjeremylt   // Create the mass or diff operator
6620c59ef15SJeremy L Thompson   CeedOperatorCreate(ceed, qf_apply, NULL, NULL, &op_apply);
6630c59ef15SJeremy L Thompson   CeedOperatorSetField(op_apply, "u", Erestrictu, CEED_TRANSPOSE,
6640c59ef15SJeremy L Thompson                        basisu, CEED_VECTOR_ACTIVE);
6650c59ef15SJeremy L Thompson   CeedOperatorSetField(op_apply, "rho", Erestrictqdi, CEED_NOTRANSPOSE,
6660c59ef15SJeremy L Thompson                        CEED_BASIS_COLLOCATED, rho);
6670c59ef15SJeremy L Thompson   CeedOperatorSetField(op_apply, "v", Erestrictu, CEED_TRANSPOSE,
6680c59ef15SJeremy L Thompson                        basisu, CEED_VECTOR_ACTIVE);
6690c59ef15SJeremy L Thompson 
6700c59ef15SJeremy L Thompson   // Create the error operator
6710c59ef15SJeremy L Thompson   CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error);
6720c59ef15SJeremy L Thompson   CeedOperatorSetField(op_error, "u", Erestrictu, CEED_TRANSPOSE,
6730c59ef15SJeremy L Thompson                        basisu, CEED_VECTOR_ACTIVE);
6740c59ef15SJeremy L Thompson   CeedOperatorSetField(op_error, "true_soln", Erestrictui, CEED_NOTRANSPOSE,
6750c59ef15SJeremy L Thompson                        CEED_BASIS_COLLOCATED, target);
6760c59ef15SJeremy L Thompson   CeedOperatorSetField(op_error, "error", Erestrictui, CEED_NOTRANSPOSE,
6770c59ef15SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
6780c59ef15SJeremy L Thompson 
6790c59ef15SJeremy L Thompson 
6800c59ef15SJeremy L Thompson   // Set up Mat
6810c59ef15SJeremy L Thompson   ierr = PetscMalloc1(1, &user); CHKERRQ(ierr);
6820c59ef15SJeremy L Thompson   user->comm = comm;
6830c59ef15SJeremy L Thompson   user->ltog = ltog;
6840c59ef15SJeremy L Thompson   if (bpChoice != CEED_BP1 && bpChoice != CEED_BP2) {
6850c59ef15SJeremy L Thompson     user->ltog0 = ltog0;
6860c59ef15SJeremy L Thompson     user->gtogD = gtogD;
6870c59ef15SJeremy L Thompson   }
6880c59ef15SJeremy L Thompson   user->Xloc = Xloc;
6890c59ef15SJeremy L Thompson   ierr = VecDuplicate(Xloc, &user->Yloc); CHKERRQ(ierr);
6900c59ef15SJeremy L Thompson   CeedVectorCreate(ceed, lsize*vscale, &user->xceed);
6910c59ef15SJeremy L Thompson   CeedVectorCreate(ceed, lsize*vscale, &user->yceed);
6920c59ef15SJeremy L Thompson   user->op = op_apply;
6930c59ef15SJeremy L Thompson   user->rho = rho;
6940c59ef15SJeremy L Thompson   user->ceed = ceed;
6950c59ef15SJeremy L Thompson 
6960c59ef15SJeremy L Thompson   ierr = MatCreateShell(comm, mdof[0]*mdof[1]*mdof[2]*vscale,
6970c59ef15SJeremy L Thompson                         mdof[0]*mdof[1]*mdof[2]*vscale,
6980c59ef15SJeremy L Thompson                         PETSC_DECIDE, PETSC_DECIDE, user, &mat); CHKERRQ(ierr);
6990c59ef15SJeremy L Thompson   if (bpChoice == CEED_BP1 || bpChoice == CEED_BP2) {
7000c59ef15SJeremy L Thompson     ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Mass);
7010c59ef15SJeremy L Thompson     CHKERRQ(ierr);
7020c59ef15SJeremy L Thompson   } else {
7030c59ef15SJeremy L Thompson     ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Diff);
7040c59ef15SJeremy L Thompson     CHKERRQ(ierr);
7050c59ef15SJeremy L Thompson   }
7060c59ef15SJeremy L Thompson   ierr = MatCreateVecs(mat, &rhs, NULL); CHKERRQ(ierr);
7070c59ef15SJeremy L Thompson 
7080c59ef15SJeremy L Thompson   // Get RHS vector
7090c59ef15SJeremy L Thompson   ierr = VecDuplicate(Xloc, &rhsloc); CHKERRQ(ierr);
7100c59ef15SJeremy L Thompson   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
7110c59ef15SJeremy L Thompson   ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr);
7120c59ef15SJeremy L Thompson   CeedVectorSetArray(rhsceed, CEED_MEM_HOST, CEED_USE_POINTER, r);
7130c59ef15SJeremy L Thompson 
7140c59ef15SJeremy L Thompson   // Setup rho, rhs, and target
7150c59ef15SJeremy L Thompson   CeedOperatorApply(op_setup, xcoord, rho, CEED_REQUEST_IMMEDIATE);
7160c59ef15SJeremy L Thompson   ierr = CeedVectorSyncArray(rhsceed, CEED_MEM_HOST); CHKERRQ(ierr);
7170c59ef15SJeremy L Thompson   CeedVectorDestroy(&xcoord);
7180c59ef15SJeremy L Thompson 
7190c59ef15SJeremy L Thompson   // Gather RHS
7200c59ef15SJeremy L Thompson   ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr);
7210c59ef15SJeremy L Thompson   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
7220c59ef15SJeremy L Thompson   ierr = VecScatterBegin(ltog, rhsloc, rhs, ADD_VALUES, SCATTER_FORWARD);
7230c59ef15SJeremy L Thompson   CHKERRQ(ierr);
7240c59ef15SJeremy L Thompson   ierr = VecScatterEnd(ltog, rhsloc, rhs, ADD_VALUES, SCATTER_FORWARD);
7250c59ef15SJeremy L Thompson   CHKERRQ(ierr);
7260c59ef15SJeremy L Thompson   CeedVectorDestroy(&rhsceed);
7270c59ef15SJeremy L Thompson 
7280c59ef15SJeremy L Thompson   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
7290c59ef15SJeremy L Thompson   {
7300c59ef15SJeremy L Thompson     PC pc;
7310c59ef15SJeremy L Thompson     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
7320c59ef15SJeremy L Thompson     if (bpChoice == CEED_BP1 || bpChoice == CEED_BP2) {
7330c59ef15SJeremy L Thompson       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
7340c59ef15SJeremy L Thompson       ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr);
7350c59ef15SJeremy L Thompson     } else {
7360c59ef15SJeremy L Thompson       ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr);
7370c59ef15SJeremy L Thompson     }
7380c59ef15SJeremy L Thompson     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
7390c59ef15SJeremy L Thompson     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
7400c59ef15SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
7410c59ef15SJeremy L Thompson                             PETSC_DEFAULT); CHKERRQ(ierr);
7420c59ef15SJeremy L Thompson   }
7430c59ef15SJeremy L Thompson   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
7440c59ef15SJeremy L Thompson   ierr = KSPSetOperators(ksp, mat, mat); CHKERRQ(ierr);
7450c59ef15SJeremy L Thompson   // First run, if benchmarking
7460c59ef15SJeremy L Thompson   if (benchmark_mode) {
7470c59ef15SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
7480c59ef15SJeremy L Thompson     CHKERRQ(ierr);
7490c59ef15SJeremy L Thompson     my_rt_start = MPI_Wtime();
7500c59ef15SJeremy L Thompson     ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
7510c59ef15SJeremy L Thompson     my_rt = MPI_Wtime() - my_rt_start;
7520c59ef15SJeremy L Thompson     // Set maxits based on first iteration timing
7530c59ef15SJeremy L Thompson     if (my_rt > 0.02) {
7540c59ef15SJeremy L Thompson       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
7550c59ef15SJeremy L Thompson       CHKERRQ(ierr);
7560c59ef15SJeremy L Thompson     } else {
7570c59ef15SJeremy L Thompson       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
7580c59ef15SJeremy L Thompson       CHKERRQ(ierr);
7590c59ef15SJeremy L Thompson     }
7600c59ef15SJeremy L Thompson   }
7610c59ef15SJeremy L Thompson   // Timed solve
7620c59ef15SJeremy L Thompson   my_rt_start = MPI_Wtime();
7630c59ef15SJeremy L Thompson   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
7640c59ef15SJeremy L Thompson   my_rt = MPI_Wtime() - my_rt_start;
7650c59ef15SJeremy L Thompson   {
7660c59ef15SJeremy L Thompson     KSPType ksptype;
7670c59ef15SJeremy L Thompson     KSPConvergedReason reason;
7680c59ef15SJeremy L Thompson     PetscReal rnorm;
7690c59ef15SJeremy L Thompson     PetscInt its;
7700c59ef15SJeremy L Thompson     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
7710c59ef15SJeremy L Thompson     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
7720c59ef15SJeremy L Thompson     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
7730c59ef15SJeremy L Thompson     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
7740c59ef15SJeremy L Thompson     if (!test_mode || reason < 0 || rnorm > 1e-8) {
7753ce86546SJeremy L Thompson       ierr = PetscPrintf(comm,
7763ce86546SJeremy L Thompson                          "  KSP:\n"
7773ce86546SJeremy L Thompson                          "    KSP Type                           : %s\n"
7783ce86546SJeremy L Thompson                          "    KSP Convergence                    : %s\n"
7793ce86546SJeremy L Thompson                          "    Total KSP Iterations               : %D\n"
7803ce86546SJeremy L Thompson                          "    Final rnorm                        : %e\n",
7813ce86546SJeremy L Thompson                          ksptype, KSPConvergedReasons[reason], its,
7823ce86546SJeremy L Thompson                          (double)rnorm); CHKERRQ(ierr);
7830c59ef15SJeremy L Thompson     }
7840c59ef15SJeremy L Thompson     if (benchmark_mode && (!test_mode)) {
7850c59ef15SJeremy L Thompson       CeedInt gsize;
7860c59ef15SJeremy L Thompson       ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
7870c59ef15SJeremy L Thompson       MPI_Reduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, 0, comm);
7880c59ef15SJeremy L Thompson       MPI_Reduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, 0, comm);
7890c59ef15SJeremy L Thompson       ierr = PetscPrintf(comm,
7903ce86546SJeremy L Thompson                          "  Performance:\n"
7913ce86546SJeremy L Thompson                          "    CG Solve Time                      : %g (%g) sec\n"
7923ce86546SJeremy L Thompson                          "    DOFs/Sec in CG                     : %g (%g) million\n",
7933ce86546SJeremy L Thompson                          rt_max, rt_min, 1e-6*gsize*its/rt_max,
7943ce86546SJeremy L Thompson                          1e-6*gsize*its/rt_min); CHKERRQ(ierr);
7950c59ef15SJeremy L Thompson     }
7960c59ef15SJeremy L Thompson   }
7970c59ef15SJeremy L Thompson 
7980c59ef15SJeremy L Thompson   {
7990c59ef15SJeremy L Thompson     PetscReal maxerror;
8000c59ef15SJeremy L Thompson     ierr = ComputeErrorMax(user, op_error, X, target, &maxerror); CHKERRQ(ierr);
8010c59ef15SJeremy L Thompson     PetscReal tol = (bpChoice == CEED_BP1 || bpChoice == CEED_BP2) ? 5e-3 : 5e-2;
8020c59ef15SJeremy L Thompson     if (!test_mode || maxerror > tol) {
8033ce86546SJeremy L Thompson       ierr = PetscPrintf(comm,
8043ce86546SJeremy L Thompson                          "    Pointwise Error (max)              : %e\n",
8053ce86546SJeremy L Thompson                          (double)maxerror); CHKERRQ(ierr);
8060c59ef15SJeremy L Thompson     }
8070c59ef15SJeremy L Thompson   }
8080c59ef15SJeremy L Thompson 
809*6bd9afcaSjeremylt   if (write_solution) {
810*6bd9afcaSjeremylt     PetscViewer vtkviewersoln;
811*6bd9afcaSjeremylt 
812*6bd9afcaSjeremylt     ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr);
813*6bd9afcaSjeremylt     ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr);
814*6bd9afcaSjeremylt     ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr);
815*6bd9afcaSjeremylt     ierr = VecView(X, vtkviewersoln); CHKERRQ(ierr);
816*6bd9afcaSjeremylt     ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr);
817*6bd9afcaSjeremylt   }
818*6bd9afcaSjeremylt 
8190c59ef15SJeremy L Thompson   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
8200c59ef15SJeremy L Thompson   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
8210c59ef15SJeremy L Thompson   ierr = VecDestroy(&X); CHKERRQ(ierr);
8220c59ef15SJeremy L Thompson   ierr = VecDestroy(&user->Xloc); CHKERRQ(ierr);
8230c59ef15SJeremy L Thompson   ierr = VecDestroy(&user->Yloc); CHKERRQ(ierr);
8240c59ef15SJeremy L Thompson   ierr = VecScatterDestroy(&ltog); CHKERRQ(ierr);
8250c59ef15SJeremy L Thompson   ierr = VecScatterDestroy(&ltog0); CHKERRQ(ierr);
8260c59ef15SJeremy L Thompson   ierr = VecScatterDestroy(&gtogD); CHKERRQ(ierr);
8270c59ef15SJeremy L Thompson   ierr = MatDestroy(&mat); CHKERRQ(ierr);
8280c59ef15SJeremy L Thompson   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
8290c59ef15SJeremy L Thompson 
8300c59ef15SJeremy L Thompson   CeedVectorDestroy(&user->xceed);
8310c59ef15SJeremy L Thompson   CeedVectorDestroy(&user->yceed);
8320c59ef15SJeremy L Thompson   CeedVectorDestroy(&user->rho);
8330c59ef15SJeremy L Thompson   CeedVectorDestroy(&target);
8340c59ef15SJeremy L Thompson   CeedOperatorDestroy(&op_setup);
8350c59ef15SJeremy L Thompson   CeedOperatorDestroy(&op_apply);
8360c59ef15SJeremy L Thompson   CeedOperatorDestroy(&op_error);
8370c59ef15SJeremy L Thompson   CeedElemRestrictionDestroy(&Erestrictu);
8380c59ef15SJeremy L Thompson   CeedElemRestrictionDestroy(&Erestrictx);
8390c59ef15SJeremy L Thompson   CeedElemRestrictionDestroy(&Erestrictui);
8400c59ef15SJeremy L Thompson   CeedElemRestrictionDestroy(&Erestrictxi);
8410c59ef15SJeremy L Thompson   CeedElemRestrictionDestroy(&Erestrictqdi);
8420c59ef15SJeremy L Thompson   CeedQFunctionDestroy(&qf_setup);
8430c59ef15SJeremy L Thompson   CeedQFunctionDestroy(&qf_apply);
8440c59ef15SJeremy L Thompson   CeedQFunctionDestroy(&qf_error);
8450c59ef15SJeremy L Thompson   CeedBasisDestroy(&basisu);
8460c59ef15SJeremy L Thompson   CeedBasisDestroy(&basisx);
8470c59ef15SJeremy L Thompson   CeedDestroy(&ceed);
8480c59ef15SJeremy L Thompson   ierr = PetscFree(user); CHKERRQ(ierr);
8490c59ef15SJeremy L Thompson   return PetscFinalize();
8500c59ef15SJeremy L Thompson }
851