xref: /libCEED/examples/petsc/bps.c (revision 2d03409ca85173e0565dadbac5c9383b99bac8d1)
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 //
389e168476Sjeremylt //TESTARGS -ceed {ceed_resource} -test -problem bp2 -degree 5 -qextra 5
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>
474d537eeaSYohann #include <petscksp.h>
484d537eeaSYohann #include <ceed.h>
496c5df90dSjeremylt #include "qfunctions/common.h"
506c5df90dSjeremylt #include "qfunctions/bp1.h"
516c5df90dSjeremylt #include "qfunctions/bp2.h"
526c5df90dSjeremylt #include "qfunctions/bp3.h"
536c5df90dSjeremylt #include "qfunctions/bp4.h"
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 }
70b52ddc88Sjeremylt static void GlobalNodes(const PetscInt p[3], const PetscInt irank[3],
710c59ef15SJeremy L Thompson                         PetscInt degree, const PetscInt melem[3],
72b52ddc88Sjeremylt                         PetscInt mnodes[3]) {
730c59ef15SJeremy L Thompson   for (int d=0; d<3; d++)
74b52ddc88Sjeremylt     mnodes[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++) {
83b52ddc88Sjeremylt         PetscInt mnodes[3], ijkrank[] = {i,j,k};
840c59ef15SJeremy L Thompson         if (i == irank[0] && j == irank[1] && k == irank[2]) return start;
85b52ddc88Sjeremylt         GlobalNodes(p, ijkrank, degree, melem, mnodes);
86b52ddc88Sjeremylt         start += mnodes[0] * mnodes[1] * mnodes[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];
96b52ddc88Sjeremylt   PetscInt mnodes[3], *idx, *idxp;
970c59ef15SJeremy L Thompson 
98819eb1b3Sjeremylt   // Get indicies
99b52ddc88Sjeremylt   for (int d=0; d<3; d++) mnodes[d] = melem[d]*(P-1) + 1;
100819eb1b3Sjeremylt   idxp = idx = malloc(nelem*P*P*P*sizeof idx[0]);
10142030c05Sjeremylt   for (CeedInt i=0; i<melem[0]; i++)
10242030c05Sjeremylt     for (CeedInt j=0; j<melem[1]; j++)
10342030c05Sjeremylt       for (CeedInt k=0; k<melem[2]; k++,idxp += P*P*P)
10442030c05Sjeremylt         for (CeedInt ii=0; ii<P; ii++)
10542030c05Sjeremylt           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
108b52ddc88Sjeremylt                 idxp[(ii*P+jj)*P+kk] = (((i*(P-1)+ii)*mnodes[1]
109b52ddc88Sjeremylt                                          + (j*(P-1)+jj))*mnodes[2]
1100c59ef15SJeremy L Thompson                                         + (k*(P-1)+kk));
1110c59ef15SJeremy L Thompson               } else { // (k,j,i) ordering for consistency with MFEM example
112b52ddc88Sjeremylt                 idxp[ii+P*(jj+P*kk)] = (((i*(P-1)+ii)*mnodes[1]
113b52ddc88Sjeremylt                                          + (j*(P-1)+jj))*mnodes[2]
1140c59ef15SJeremy L Thompson                                         + (k*(P-1)+kk));
1150c59ef15SJeremy L Thompson               }
1160c59ef15SJeremy L Thompson             }
117819eb1b3Sjeremylt 
118819eb1b3Sjeremylt   // Setup CEED restriction
119288c0443SJeremy L Thompson   CeedElemRestrictionCreate(ceed, nelem, P*P*P, mnodes[0]*mnodes[1]*mnodes[2],
12042030c05Sjeremylt                             ncomp, CEED_MEM_HOST, CEED_OWN_POINTER, idx,
12142030c05Sjeremylt                             Erestrict);
122819eb1b3Sjeremylt 
1230c59ef15SJeremy L Thompson   PetscFunctionReturn(0);
1240c59ef15SJeremy L Thompson }
1250c59ef15SJeremy L Thompson 
1260c59ef15SJeremy L Thompson // Data for PETSc
1270c59ef15SJeremy L Thompson typedef struct User_ *User;
1280c59ef15SJeremy L Thompson struct User_ {
1290c59ef15SJeremy L Thompson   MPI_Comm comm;
1300c59ef15SJeremy L Thompson   VecScatter ltog;              // Scatter for all entries
1310c59ef15SJeremy L Thompson   VecScatter ltog0;             // Skip Dirichlet values
1320c59ef15SJeremy L Thompson   VecScatter gtogD;             // global-to-global; only Dirichlet values
1330c59ef15SJeremy L Thompson   Vec Xloc, Yloc;
1340c59ef15SJeremy L Thompson   CeedVector xceed, yceed;
1350c59ef15SJeremy L Thompson   CeedOperator op;
136a2fa7910SValeria Barra   CeedVector qdata;
1370c59ef15SJeremy L Thompson   Ceed ceed;
1380c59ef15SJeremy L Thompson };
1390c59ef15SJeremy L Thompson 
1400c59ef15SJeremy L Thompson // BP Options
1410c59ef15SJeremy L Thompson typedef enum {
1420c59ef15SJeremy L Thompson   CEED_BP1 = 0, CEED_BP2 = 1, CEED_BP3 = 2,
1430c59ef15SJeremy L Thompson   CEED_BP4 = 3, CEED_BP5 = 4, CEED_BP6 = 5
1440c59ef15SJeremy L Thompson } bpType;
1450c59ef15SJeremy L Thompson static const char *const bpTypes[] = {"bp1","bp2","bp3","bp4","bp5","bp6",
1464d537eeaSYohann                                       "bpType","CEED_BP",0
1474d537eeaSYohann                                      };
1480c59ef15SJeremy L Thompson 
1490c59ef15SJeremy L Thompson // BP specific data
1500c59ef15SJeremy L Thompson typedef struct {
15134d77899SValeria Barra   CeedInt ncompu, qdatasize, qextra;
1526c5df90dSjeremylt   CeedQFunctionUser setupgeo, setuprhs, apply, error;
1536c5df90dSjeremylt   const char *setupgeofname, *setuprhsfname, *applyfname, *errorfname;
1540c59ef15SJeremy L Thompson   CeedEvalMode inmode, outmode;
155e07e9ddcSjeremylt   CeedQuadMode qmode;
1560c59ef15SJeremy L Thompson } bpData;
1570c59ef15SJeremy L Thompson 
1580c59ef15SJeremy L Thompson bpData bpOptions[6] = {
1590c59ef15SJeremy L Thompson   [CEED_BP1] = {
16034d77899SValeria Barra     .ncompu = 1,
1610c59ef15SJeremy L Thompson     .qdatasize = 1,
162194c25f7Sjeremylt     .qextra = 1,
1636c5df90dSjeremylt     .setupgeo = SetupMassGeo,
1646c5df90dSjeremylt     .setuprhs = SetupMassRhs,
1650c59ef15SJeremy L Thompson     .apply = Mass,
1660c59ef15SJeremy L Thompson     .error = Error,
1676c5df90dSjeremylt     .setupgeofname = SetupMassGeo_loc,
1686c5df90dSjeremylt     .setuprhsfname = SetupMassRhs_loc,
1694d537eeaSYohann     .applyfname = Mass_loc,
1704d537eeaSYohann     .errorfname = Error_loc,
1710c59ef15SJeremy L Thompson     .inmode = CEED_EVAL_INTERP,
172e07e9ddcSjeremylt     .outmode = CEED_EVAL_INTERP,
173f90c8643Sjeremylt     .qmode = CEED_GAUSS
174f90c8643Sjeremylt   },
1750c59ef15SJeremy L Thompson   [CEED_BP2] = {
17634d77899SValeria Barra     .ncompu = 3,
1770c59ef15SJeremy L Thompson     .qdatasize = 1,
178194c25f7Sjeremylt     .qextra = 1,
1796c5df90dSjeremylt     .setupgeo = SetupMassGeo,
1806c5df90dSjeremylt     .setuprhs = SetupMassRhs3,
1810c59ef15SJeremy L Thompson     .apply = Mass3,
1820c59ef15SJeremy L Thompson     .error = Error3,
1836c5df90dSjeremylt     .setupgeofname = SetupMassGeo_loc,
1846c5df90dSjeremylt     .setuprhsfname = SetupMassRhs3_loc,
1854d537eeaSYohann     .applyfname = Mass3_loc,
1864d537eeaSYohann     .errorfname = Error3_loc,
1870c59ef15SJeremy L Thompson     .inmode = CEED_EVAL_INTERP,
188e07e9ddcSjeremylt     .outmode = CEED_EVAL_INTERP,
189f90c8643Sjeremylt     .qmode = CEED_GAUSS
190f90c8643Sjeremylt   },
1910c59ef15SJeremy L Thompson   [CEED_BP3] = {
19234d77899SValeria Barra     .ncompu = 1,
1930c59ef15SJeremy L Thompson     .qdatasize = 6,
194194c25f7Sjeremylt     .qextra = 1,
1956c5df90dSjeremylt     .setupgeo = SetupDiffGeo,
1966c5df90dSjeremylt     .setuprhs = SetupDiffRhs,
1970c59ef15SJeremy L Thompson     .apply = Diff,
1980c59ef15SJeremy L Thompson     .error = Error,
1996c5df90dSjeremylt     .setupgeofname = SetupDiffGeo_loc,
2006c5df90dSjeremylt     .setuprhsfname = SetupDiffRhs_loc,
2014d537eeaSYohann     .applyfname = Diff_loc,
2024d537eeaSYohann     .errorfname = Error_loc,
2030c59ef15SJeremy L Thompson     .inmode = CEED_EVAL_GRAD,
204e07e9ddcSjeremylt     .outmode = CEED_EVAL_GRAD,
205f90c8643Sjeremylt     .qmode = CEED_GAUSS
206f90c8643Sjeremylt   },
2070c59ef15SJeremy L Thompson   [CEED_BP4] = {
20834d77899SValeria Barra     .ncompu = 3,
2090c59ef15SJeremy L Thompson     .qdatasize = 6,
210194c25f7Sjeremylt     .qextra = 1,
2116c5df90dSjeremylt     .setupgeo = SetupDiffGeo,
2126c5df90dSjeremylt     .setuprhs = SetupDiffRhs3,
2130c59ef15SJeremy L Thompson     .apply = Diff3,
2144b5b4ec1Sjeremylt     .error = Error3,
2156c5df90dSjeremylt     .setupgeofname = SetupDiffGeo_loc,
2166c5df90dSjeremylt     .setuprhsfname = SetupDiffRhs3_loc,
2176c5df90dSjeremylt     .applyfname = Diff_loc,
2184d537eeaSYohann     .errorfname = Error3_loc,
2190c59ef15SJeremy L Thompson     .inmode = CEED_EVAL_GRAD,
220e07e9ddcSjeremylt     .outmode = CEED_EVAL_GRAD,
221f90c8643Sjeremylt     .qmode = CEED_GAUSS
222f90c8643Sjeremylt   },
2230c59ef15SJeremy L Thompson   [CEED_BP5] = {
22434d77899SValeria Barra     .ncompu = 1,
2250c59ef15SJeremy L Thompson     .qdatasize = 6,
226194c25f7Sjeremylt     .qextra = 0,
2276c5df90dSjeremylt     .setupgeo = SetupDiffGeo,
2286c5df90dSjeremylt     .setuprhs = SetupDiffRhs,
2290c59ef15SJeremy L Thompson     .apply = Diff,
2300c59ef15SJeremy L Thompson     .error = Error,
2316c5df90dSjeremylt     .setupgeofname = SetupDiffGeo_loc,
2326c5df90dSjeremylt     .setuprhsfname = SetupDiffRhs_loc,
2334d537eeaSYohann     .applyfname = Diff_loc,
2344d537eeaSYohann     .errorfname = Error_loc,
2350c59ef15SJeremy L Thompson     .inmode = CEED_EVAL_GRAD,
236e07e9ddcSjeremylt     .outmode = CEED_EVAL_GRAD,
237f90c8643Sjeremylt     .qmode = CEED_GAUSS_LOBATTO
238f90c8643Sjeremylt   },
2390c59ef15SJeremy L Thompson   [CEED_BP6] = {
24034d77899SValeria Barra     .ncompu = 3,
2410c59ef15SJeremy L Thompson     .qdatasize = 6,
242194c25f7Sjeremylt     .qextra = 0,
2436c5df90dSjeremylt     .setupgeo = SetupDiffGeo,
2446c5df90dSjeremylt     .setuprhs = SetupDiffRhs3,
2450c59ef15SJeremy L Thompson     .apply = Diff3,
2464b5b4ec1Sjeremylt     .error = Error3,
2476c5df90dSjeremylt     .setupgeofname = SetupDiffGeo_loc,
2486c5df90dSjeremylt     .setuprhsfname = SetupDiffRhs3_loc,
2496c5df90dSjeremylt     .applyfname = Diff_loc,
2504d537eeaSYohann     .errorfname = Error3_loc,
2510c59ef15SJeremy L Thompson     .inmode = CEED_EVAL_GRAD,
252e07e9ddcSjeremylt     .outmode = CEED_EVAL_GRAD,
2537f823360Sjeremylt     .qmode = CEED_GAUSS_LOBATTO
2547f823360Sjeremylt   }
2550c59ef15SJeremy L Thompson };
2560c59ef15SJeremy L Thompson 
2570c59ef15SJeremy L Thompson // This function uses libCEED to compute the action of the mass matrix
2580c59ef15SJeremy L Thompson static PetscErrorCode MatMult_Mass(Mat A, Vec X, Vec Y) {
2590c59ef15SJeremy L Thompson   PetscErrorCode ierr;
2600c59ef15SJeremy L Thompson   User user;
2610c59ef15SJeremy L Thompson   PetscScalar *x, *y;
2620c59ef15SJeremy L Thompson 
2630c59ef15SJeremy L Thompson   PetscFunctionBeginUser;
2640c59ef15SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
2650c59ef15SJeremy L Thompson   ierr = VecScatterBegin(user->ltog, X, user->Xloc, INSERT_VALUES,
2660c59ef15SJeremy L Thompson                          SCATTER_REVERSE); CHKERRQ(ierr);
2670c59ef15SJeremy L Thompson   ierr = VecScatterEnd(user->ltog, X, user->Xloc, INSERT_VALUES, SCATTER_REVERSE);
2680c59ef15SJeremy L Thompson   CHKERRQ(ierr);
2690c59ef15SJeremy L Thompson   ierr = VecZeroEntries(user->Yloc); CHKERRQ(ierr);
2700c59ef15SJeremy L Thompson 
2710c59ef15SJeremy L Thompson   ierr = VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr);
2720c59ef15SJeremy L Thompson   ierr = VecGetArray(user->Yloc, &y); CHKERRQ(ierr);
2730c59ef15SJeremy L Thompson   CeedVectorSetArray(user->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
2740c59ef15SJeremy L Thompson   CeedVectorSetArray(user->yceed, CEED_MEM_HOST, CEED_USE_POINTER, y);
2750c59ef15SJeremy L Thompson 
2760c59ef15SJeremy L Thompson   CeedOperatorApply(user->op, user->xceed, user->yceed,
2770c59ef15SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
2780c59ef15SJeremy L Thompson   ierr = CeedVectorSyncArray(user->yceed, CEED_MEM_HOST); CHKERRQ(ierr);
2790c59ef15SJeremy L Thompson 
2800c59ef15SJeremy L Thompson   ierr = VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr);
2810c59ef15SJeremy L Thompson   ierr = VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr);
2820c59ef15SJeremy L Thompson 
2830c59ef15SJeremy L Thompson   if (Y) {
2840c59ef15SJeremy L Thompson     ierr = VecZeroEntries(Y); CHKERRQ(ierr);
2850c59ef15SJeremy L Thompson     ierr = VecScatterBegin(user->ltog, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD);
2860c59ef15SJeremy L Thompson     CHKERRQ(ierr);
2870c59ef15SJeremy L Thompson     ierr = VecScatterEnd(user->ltog, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD);
2880c59ef15SJeremy L Thompson     CHKERRQ(ierr);
2890c59ef15SJeremy L Thompson   }
2900c59ef15SJeremy L Thompson   PetscFunctionReturn(0);
2910c59ef15SJeremy L Thompson }
2920c59ef15SJeremy L Thompson 
2930c59ef15SJeremy L Thompson // This function uses libCEED to compute the action of the Laplacian with
2940c59ef15SJeremy L Thompson // Dirichlet boundary conditions
2950c59ef15SJeremy L Thompson static PetscErrorCode MatMult_Diff(Mat A, Vec X, Vec Y) {
2960c59ef15SJeremy L Thompson   PetscErrorCode ierr;
2970c59ef15SJeremy L Thompson   User user;
2980c59ef15SJeremy L Thompson   PetscScalar *x, *y;
2990c59ef15SJeremy L Thompson 
3000c59ef15SJeremy L Thompson   PetscFunctionBeginUser;
3010c59ef15SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
302819eb1b3Sjeremylt 
303819eb1b3Sjeremylt   // Global-to-local
3040c59ef15SJeremy L Thompson   ierr = VecScatterBegin(user->ltog0, X, user->Xloc, INSERT_VALUES,
3050c59ef15SJeremy L Thompson                          SCATTER_REVERSE); CHKERRQ(ierr);
3060c59ef15SJeremy L Thompson   ierr = VecScatterEnd(user->ltog0, X, user->Xloc, INSERT_VALUES,
3070c59ef15SJeremy L Thompson                        SCATTER_REVERSE);
3080c59ef15SJeremy L Thompson   CHKERRQ(ierr);
3090c59ef15SJeremy L Thompson   ierr = VecZeroEntries(user->Yloc); CHKERRQ(ierr);
3100c59ef15SJeremy L Thompson 
311819eb1b3Sjeremylt   // Setup CEED vectors
3120c59ef15SJeremy L Thompson   ierr = VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr);
3130c59ef15SJeremy L Thompson   ierr = VecGetArray(user->Yloc, &y); CHKERRQ(ierr);
3140c59ef15SJeremy L Thompson   CeedVectorSetArray(user->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
3150c59ef15SJeremy L Thompson   CeedVectorSetArray(user->yceed, CEED_MEM_HOST, CEED_USE_POINTER, y);
3160c59ef15SJeremy L Thompson 
317819eb1b3Sjeremylt   // Apply CEED operator
3180c59ef15SJeremy L Thompson   CeedOperatorApply(user->op, user->xceed, user->yceed,
3190c59ef15SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
3200c59ef15SJeremy L Thompson   ierr = CeedVectorSyncArray(user->yceed, CEED_MEM_HOST); CHKERRQ(ierr);
3210c59ef15SJeremy L Thompson 
322819eb1b3Sjeremylt   // Restore PETSc vectors
3230c59ef15SJeremy L Thompson   ierr = VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr);
3240c59ef15SJeremy L Thompson   ierr = VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr);
3250c59ef15SJeremy L Thompson 
326819eb1b3Sjeremylt   // Local-to-global
3270c59ef15SJeremy L Thompson   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
3280c59ef15SJeremy L Thompson   ierr = VecScatterBegin(user->gtogD, X, Y, INSERT_VALUES, SCATTER_FORWARD);
3290c59ef15SJeremy L Thompson   CHKERRQ(ierr);
3300c59ef15SJeremy L Thompson   ierr = VecScatterEnd(user->gtogD, X, Y, INSERT_VALUES, SCATTER_FORWARD);
3310c59ef15SJeremy L Thompson   CHKERRQ(ierr);
3320c59ef15SJeremy L Thompson   ierr = VecScatterBegin(user->ltog0, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD);
3330c59ef15SJeremy L Thompson   CHKERRQ(ierr);
3340c59ef15SJeremy L Thompson   ierr = VecScatterEnd(user->ltog0, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD);
3350c59ef15SJeremy L Thompson   CHKERRQ(ierr);
336819eb1b3Sjeremylt 
3370c59ef15SJeremy L Thompson   PetscFunctionReturn(0);
3380c59ef15SJeremy L Thompson }
3390c59ef15SJeremy L Thompson 
3400c59ef15SJeremy L Thompson // This function calculates the error in the final solution
3410c59ef15SJeremy L Thompson static PetscErrorCode ComputeErrorMax(User user, CeedOperator op_error, Vec X,
3420c59ef15SJeremy L Thompson                                       CeedVector target, PetscReal *maxerror) {
3430c59ef15SJeremy L Thompson   PetscErrorCode ierr;
3440c59ef15SJeremy L Thompson   PetscScalar *x;
3450c59ef15SJeremy L Thompson   CeedVector collocated_error;
3460c59ef15SJeremy L Thompson   CeedInt length;
3470c59ef15SJeremy L Thompson 
3480c59ef15SJeremy L Thompson   PetscFunctionBeginUser;
3490c59ef15SJeremy L Thompson   CeedVectorGetLength(target, &length);
3500c59ef15SJeremy L Thompson   CeedVectorCreate(user->ceed, length, &collocated_error);
351819eb1b3Sjeremylt 
352819eb1b3Sjeremylt   // Global-to-local
3530c59ef15SJeremy L Thompson   ierr = VecScatterBegin(user->ltog, X, user->Xloc, INSERT_VALUES,
3540c59ef15SJeremy L Thompson                          SCATTER_REVERSE); CHKERRQ(ierr);
3550c59ef15SJeremy L Thompson   ierr = VecScatterEnd(user->ltog, X, user->Xloc, INSERT_VALUES, SCATTER_REVERSE);
3560c59ef15SJeremy L Thompson   CHKERRQ(ierr);
357819eb1b3Sjeremylt 
358819eb1b3Sjeremylt   // Setup CEED vector
3590c59ef15SJeremy L Thompson   ierr = VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr);
3600c59ef15SJeremy L Thompson   CeedVectorSetArray(user->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
361819eb1b3Sjeremylt 
362819eb1b3Sjeremylt   // Apply CEED operator
3630c59ef15SJeremy L Thompson   CeedOperatorApply(op_error, user->xceed, collocated_error,
3640c59ef15SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
365819eb1b3Sjeremylt 
366819eb1b3Sjeremylt   // Restore PETSc vector
3670c59ef15SJeremy L Thompson   VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr);
3680c59ef15SJeremy L Thompson 
369819eb1b3Sjeremylt   // Reduce max error
3700c59ef15SJeremy L Thompson   *maxerror = 0;
3710c59ef15SJeremy L Thompson   const CeedScalar *e;
3720c59ef15SJeremy L Thompson   CeedVectorGetArrayRead(collocated_error, CEED_MEM_HOST, &e);
3730c59ef15SJeremy L Thompson   for (CeedInt i=0; i<length; i++) {
3740c59ef15SJeremy L Thompson     *maxerror = PetscMax(*maxerror, PetscAbsScalar(e[i]));
3750c59ef15SJeremy L Thompson   }
3760c59ef15SJeremy L Thompson   CeedVectorRestoreArrayRead(collocated_error, &e);
37700742fa8SJed Brown   ierr = MPI_Allreduce(MPI_IN_PLACE, maxerror,
37800742fa8SJed Brown                        1, MPIU_REAL, MPIU_MAX, user->comm); CHKERRQ(ierr);
379819eb1b3Sjeremylt 
380819eb1b3Sjeremylt   // Cleanup
3810c59ef15SJeremy L Thompson   CeedVectorDestroy(&collocated_error);
382819eb1b3Sjeremylt 
3830c59ef15SJeremy L Thompson   PetscFunctionReturn(0);
3840c59ef15SJeremy L Thompson }
3850c59ef15SJeremy L Thompson 
3860c59ef15SJeremy L Thompson int main(int argc, char **argv) {
3870c59ef15SJeremy L Thompson   PetscInt ierr;
3880c59ef15SJeremy L Thompson   MPI_Comm comm;
3890c59ef15SJeremy L Thompson   char ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self";
390819eb1b3Sjeremylt   double my_rt_start, my_rt, rt_min, rt_max;
391b52ddc88Sjeremylt   PetscInt degree, qextra, localnodes, localelem, melem[3], mnodes[3], p[3],
39234d77899SValeria Barra            irank[3], lnodes[3], lsize, ncompu = 1;
3930c59ef15SJeremy L Thompson   PetscScalar *r;
3946bd9afcaSjeremylt   PetscBool test_mode, benchmark_mode, write_solution;
3950c59ef15SJeremy L Thompson   PetscMPIInt size, rank;
396819eb1b3Sjeremylt   Vec X, Xloc, rhs, rhsloc;
397819eb1b3Sjeremylt   Mat mat;
398819eb1b3Sjeremylt   KSP ksp;
3994b5b4ec1Sjeremylt   VecScatter ltog, ltog0, gtogD;
400819eb1b3Sjeremylt   User user;
4010c59ef15SJeremy L Thompson   Ceed ceed;
4020c59ef15SJeremy L Thompson   CeedBasis basisx, basisu;
4030c59ef15SJeremy L Thompson   CeedElemRestriction Erestrictx, Erestrictu, Erestrictxi, Erestrictui,
4040c59ef15SJeremy L Thompson                       Erestrictqdi;
4056c5df90dSjeremylt   CeedQFunction qf_setupgeo, qf_setuprhs, qf_apply, qf_error;
4066c5df90dSjeremylt   CeedOperator op_setupgeo, op_setuprhs, op_apply, op_error;
407a2fa7910SValeria Barra   CeedVector xcoord, qdata, rhsceed, target;
4080c59ef15SJeremy L Thompson   CeedInt P, Q;
40934d77899SValeria Barra   const CeedInt dim = 3, ncompx = 3;
4100c59ef15SJeremy L Thompson   bpType bpChoice;
4110c59ef15SJeremy L Thompson 
4120c59ef15SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
4130c59ef15SJeremy L Thompson   if (ierr) return ierr;
4140c59ef15SJeremy L Thompson   comm = PETSC_COMM_WORLD;
4150c59ef15SJeremy L Thompson   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr);
4160c59ef15SJeremy L Thompson   bpChoice = CEED_BP1;
4170c59ef15SJeremy L Thompson   ierr = PetscOptionsEnum("-problem",
4180c59ef15SJeremy L Thompson                           "CEED benchmark problem to solve", NULL,
4190c59ef15SJeremy L Thompson                           bpTypes, (PetscEnum)bpChoice, (PetscEnum *)&bpChoice,
4200c59ef15SJeremy L Thompson                           NULL); CHKERRQ(ierr);
42134d77899SValeria Barra   ncompu = bpOptions[bpChoice].ncompu;
4220c59ef15SJeremy L Thompson   test_mode = PETSC_FALSE;
4230c59ef15SJeremy L Thompson   ierr = PetscOptionsBool("-test",
4240c59ef15SJeremy L Thompson                           "Testing mode (do not print unless error is large)",
4250c59ef15SJeremy L Thompson                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
4260c59ef15SJeremy L Thompson   benchmark_mode = PETSC_FALSE;
4270c59ef15SJeremy L Thompson   ierr = PetscOptionsBool("-benchmark",
4280c59ef15SJeremy L Thompson                           "Benchmarking mode (prints benchmark statistics)",
4290c59ef15SJeremy L Thompson                           NULL, benchmark_mode, &benchmark_mode, NULL);
4300c59ef15SJeremy L Thompson   CHKERRQ(ierr);
4316bd9afcaSjeremylt   write_solution = PETSC_FALSE;
4326bd9afcaSjeremylt   ierr = PetscOptionsBool("-write_solution",
4336bd9afcaSjeremylt                           "Write solution for visualization",
4346bd9afcaSjeremylt                           NULL, write_solution, &write_solution, NULL);
4356bd9afcaSjeremylt   CHKERRQ(ierr);
4360c59ef15SJeremy L Thompson   degree = test_mode ? 3 : 1;
4370c59ef15SJeremy L Thompson   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
4380c59ef15SJeremy L Thompson                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
4390c59ef15SJeremy L Thompson   qextra = bpOptions[bpChoice].qextra;
4400c59ef15SJeremy L Thompson   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
4410c59ef15SJeremy L Thompson                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
4420c59ef15SJeremy L Thompson   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
4430c59ef15SJeremy L Thompson                             NULL, ceedresource, ceedresource,
4440c59ef15SJeremy L Thompson                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
445b52ddc88Sjeremylt   localnodes = 1000;
4460c59ef15SJeremy L Thompson   ierr = PetscOptionsInt("-local",
447432a1099Sjeremylt                          "Target number of locally owned nodes per process",
448b52ddc88Sjeremylt                          NULL, localnodes, &localnodes, NULL); CHKERRQ(ierr);
4490c59ef15SJeremy L Thompson   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
4503ce86546SJeremy L Thompson   P = degree + 1;
4513ce86546SJeremy L Thompson   Q = P + qextra;
4520c59ef15SJeremy L Thompson 
4530c59ef15SJeremy L Thompson   // Determine size of process grid
4540c59ef15SJeremy L Thompson   ierr = MPI_Comm_size(comm, &size); CHKERRQ(ierr);
4550c59ef15SJeremy L Thompson   Split3(size, p, false);
4560c59ef15SJeremy L Thompson 
457b52ddc88Sjeremylt   // Find a nicely composite number of elements no less than localnodes
458b52ddc88Sjeremylt   for (localelem = PetscMax(1, localnodes / (degree*degree*degree)); ;
4590c59ef15SJeremy L Thompson        localelem++) {
4600c59ef15SJeremy L Thompson     Split3(localelem, melem, true);
4610c59ef15SJeremy L Thompson     if (Max3(melem) / Min3(melem) <= 2) break;
4620c59ef15SJeremy L Thompson   }
4630c59ef15SJeremy L Thompson 
4640c59ef15SJeremy L Thompson   // Find my location in the process grid
4650c59ef15SJeremy L Thompson   ierr = MPI_Comm_rank(comm, &rank); CHKERRQ(ierr);
46634d77899SValeria Barra   for (int d=0,rankleft=rank; d<dim; d++) {
4670c59ef15SJeremy L Thompson     const int pstride[3] = {p[1] *p[2], p[2], 1};
4680c59ef15SJeremy L Thompson     irank[d] = rankleft / pstride[d];
4690c59ef15SJeremy L Thompson     rankleft -= irank[d] * pstride[d];
4700c59ef15SJeremy L Thompson   }
4710c59ef15SJeremy L Thompson 
472b52ddc88Sjeremylt   GlobalNodes(p, irank, degree, melem, mnodes);
4730c59ef15SJeremy L Thompson 
474819eb1b3Sjeremylt   // Setup global vector
4750c59ef15SJeremy L Thompson   ierr = VecCreate(comm, &X); CHKERRQ(ierr);
47634d77899SValeria Barra   ierr = VecSetSizes(X, mnodes[0]*mnodes[1]*mnodes[2]*ncompu, PETSC_DECIDE);
4770c59ef15SJeremy L Thompson   CHKERRQ(ierr);
4780c59ef15SJeremy L Thompson   ierr = VecSetUp(X); CHKERRQ(ierr);
4790c59ef15SJeremy L Thompson 
480*2d03409cSjeremylt   // Set up libCEED
481*2d03409cSjeremylt   CeedInit(ceedresource, &ceed);
482*2d03409cSjeremylt 
483819eb1b3Sjeremylt   // Print summary
4840c59ef15SJeremy L Thompson   if (!test_mode) {
4850c59ef15SJeremy L Thompson     CeedInt gsize;
4860c59ef15SJeremy L Thompson     ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
487*2d03409cSjeremylt     const char *usedresource;
488*2d03409cSjeremylt     CeedGetResource(ceed, &usedresource);
4893ce86546SJeremy L Thompson     ierr = PetscPrintf(comm,
4903ce86546SJeremy L Thompson                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n"
4913ce86546SJeremy L Thompson                        "  libCEED:\n"
4923ce86546SJeremy L Thompson                        "    libCEED Backend                    : %s\n"
4933ce86546SJeremy L Thompson                        "  Mesh:\n"
4943ce86546SJeremy L Thompson                        "    Number of 1D Basis Nodes (p)       : %d\n"
4953ce86546SJeremy L Thompson                        "    Number of 1D Quadrature Points (q) : %d\n"
496b52ddc88Sjeremylt                        "    Global nodes                       : %D\n"
4973ce86546SJeremy L Thompson                        "    Process Decomposition              : %D %D %D\n"
4983ce86546SJeremy L Thompson                        "    Local Elements                     : %D = %D %D %D\n"
499b52ddc88Sjeremylt                        "    Owned nodes                        : %D = %D %D %D\n",
500*2d03409cSjeremylt                        bpChoice+1, usedresource, P, Q,  gsize/ncompu, p[0],
5013ce86546SJeremy L Thompson                        p[1], p[2], localelem, melem[0], melem[1], melem[2],
502*2d03409cSjeremylt                        mnodes[0]*mnodes[1]*mnodes[2], mnodes[0], mnodes[1],
503*2d03409cSjeremylt                        mnodes[2]); CHKERRQ(ierr);
5040c59ef15SJeremy L Thompson   }
5050c59ef15SJeremy L Thompson 
5060c59ef15SJeremy L Thompson   {
5070c59ef15SJeremy L Thompson     lsize = 1;
50834d77899SValeria Barra     for (int d=0; d<dim; d++) {
509b52ddc88Sjeremylt       lnodes[d] = melem[d]*degree + 1;
510b52ddc88Sjeremylt       lsize *= lnodes[d];
5110c59ef15SJeremy L Thompson     }
5120c59ef15SJeremy L Thompson     ierr = VecCreate(PETSC_COMM_SELF, &Xloc); CHKERRQ(ierr);
51334d77899SValeria Barra     ierr = VecSetSizes(Xloc, lsize*ncompu, PETSC_DECIDE); CHKERRQ(ierr);
5140c59ef15SJeremy L Thompson     ierr = VecSetUp(Xloc); CHKERRQ(ierr);
5150c59ef15SJeremy L Thompson 
5160c59ef15SJeremy L Thompson     // Create local-to-global scatter
5170c59ef15SJeremy L Thompson     PetscInt *ltogind, *ltogind0, *locind, l0count;
5180c59ef15SJeremy L Thompson     IS ltogis, ltogis0, locis;
51934d77899SValeria Barra     PetscInt gstart[2][2][2], gmnodes[2][2][2][dim];
5200c59ef15SJeremy L Thompson 
5210c59ef15SJeremy L Thompson     for (int i=0; i<2; i++) {
5220c59ef15SJeremy L Thompson       for (int j=0; j<2; j++) {
5230c59ef15SJeremy L Thompson         for (int k=0; k<2; k++) {
5240c59ef15SJeremy L Thompson           PetscInt ijkrank[3] = {irank[0]+i, irank[1]+j, irank[2]+k};
5250c59ef15SJeremy L Thompson           gstart[i][j][k] = GlobalStart(p, ijkrank, degree, melem);
526b52ddc88Sjeremylt           GlobalNodes(p, ijkrank, degree, melem, gmnodes[i][j][k]);
5270c59ef15SJeremy L Thompson         }
5280c59ef15SJeremy L Thompson       }
5290c59ef15SJeremy L Thompson     }
5300c59ef15SJeremy L Thompson 
5310c59ef15SJeremy L Thompson     ierr = PetscMalloc1(lsize, &ltogind); CHKERRQ(ierr);
5320c59ef15SJeremy L Thompson     ierr = PetscMalloc1(lsize, &ltogind0); CHKERRQ(ierr);
5330c59ef15SJeremy L Thompson     ierr = PetscMalloc1(lsize, &locind); CHKERRQ(ierr);
5340c59ef15SJeremy L Thompson     l0count = 0;
535288c0443SJeremy L Thompson     for (PetscInt i=0,ir,ii; ir=i>=mnodes[0], ii=i-ir*mnodes[0], i<lnodes[0]; i++)
536288c0443SJeremy L Thompson       for (PetscInt j=0,jr,jj; jr=j>=mnodes[1], jj=j-jr*mnodes[1], j<lnodes[1]; j++)
537b52ddc88Sjeremylt         for (PetscInt k=0,kr,kk; kr=k>=mnodes[2], kk=k-kr*mnodes[2], k<lnodes[2]; k++) {
538b52ddc88Sjeremylt           PetscInt here = (i*lnodes[1]+j)*lnodes[2]+k;
5390c59ef15SJeremy L Thompson           ltogind[here] =
540b52ddc88Sjeremylt             gstart[ir][jr][kr] + (ii*gmnodes[ir][jr][kr][1]+jj)*gmnodes[ir][jr][kr][2]+kk;
5410c59ef15SJeremy L Thompson           if ((irank[0] == 0 && i == 0)
5420c59ef15SJeremy L Thompson               || (irank[1] == 0 && j == 0)
5430c59ef15SJeremy L Thompson               || (irank[2] == 0 && k == 0)
544b52ddc88Sjeremylt               || (irank[0]+1 == p[0] && i+1 == lnodes[0])
545b52ddc88Sjeremylt               || (irank[1]+1 == p[1] && j+1 == lnodes[1])
546b52ddc88Sjeremylt               || (irank[2]+1 == p[2] && k+1 == lnodes[2]))
5470c59ef15SJeremy L Thompson             continue;
5480c59ef15SJeremy L Thompson           ltogind0[l0count] = ltogind[here];
5490c59ef15SJeremy L Thompson           locind[l0count++] = here;
5500c59ef15SJeremy L Thompson         }
55134d77899SValeria Barra     ierr = ISCreateBlock(comm, ncompu, lsize, ltogind, PETSC_OWN_POINTER,
5520c59ef15SJeremy L Thompson                          &ltogis); CHKERRQ(ierr);
5530c59ef15SJeremy L Thompson     ierr = VecScatterCreate(Xloc, NULL, X, ltogis, &ltog); CHKERRQ(ierr);
5540c59ef15SJeremy L Thompson     CHKERRQ(ierr);
55534d77899SValeria Barra     ierr = ISCreateBlock(comm, ncompu, l0count, ltogind0, PETSC_OWN_POINTER,
5560c59ef15SJeremy L Thompson                          &ltogis0); CHKERRQ(ierr);
55734d77899SValeria Barra     ierr = ISCreateBlock(comm, ncompu, l0count, locind, PETSC_OWN_POINTER,
5580c59ef15SJeremy L Thompson                          &locis); CHKERRQ(ierr);
5590c59ef15SJeremy L Thompson     ierr = VecScatterCreate(Xloc, locis, X, ltogis0, &ltog0); CHKERRQ(ierr);
5600c59ef15SJeremy L Thompson     {
5610c59ef15SJeremy L Thompson       // Create global-to-global scatter for Dirichlet values (everything not in
5620c59ef15SJeremy L Thompson       // ltogis0, which is the range of ltog0)
5630c59ef15SJeremy L Thompson       PetscInt xstart, xend, *indD, countD = 0;
5640c59ef15SJeremy L Thompson       IS isD;
5650c59ef15SJeremy L Thompson       const PetscScalar *x;
5660c59ef15SJeremy L Thompson       ierr = VecZeroEntries(Xloc); CHKERRQ(ierr);
5670c59ef15SJeremy L Thompson       ierr = VecSet(X, 1.0); CHKERRQ(ierr);
5680c59ef15SJeremy L Thompson       ierr = VecScatterBegin(ltog0, Xloc, X, INSERT_VALUES, SCATTER_FORWARD);
5690c59ef15SJeremy L Thompson       CHKERRQ(ierr);
5700c59ef15SJeremy L Thompson       ierr = VecScatterEnd(ltog0, Xloc, X, INSERT_VALUES, SCATTER_FORWARD);
5710c59ef15SJeremy L Thompson       CHKERRQ(ierr);
5720c59ef15SJeremy L Thompson       ierr = VecGetOwnershipRange(X, &xstart, &xend); CHKERRQ(ierr);
5730c59ef15SJeremy L Thompson       ierr = PetscMalloc1(xend-xstart, &indD); CHKERRQ(ierr);
5740c59ef15SJeremy L Thompson       ierr = VecGetArrayRead(X, &x); CHKERRQ(ierr);
5750c59ef15SJeremy L Thompson       for (PetscInt i=0; i<xend-xstart; i++) {
5760c59ef15SJeremy L Thompson         if (x[i] == 1.) indD[countD++] = xstart + i;
5770c59ef15SJeremy L Thompson       }
5780c59ef15SJeremy L Thompson       ierr = VecRestoreArrayRead(X, &x); CHKERRQ(ierr);
5790c59ef15SJeremy L Thompson       ierr = ISCreateGeneral(comm, countD, indD, PETSC_COPY_VALUES, &isD);
5800c59ef15SJeremy L Thompson       CHKERRQ(ierr);
5810c59ef15SJeremy L Thompson       ierr = PetscFree(indD); CHKERRQ(ierr);
5820c59ef15SJeremy L Thompson       ierr = VecScatterCreate(X, isD, X, isD, &gtogD); CHKERRQ(ierr);
5830c59ef15SJeremy L Thompson       ierr = ISDestroy(&isD); CHKERRQ(ierr);
5840c59ef15SJeremy L Thompson     }
5850c59ef15SJeremy L Thompson     ierr = ISDestroy(&ltogis); CHKERRQ(ierr);
5860c59ef15SJeremy L Thompson     ierr = ISDestroy(&ltogis0); CHKERRQ(ierr);
5870c59ef15SJeremy L Thompson     ierr = ISDestroy(&locis); CHKERRQ(ierr);
5880c59ef15SJeremy L Thompson   }
5890c59ef15SJeremy L Thompson 
59034d77899SValeria Barra   // CEED bases
59134d77899SValeria Barra   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompu, P, Q,
592e07e9ddcSjeremylt                                   bpOptions[bpChoice].qmode, &basisu);
59334d77899SValeria Barra   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, Q,
594e07e9ddcSjeremylt                                   bpOptions[bpChoice].qmode, &basisx);
5950c59ef15SJeremy L Thompson 
59634d77899SValeria Barra   // CEED restrictions
59734d77899SValeria Barra   CreateRestriction(ceed, melem, P, ncompu, &Erestrictu);
59834d77899SValeria Barra   CreateRestriction(ceed, melem, 2, dim, &Erestrictx);
5990c59ef15SJeremy L Thompson   CeedInt nelem = melem[0]*melem[1]*melem[2];
60034d77899SValeria Barra   CeedElemRestrictionCreateIdentity(ceed, nelem, Q*Q*Q, nelem*Q*Q*Q, ncompu,
6010c59ef15SJeremy L Thompson                                     &Erestrictui);
6020c59ef15SJeremy L Thompson   CeedElemRestrictionCreateIdentity(ceed, nelem,
603ea339264SYohann                                     Q*Q*Q,
604ea339264SYohann                                     nelem*Q*Q*Q,
605ea339264SYohann                                     bpOptions[bpChoice].qdatasize, &Erestrictqdi);
6060c59ef15SJeremy L Thompson   CeedElemRestrictionCreateIdentity(ceed, nelem, Q*Q*Q, nelem*Q*Q*Q, 1,
6070c59ef15SJeremy L Thompson                                     &Erestrictxi);
6080c59ef15SJeremy L Thompson   {
6090c59ef15SJeremy L Thompson     CeedScalar *xloc;
6100c59ef15SJeremy L Thompson     CeedInt shape[3] = {melem[0]+1, melem[1]+1, melem[2]+1}, len =
6110c59ef15SJeremy L Thompson                          shape[0]*shape[1]*shape[2];
61234d77899SValeria Barra     xloc = malloc(len*ncompx*sizeof xloc[0]);
6130c59ef15SJeremy L Thompson     for (CeedInt i=0; i<shape[0]; i++) {
6140c59ef15SJeremy L Thompson       for (CeedInt j=0; j<shape[1]; j++) {
6150c59ef15SJeremy L Thompson         for (CeedInt k=0; k<shape[2]; k++) {
6160c59ef15SJeremy L Thompson           xloc[((i*shape[1]+j)*shape[2]+k) + 0*len] = 1.*(irank[0]*melem[0]+i) /
6170c59ef15SJeremy L Thompson               (p[0]*melem[0]);
6180c59ef15SJeremy L Thompson           xloc[((i*shape[1]+j)*shape[2]+k) + 1*len] = 1.*(irank[1]*melem[1]+j) /
6190c59ef15SJeremy L Thompson               (p[1]*melem[1]);
6200c59ef15SJeremy L Thompson           xloc[((i*shape[1]+j)*shape[2]+k) + 2*len] = 1.*(irank[2]*melem[2]+k) /
6210c59ef15SJeremy L Thompson               (p[2]*melem[2]);
6220c59ef15SJeremy L Thompson         }
6230c59ef15SJeremy L Thompson       }
6240c59ef15SJeremy L Thompson     }
62534d77899SValeria Barra     CeedVectorCreate(ceed, len*ncompx, &xcoord);
6260c59ef15SJeremy L Thompson     CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_OWN_POINTER, xloc);
6270c59ef15SJeremy L Thompson   }
6280c59ef15SJeremy L Thompson 
6296c5df90dSjeremylt   // Create the Qfunction that builds the operator quadrature data
6306c5df90dSjeremylt   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].setupgeo,
6316c5df90dSjeremylt                               bpOptions[bpChoice].setupgeofname, &qf_setupgeo);
6326c5df90dSjeremylt   CeedQFunctionAddInput(qf_setupgeo, "dx", ncompx*dim, CEED_EVAL_GRAD);
6336c5df90dSjeremylt   CeedQFunctionAddInput(qf_setupgeo, "weight", 1, CEED_EVAL_WEIGHT);
6346c5df90dSjeremylt   CeedQFunctionAddOutput(qf_setupgeo, "qdata", bpOptions[bpChoice].qdatasize,
6350c59ef15SJeremy L Thompson                          CEED_EVAL_NONE);
6366c5df90dSjeremylt 
6376c5df90dSjeremylt   // Create the Qfunction that sets up the RHS and true solution
6386c5df90dSjeremylt   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].setuprhs,
6396c5df90dSjeremylt                               bpOptions[bpChoice].setuprhsfname, &qf_setuprhs);
6406c5df90dSjeremylt   CeedQFunctionAddInput(qf_setuprhs, "x", ncompx, CEED_EVAL_INTERP);
6416c5df90dSjeremylt   CeedQFunctionAddInput(qf_setuprhs, "dx", ncompx*dim, CEED_EVAL_GRAD);
6426c5df90dSjeremylt   CeedQFunctionAddInput(qf_setuprhs, "weight", 1, CEED_EVAL_WEIGHT);
6436c5df90dSjeremylt   CeedQFunctionAddOutput(qf_setuprhs, "true_soln", ncompu, CEED_EVAL_NONE);
6446c5df90dSjeremylt   CeedQFunctionAddOutput(qf_setuprhs, "rhs", ncompu, CEED_EVAL_INTERP);
6450c59ef15SJeremy L Thompson 
6460c59ef15SJeremy L Thompson   // Set up PDE operator
6470c59ef15SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].apply,
6480c59ef15SJeremy L Thompson                               bpOptions[bpChoice].applyfname, &qf_apply);
6490c59ef15SJeremy L Thompson   // Add inputs and outputs
6506c5df90dSjeremylt   CeedInt inscale = bpOptions[bpChoice].inmode==CEED_EVAL_GRAD ? 3 : 1;
6516c5df90dSjeremylt   CeedInt outscale = bpOptions[bpChoice].outmode==CEED_EVAL_GRAD ? 3 : 1;
6526c5df90dSjeremylt   CeedQFunctionAddInput(qf_apply, "u", ncompu*inscale,
6534d537eeaSYohann                         bpOptions[bpChoice].inmode);
654a2fa7910SValeria Barra   CeedQFunctionAddInput(qf_apply, "qdata", bpOptions[bpChoice].qdatasize,
6550c59ef15SJeremy L Thompson                         CEED_EVAL_NONE);
6566c5df90dSjeremylt   CeedQFunctionAddOutput(qf_apply, "v", ncompu*outscale,
6574d537eeaSYohann                          bpOptions[bpChoice].outmode);
6580c59ef15SJeremy L Thompson 
6590c59ef15SJeremy L Thompson   // Create the error qfunction
6600c59ef15SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].error,
6610c59ef15SJeremy L Thompson                               bpOptions[bpChoice].errorfname, &qf_error);
66234d77899SValeria Barra   CeedQFunctionAddInput(qf_error, "u", ncompu, CEED_EVAL_INTERP);
66334d77899SValeria Barra   CeedQFunctionAddInput(qf_error, "true_soln", ncompu, CEED_EVAL_NONE);
66434d77899SValeria Barra   CeedQFunctionAddOutput(qf_error, "error", ncompu, CEED_EVAL_NONE);
6650c59ef15SJeremy L Thompson 
6660c59ef15SJeremy L Thompson   // Create the persistent vectors that will be needed in setup
667819eb1b3Sjeremylt   CeedInt nqpts;
668819eb1b3Sjeremylt   CeedBasisGetNumQuadraturePoints(basisu, &nqpts);
669a2fa7910SValeria Barra   CeedVectorCreate(ceed, bpOptions[bpChoice].qdatasize*nelem*nqpts, &qdata);
67034d77899SValeria Barra   CeedVectorCreate(ceed, nelem*nqpts*ncompu, &target);
67134d77899SValeria Barra   CeedVectorCreate(ceed, lsize*ncompu, &rhsceed);
6720c59ef15SJeremy L Thompson 
6734b5b4ec1Sjeremylt   // Create the operator that builds the quadrature data for the ceed operator
6746c5df90dSjeremylt   CeedOperatorCreate(ceed, qf_setupgeo, NULL, NULL, &op_setupgeo);
6756c5df90dSjeremylt   CeedOperatorSetField(op_setupgeo, "dx", Erestrictx, CEED_NOTRANSPOSE,
6760c59ef15SJeremy L Thompson                        basisx, CEED_VECTOR_ACTIVE);
6776c5df90dSjeremylt   CeedOperatorSetField(op_setupgeo, "weight", Erestrictxi, CEED_NOTRANSPOSE,
6780c59ef15SJeremy L Thompson                        basisx, CEED_VECTOR_NONE);
6796c5df90dSjeremylt   CeedOperatorSetField(op_setupgeo, "qdata", Erestrictqdi, CEED_NOTRANSPOSE,
6800c59ef15SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
6816c5df90dSjeremylt 
6826c5df90dSjeremylt   // Create the operator that builds the RHS and true solution
6836c5df90dSjeremylt   CeedOperatorCreate(ceed, qf_setuprhs, NULL, NULL, &op_setuprhs);
6846c5df90dSjeremylt   CeedOperatorSetField(op_setuprhs, "x", Erestrictx, CEED_NOTRANSPOSE,
6856c5df90dSjeremylt                        basisx, CEED_VECTOR_ACTIVE);
6866c5df90dSjeremylt   CeedOperatorSetField(op_setuprhs, "dx", Erestrictx, CEED_NOTRANSPOSE,
6876c5df90dSjeremylt                        basisx, CEED_VECTOR_ACTIVE);
6886c5df90dSjeremylt   CeedOperatorSetField(op_setuprhs, "weight", Erestrictxi, CEED_NOTRANSPOSE,
6896c5df90dSjeremylt                        basisx, CEED_VECTOR_NONE);
6906c5df90dSjeremylt   CeedOperatorSetField(op_setuprhs, "true_soln", Erestrictui, CEED_NOTRANSPOSE,
6910c59ef15SJeremy L Thompson                        CEED_BASIS_COLLOCATED, target);
6926c5df90dSjeremylt   CeedOperatorSetField(op_setuprhs, "rhs", Erestrictu, CEED_TRANSPOSE,
6936c5df90dSjeremylt                        basisu, CEED_VECTOR_ACTIVE);
6940c59ef15SJeremy L Thompson 
6954b5b4ec1Sjeremylt   // Create the mass or diff operator
6960c59ef15SJeremy L Thompson   CeedOperatorCreate(ceed, qf_apply, NULL, NULL, &op_apply);
6970c59ef15SJeremy L Thompson   CeedOperatorSetField(op_apply, "u", Erestrictu, CEED_TRANSPOSE,
6980c59ef15SJeremy L Thompson                        basisu, CEED_VECTOR_ACTIVE);
699a2fa7910SValeria Barra   CeedOperatorSetField(op_apply, "qdata", Erestrictqdi, CEED_NOTRANSPOSE,
700a2fa7910SValeria Barra                        CEED_BASIS_COLLOCATED, qdata);
7010c59ef15SJeremy L Thompson   CeedOperatorSetField(op_apply, "v", Erestrictu, CEED_TRANSPOSE,
7020c59ef15SJeremy L Thompson                        basisu, CEED_VECTOR_ACTIVE);
7030c59ef15SJeremy L Thompson 
7040c59ef15SJeremy L Thompson   // Create the error operator
7050c59ef15SJeremy L Thompson   CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error);
7060c59ef15SJeremy L Thompson   CeedOperatorSetField(op_error, "u", Erestrictu, CEED_TRANSPOSE,
7070c59ef15SJeremy L Thompson                        basisu, CEED_VECTOR_ACTIVE);
7080c59ef15SJeremy L Thompson   CeedOperatorSetField(op_error, "true_soln", Erestrictui, CEED_NOTRANSPOSE,
7090c59ef15SJeremy L Thompson                        CEED_BASIS_COLLOCATED, target);
7100c59ef15SJeremy L Thompson   CeedOperatorSetField(op_error, "error", Erestrictui, CEED_NOTRANSPOSE,
7110c59ef15SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
7120c59ef15SJeremy L Thompson 
7130c59ef15SJeremy L Thompson   // Set up Mat
7140c59ef15SJeremy L Thompson   ierr = PetscMalloc1(1, &user); CHKERRQ(ierr);
7150c59ef15SJeremy L Thompson   user->comm = comm;
7160c59ef15SJeremy L Thompson   user->ltog = ltog;
7170c59ef15SJeremy L Thompson   if (bpChoice != CEED_BP1 && bpChoice != CEED_BP2) {
7180c59ef15SJeremy L Thompson     user->ltog0 = ltog0;
7190c59ef15SJeremy L Thompson     user->gtogD = gtogD;
7200c59ef15SJeremy L Thompson   }
7210c59ef15SJeremy L Thompson   user->Xloc = Xloc;
7220c59ef15SJeremy L Thompson   ierr = VecDuplicate(Xloc, &user->Yloc); CHKERRQ(ierr);
72334d77899SValeria Barra   CeedVectorCreate(ceed, lsize*ncompu, &user->xceed);
72434d77899SValeria Barra   CeedVectorCreate(ceed, lsize*ncompu, &user->yceed);
7250c59ef15SJeremy L Thompson   user->op = op_apply;
726a2fa7910SValeria Barra   user->qdata = qdata;
7270c59ef15SJeremy L Thompson   user->ceed = ceed;
7280c59ef15SJeremy L Thompson 
72934d77899SValeria Barra   ierr = MatCreateShell(comm, mnodes[0]*mnodes[1]*mnodes[2]*ncompu,
73034d77899SValeria Barra                         mnodes[0]*mnodes[1]*mnodes[2]*ncompu,
7310c59ef15SJeremy L Thompson                         PETSC_DECIDE, PETSC_DECIDE, user, &mat); CHKERRQ(ierr);
7320c59ef15SJeremy L Thompson   if (bpChoice == CEED_BP1 || bpChoice == CEED_BP2) {
7330c59ef15SJeremy L Thompson     ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Mass);
7340c59ef15SJeremy L Thompson     CHKERRQ(ierr);
7350c59ef15SJeremy L Thompson   } else {
7360c59ef15SJeremy L Thompson     ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Diff);
7370c59ef15SJeremy L Thompson     CHKERRQ(ierr);
7380c59ef15SJeremy L Thompson   }
7390c59ef15SJeremy L Thompson   ierr = MatCreateVecs(mat, &rhs, NULL); CHKERRQ(ierr);
7400c59ef15SJeremy L Thompson 
7410c59ef15SJeremy L Thompson   // Get RHS vector
7420c59ef15SJeremy L Thompson   ierr = VecDuplicate(Xloc, &rhsloc); CHKERRQ(ierr);
7430c59ef15SJeremy L Thompson   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
7440c59ef15SJeremy L Thompson   ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr);
7450c59ef15SJeremy L Thompson   CeedVectorSetArray(rhsceed, CEED_MEM_HOST, CEED_USE_POINTER, r);
7460c59ef15SJeremy L Thompson 
747a2fa7910SValeria Barra   // Setup qdata, rhs, and target
7486c5df90dSjeremylt   CeedOperatorApply(op_setupgeo, xcoord, qdata, CEED_REQUEST_IMMEDIATE);
7496c5df90dSjeremylt   CeedOperatorApply(op_setuprhs, xcoord, rhsceed, CEED_REQUEST_IMMEDIATE);
7500c59ef15SJeremy L Thompson   ierr = CeedVectorSyncArray(rhsceed, CEED_MEM_HOST); CHKERRQ(ierr);
7510c59ef15SJeremy L Thompson   CeedVectorDestroy(&xcoord);
7520c59ef15SJeremy L Thompson 
7530c59ef15SJeremy L Thompson   // Gather RHS
7540c59ef15SJeremy L Thompson   ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr);
7550c59ef15SJeremy L Thompson   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
7560c59ef15SJeremy L Thompson   ierr = VecScatterBegin(ltog, rhsloc, rhs, ADD_VALUES, SCATTER_FORWARD);
7570c59ef15SJeremy L Thompson   CHKERRQ(ierr);
7580c59ef15SJeremy L Thompson   ierr = VecScatterEnd(ltog, rhsloc, rhs, ADD_VALUES, SCATTER_FORWARD);
7590c59ef15SJeremy L Thompson   CHKERRQ(ierr);
7600c59ef15SJeremy L Thompson   CeedVectorDestroy(&rhsceed);
7610c59ef15SJeremy L Thompson 
7620c59ef15SJeremy L Thompson   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
7630c59ef15SJeremy L Thompson   {
7640c59ef15SJeremy L Thompson     PC pc;
7650c59ef15SJeremy L Thompson     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
7660c59ef15SJeremy L Thompson     if (bpChoice == CEED_BP1 || bpChoice == CEED_BP2) {
7670c59ef15SJeremy L Thompson       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
7680c59ef15SJeremy L Thompson       ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr);
7690c59ef15SJeremy L Thompson     } else {
7700c59ef15SJeremy L Thompson       ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr);
7710c59ef15SJeremy L Thompson     }
7720c59ef15SJeremy L Thompson     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
7730c59ef15SJeremy L Thompson     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
7740c59ef15SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
7750c59ef15SJeremy L Thompson                             PETSC_DEFAULT); CHKERRQ(ierr);
7760c59ef15SJeremy L Thompson   }
7770c59ef15SJeremy L Thompson   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
7780c59ef15SJeremy L Thompson   ierr = KSPSetOperators(ksp, mat, mat); CHKERRQ(ierr);
7790c59ef15SJeremy L Thompson   // First run, if benchmarking
7800c59ef15SJeremy L Thompson   if (benchmark_mode) {
7810c59ef15SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
7820c59ef15SJeremy L Thompson     CHKERRQ(ierr);
7830c59ef15SJeremy L Thompson     my_rt_start = MPI_Wtime();
7840c59ef15SJeremy L Thompson     ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
7850c59ef15SJeremy L Thompson     my_rt = MPI_Wtime() - my_rt_start;
786b5009ec9Sjeremylt     ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
787b5009ec9Sjeremylt     CHKERRQ(ierr);
7880c59ef15SJeremy L Thompson     // Set maxits based on first iteration timing
7890c59ef15SJeremy L Thompson     if (my_rt > 0.02) {
7900c59ef15SJeremy L Thompson       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
7910c59ef15SJeremy L Thompson       CHKERRQ(ierr);
7920c59ef15SJeremy L Thompson     } else {
7930c59ef15SJeremy L Thompson       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
7940c59ef15SJeremy L Thompson       CHKERRQ(ierr);
7950c59ef15SJeremy L Thompson     }
7960c59ef15SJeremy L Thompson   }
7970c59ef15SJeremy L Thompson   // Timed solve
7983b3d4a15SJeremy L Thompson   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
7990c59ef15SJeremy L Thompson   my_rt_start = MPI_Wtime();
8000c59ef15SJeremy L Thompson   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
8010c59ef15SJeremy L Thompson   my_rt = MPI_Wtime() - my_rt_start;
8020c59ef15SJeremy L Thompson   {
8030c59ef15SJeremy L Thompson     KSPType ksptype;
8040c59ef15SJeremy L Thompson     KSPConvergedReason reason;
8050c59ef15SJeremy L Thompson     PetscReal rnorm;
8060c59ef15SJeremy L Thompson     PetscInt its;
8070c59ef15SJeremy L Thompson     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
8080c59ef15SJeremy L Thompson     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
8090c59ef15SJeremy L Thompson     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
8100c59ef15SJeremy L Thompson     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
8110c59ef15SJeremy L Thompson     if (!test_mode || reason < 0 || rnorm > 1e-8) {
8123ce86546SJeremy L Thompson       ierr = PetscPrintf(comm,
8133ce86546SJeremy L Thompson                          "  KSP:\n"
8143ce86546SJeremy L Thompson                          "    KSP Type                           : %s\n"
8153ce86546SJeremy L Thompson                          "    KSP Convergence                    : %s\n"
8163ce86546SJeremy L Thompson                          "    Total KSP Iterations               : %D\n"
8173ce86546SJeremy L Thompson                          "    Final rnorm                        : %e\n",
8183ce86546SJeremy L Thompson                          ksptype, KSPConvergedReasons[reason], its,
8193ce86546SJeremy L Thompson                          (double)rnorm); CHKERRQ(ierr);
8200c59ef15SJeremy L Thompson     }
8210c59ef15SJeremy L Thompson     if (benchmark_mode && (!test_mode)) {
8220c59ef15SJeremy L Thompson       CeedInt gsize;
8230c59ef15SJeremy L Thompson       ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
824b5009ec9Sjeremylt       ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
825b5009ec9Sjeremylt       CHKERRQ(ierr);
826b5009ec9Sjeremylt       ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
827b5009ec9Sjeremylt       CHKERRQ(ierr);
8280c59ef15SJeremy L Thompson       ierr = PetscPrintf(comm,
8293ce86546SJeremy L Thompson                          "  Performance:\n"
8303ce86546SJeremy L Thompson                          "    CG Solve Time                      : %g (%g) sec\n"
831b52ddc88Sjeremylt                          "    DoFs/Sec in CG                     : %g (%g) million\n",
8323ce86546SJeremy L Thompson                          rt_max, rt_min, 1e-6*gsize*its/rt_max,
8333ce86546SJeremy L Thompson                          1e-6*gsize*its/rt_min); CHKERRQ(ierr);
8340c59ef15SJeremy L Thompson     }
8350c59ef15SJeremy L Thompson   }
8360c59ef15SJeremy L Thompson 
8370c59ef15SJeremy L Thompson   {
8380c59ef15SJeremy L Thompson     PetscReal maxerror;
8390c59ef15SJeremy L Thompson     ierr = ComputeErrorMax(user, op_error, X, target, &maxerror); CHKERRQ(ierr);
8400c59ef15SJeremy L Thompson     PetscReal tol = (bpChoice == CEED_BP1 || bpChoice == CEED_BP2) ? 5e-3 : 5e-2;
8410c59ef15SJeremy L Thompson     if (!test_mode || maxerror > tol) {
8423ce86546SJeremy L Thompson       ierr = PetscPrintf(comm,
8433ce86546SJeremy L Thompson                          "    Pointwise Error (max)              : %e\n",
8443ce86546SJeremy L Thompson                          (double)maxerror); CHKERRQ(ierr);
8450c59ef15SJeremy L Thompson     }
8460c59ef15SJeremy L Thompson   }
8470c59ef15SJeremy L Thompson 
8486bd9afcaSjeremylt   if (write_solution) {
8496bd9afcaSjeremylt     PetscViewer vtkviewersoln;
8506bd9afcaSjeremylt 
8516bd9afcaSjeremylt     ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr);
8526bd9afcaSjeremylt     ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr);
8536bd9afcaSjeremylt     ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr);
8546bd9afcaSjeremylt     ierr = VecView(X, vtkviewersoln); CHKERRQ(ierr);
8556bd9afcaSjeremylt     ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr);
8566bd9afcaSjeremylt   }
8576bd9afcaSjeremylt 
8580c59ef15SJeremy L Thompson   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
8590c59ef15SJeremy L Thompson   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
8600c59ef15SJeremy L Thompson   ierr = VecDestroy(&X); CHKERRQ(ierr);
8610c59ef15SJeremy L Thompson   ierr = VecDestroy(&user->Xloc); CHKERRQ(ierr);
8620c59ef15SJeremy L Thompson   ierr = VecDestroy(&user->Yloc); CHKERRQ(ierr);
8630c59ef15SJeremy L Thompson   ierr = VecScatterDestroy(&ltog); CHKERRQ(ierr);
8640c59ef15SJeremy L Thompson   ierr = VecScatterDestroy(&ltog0); CHKERRQ(ierr);
8650c59ef15SJeremy L Thompson   ierr = VecScatterDestroy(&gtogD); CHKERRQ(ierr);
8660c59ef15SJeremy L Thompson   ierr = MatDestroy(&mat); CHKERRQ(ierr);
8670c59ef15SJeremy L Thompson   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
8680c59ef15SJeremy L Thompson 
8690c59ef15SJeremy L Thompson   CeedVectorDestroy(&user->xceed);
8700c59ef15SJeremy L Thompson   CeedVectorDestroy(&user->yceed);
871a2fa7910SValeria Barra   CeedVectorDestroy(&user->qdata);
8720c59ef15SJeremy L Thompson   CeedVectorDestroy(&target);
8736c5df90dSjeremylt   CeedOperatorDestroy(&op_setupgeo);
8746c5df90dSjeremylt   CeedOperatorDestroy(&op_setuprhs);
8750c59ef15SJeremy L Thompson   CeedOperatorDestroy(&op_apply);
8760c59ef15SJeremy L Thompson   CeedOperatorDestroy(&op_error);
8770c59ef15SJeremy L Thompson   CeedElemRestrictionDestroy(&Erestrictu);
8780c59ef15SJeremy L Thompson   CeedElemRestrictionDestroy(&Erestrictx);
8790c59ef15SJeremy L Thompson   CeedElemRestrictionDestroy(&Erestrictui);
8800c59ef15SJeremy L Thompson   CeedElemRestrictionDestroy(&Erestrictxi);
8810c59ef15SJeremy L Thompson   CeedElemRestrictionDestroy(&Erestrictqdi);
8826c5df90dSjeremylt   CeedQFunctionDestroy(&qf_setupgeo);
8836c5df90dSjeremylt   CeedQFunctionDestroy(&qf_setuprhs);
8840c59ef15SJeremy L Thompson   CeedQFunctionDestroy(&qf_apply);
8850c59ef15SJeremy L Thompson   CeedQFunctionDestroy(&qf_error);
8860c59ef15SJeremy L Thompson   CeedBasisDestroy(&basisu);
8870c59ef15SJeremy L Thompson   CeedBasisDestroy(&basisx);
8880c59ef15SJeremy L Thompson   CeedDestroy(&ceed);
8890c59ef15SJeremy L Thompson   ierr = PetscFree(user); CHKERRQ(ierr);
8900c59ef15SJeremy L Thompson   return PetscFinalize();
8910c59ef15SJeremy L Thompson }
892