xref: /libCEED/examples/petsc/bpsraw.c (revision 8d0bb2bb760f4e4c85ef3ddb5641d660bdf0f4d7)
1cb32e2e7SValeria Barra // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2cb32e2e7SValeria Barra // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3cb32e2e7SValeria Barra // reserved. See files LICENSE and NOTICE for details.
4cb32e2e7SValeria Barra //
5cb32e2e7SValeria Barra // This file is part of CEED, a collection of benchmarks, miniapps, software
6cb32e2e7SValeria Barra // libraries and APIs for efficient high-order finite element and spectral
7cb32e2e7SValeria Barra // element discretizations for exascale applications. For more information and
8cb32e2e7SValeria Barra // source code availability see http://github.com/ceed.
9cb32e2e7SValeria Barra //
10cb32e2e7SValeria Barra // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11cb32e2e7SValeria Barra // a collaborative effort of two U.S. Department of Energy organizations (Office
12cb32e2e7SValeria Barra // of Science and the National Nuclear Security Administration) responsible for
13cb32e2e7SValeria Barra // the planning and preparation of a capable exascale ecosystem, including
14cb32e2e7SValeria Barra // software, applications, hardware, advanced system engineering and early
15cb32e2e7SValeria Barra // testbed platforms, in support of the nation's exascale computing imperative.
16cb32e2e7SValeria Barra 
17cb32e2e7SValeria Barra //                        libCEED + PETSc Example: CEED BPs
18cb32e2e7SValeria Barra //
19cb32e2e7SValeria Barra // This example demonstrates a simple usage of libCEED with PETSc to solve the
20cb32e2e7SValeria Barra // CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps.
21cb32e2e7SValeria Barra //
22cb32e2e7SValeria Barra // The code is intentionally "raw", using only low-level communication
23cb32e2e7SValeria Barra // primitives.
24cb32e2e7SValeria Barra //
25cb32e2e7SValeria Barra // Build with:
26cb32e2e7SValeria Barra //
27cb32e2e7SValeria Barra //     make bpsraw [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
28cb32e2e7SValeria Barra //
29cb32e2e7SValeria Barra // Sample runs:
30cb32e2e7SValeria Barra //
3132d2ee49SValeria Barra //     ./bpsraw -problem bp1
3232d2ee49SValeria Barra //     ./bpsraw -problem bp2 -ceed /cpu/self
3332d2ee49SValeria Barra //     ./bpsraw -problem bp3 -ceed /gpu/occa
3432d2ee49SValeria Barra //     ./bpsraw -problem bp4 -ceed /cpu/occa
3532d2ee49SValeria Barra //     ./bpsraw -problem bp5 -ceed /omp/occa
3632d2ee49SValeria Barra //     ./bpsraw -problem bp6 -ceed /ocl/occa
37cb32e2e7SValeria Barra //
38cb32e2e7SValeria Barra //TESTARGS -ceed {ceed_resource} -test -problem bp2 -degree 5 -qextra 5
39cb32e2e7SValeria Barra 
40cb32e2e7SValeria Barra /// @file
41cb32e2e7SValeria Barra /// CEED BPs example using PETSc
42cb32e2e7SValeria Barra /// See bps.c for an implementation using DMPlex unstructured grids.
43cb32e2e7SValeria Barra const char help[] = "Solve CEED BPs using PETSc\n";
44cb32e2e7SValeria Barra 
45cb32e2e7SValeria Barra #include <stdbool.h>
46cb32e2e7SValeria Barra #include <string.h>
47cb32e2e7SValeria Barra #include <petscksp.h>
48cb32e2e7SValeria Barra #include <ceed.h>
49cb32e2e7SValeria Barra #include "qfunctions/bps/common.h"
50cb32e2e7SValeria Barra #include "qfunctions/bps/bp1.h"
51cb32e2e7SValeria Barra #include "qfunctions/bps/bp2.h"
52cb32e2e7SValeria Barra #include "qfunctions/bps/bp3.h"
53cb32e2e7SValeria Barra #include "qfunctions/bps/bp4.h"
54cb32e2e7SValeria Barra 
559396343dSjeremylt #include <petscsys.h>
569396343dSjeremylt #if PETSC_VERSION_LT(3,12,0)
579396343dSjeremylt #ifdef PETSC_HAVE_CUDA
589396343dSjeremylt #include <petsccuda.h>
599396343dSjeremylt // Note: With PETSc prior to version 3.12.0, providing the source path to
609396343dSjeremylt //       include 'cublas_v2.h' will be needed to use 'petsccuda.h'.
619396343dSjeremylt #endif
629396343dSjeremylt #endif
639396343dSjeremylt 
64cb32e2e7SValeria Barra static void Split3(PetscInt size, PetscInt m[3], bool reverse) {
65cb32e2e7SValeria Barra   for (PetscInt d=0,sizeleft=size; d<3; d++) {
66cb32e2e7SValeria Barra     PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(sizeleft, 1./(3 - d)));
67cb32e2e7SValeria Barra     while (try * (sizeleft / try) != sizeleft) try++;
68cb32e2e7SValeria Barra     m[reverse ? 2-d : d] = try;
69cb32e2e7SValeria Barra     sizeleft /= try;
70cb32e2e7SValeria Barra   }
71cb32e2e7SValeria Barra }
72cb32e2e7SValeria Barra 
73cb32e2e7SValeria Barra static PetscInt Max3(const PetscInt a[3]) {
74cb32e2e7SValeria Barra   return PetscMax(a[0], PetscMax(a[1], a[2]));
75cb32e2e7SValeria Barra }
76cb32e2e7SValeria Barra static PetscInt Min3(const PetscInt a[3]) {
77cb32e2e7SValeria Barra   return PetscMin(a[0], PetscMin(a[1], a[2]));
78cb32e2e7SValeria Barra }
79cb32e2e7SValeria Barra static void GlobalNodes(const PetscInt p[3], const PetscInt irank[3],
80cb32e2e7SValeria Barra                         PetscInt degree, const PetscInt melem[3],
81cb32e2e7SValeria Barra                         PetscInt mnodes[3]) {
82cb32e2e7SValeria Barra   for (int d=0; d<3; d++)
83cb32e2e7SValeria Barra     mnodes[d] = degree*melem[d] + (irank[d] == p[d]-1);
84cb32e2e7SValeria Barra }
85cb32e2e7SValeria Barra static PetscInt GlobalStart(const PetscInt p[3], const PetscInt irank[3],
86cb32e2e7SValeria Barra                             PetscInt degree, const PetscInt melem[3]) {
87cb32e2e7SValeria Barra   PetscInt start = 0;
88cb32e2e7SValeria Barra   // Dumb brute-force is easier to read
89cb32e2e7SValeria Barra   for (PetscInt i=0; i<p[0]; i++) {
90cb32e2e7SValeria Barra     for (PetscInt j=0; j<p[1]; j++) {
91cb32e2e7SValeria Barra       for (PetscInt k=0; k<p[2]; k++) {
92cb32e2e7SValeria Barra         PetscInt mnodes[3], ijkrank[] = {i,j,k};
93cb32e2e7SValeria Barra         if (i == irank[0] && j == irank[1] && k == irank[2]) return start;
94cb32e2e7SValeria Barra         GlobalNodes(p, ijkrank, degree, melem, mnodes);
95cb32e2e7SValeria Barra         start += mnodes[0] * mnodes[1] * mnodes[2];
96cb32e2e7SValeria Barra       }
97cb32e2e7SValeria Barra     }
98cb32e2e7SValeria Barra   }
99cb32e2e7SValeria Barra   return -1;
100cb32e2e7SValeria Barra }
101d979a051Sjeremylt static int CreateRestriction(Ceed ceed, const CeedInt melem[3], CeedInt P,
102d979a051Sjeremylt                              CeedInt ncomp, CeedElemRestriction *Erestrict) {
103cb32e2e7SValeria Barra   const PetscInt nelem = melem[0]*melem[1]*melem[2];
104cb32e2e7SValeria Barra   PetscInt mnodes[3], *idx, *idxp;
105cb32e2e7SValeria Barra 
106cb32e2e7SValeria Barra   // Get indicies
107cb32e2e7SValeria Barra   for (int d=0; d<3; d++) mnodes[d] = melem[d]*(P-1) + 1;
108cb32e2e7SValeria Barra   idxp = idx = malloc(nelem*P*P*P*sizeof idx[0]);
109cb32e2e7SValeria Barra   for (CeedInt i=0; i<melem[0]; i++)
110cb32e2e7SValeria Barra     for (CeedInt j=0; j<melem[1]; j++)
111cb32e2e7SValeria Barra       for (CeedInt k=0; k<melem[2]; k++,idxp += P*P*P)
112cb32e2e7SValeria Barra         for (CeedInt ii=0; ii<P; ii++)
113cb32e2e7SValeria Barra           for (CeedInt jj=0; jj<P; jj++)
114cb32e2e7SValeria Barra             for (CeedInt kk=0; kk<P; kk++) {
115cb32e2e7SValeria Barra               if (0) { // This is the C-style (i,j,k) ordering that I prefer
116d979a051Sjeremylt                 idxp[(ii*P+jj)*P+kk] = ncomp*(((i*(P-1)+ii)*mnodes[1]
117cb32e2e7SValeria Barra                                                + (j*(P-1)+jj))*mnodes[2]
118cb32e2e7SValeria Barra                                               + (k*(P-1)+kk));
119cb32e2e7SValeria Barra               } else { // (k,j,i) ordering for consistency with MFEM example
120d979a051Sjeremylt                 idxp[ii+P*(jj+P*kk)] = ncomp*(((i*(P-1)+ii)*mnodes[1]
121cb32e2e7SValeria Barra                                                + (j*(P-1)+jj))*mnodes[2]
122cb32e2e7SValeria Barra                                               + (k*(P-1)+kk));
123cb32e2e7SValeria Barra               }
124cb32e2e7SValeria Barra             }
125cb32e2e7SValeria Barra 
126cb32e2e7SValeria Barra   // Setup CEED restriction
127d979a051Sjeremylt   CeedElemRestrictionCreate(ceed, nelem, P*P*P, ncomp, 1,
128d979a051Sjeremylt                             mnodes[0]*mnodes[1]*mnodes[2]*ncomp,
129a8d32208Sjeremylt                             CEED_MEM_HOST, CEED_OWN_POINTER, idx, Erestrict);
130cb32e2e7SValeria Barra 
131cb32e2e7SValeria Barra   PetscFunctionReturn(0);
132cb32e2e7SValeria Barra }
133cb32e2e7SValeria Barra 
134cb32e2e7SValeria Barra // Data for PETSc
135cb32e2e7SValeria Barra typedef struct User_ *User;
136cb32e2e7SValeria Barra struct User_ {
137cb32e2e7SValeria Barra   MPI_Comm comm;
138cb32e2e7SValeria Barra   VecScatter ltog;              // Scatter for all entries
139cb32e2e7SValeria Barra   VecScatter ltog0;             // Skip Dirichlet values
140cb32e2e7SValeria Barra   VecScatter gtogD;             // global-to-global; only Dirichlet values
141cb32e2e7SValeria Barra   Vec Xloc, Yloc;
142cb32e2e7SValeria Barra   CeedVector xceed, yceed;
143cb32e2e7SValeria Barra   CeedOperator op;
144cb32e2e7SValeria Barra   CeedVector qdata;
145cb32e2e7SValeria Barra   Ceed ceed;
1469396343dSjeremylt   CeedMemType memtype;
1479396343dSjeremylt   int (*VecGetArray)(Vec, PetscScalar **);
1489396343dSjeremylt   int (*VecGetArrayRead)(Vec, const PetscScalar **);
1499396343dSjeremylt   int (*VecRestoreArray)(Vec, PetscScalar **);
1509396343dSjeremylt   int (*VecRestoreArrayRead)(Vec, const PetscScalar **);
1519396343dSjeremylt };
1529396343dSjeremylt 
1539396343dSjeremylt // MemType Options
1549396343dSjeremylt static const char *const memTypes[] = {"host","device","memType",
1559396343dSjeremylt                                        "CEED_MEM_",0
156cb32e2e7SValeria Barra                                       };
157cb32e2e7SValeria Barra 
158cb32e2e7SValeria Barra // BP Options
159cb32e2e7SValeria Barra typedef enum {
160cb32e2e7SValeria Barra   CEED_BP1 = 0, CEED_BP2 = 1, CEED_BP3 = 2,
161cb32e2e7SValeria Barra   CEED_BP4 = 3, CEED_BP5 = 4, CEED_BP6 = 5
162cb32e2e7SValeria Barra } bpType;
163cb32e2e7SValeria Barra static const char *const bpTypes[] = {"bp1","bp2","bp3","bp4","bp5","bp6",
164cb32e2e7SValeria Barra                                       "bpType","CEED_BP",0
165cb32e2e7SValeria Barra                                      };
166cb32e2e7SValeria Barra 
167cb32e2e7SValeria Barra // BP specific data
168cb32e2e7SValeria Barra typedef struct {
169cb32e2e7SValeria Barra   CeedInt ncompu, qdatasize, qextra;
170cb32e2e7SValeria Barra   CeedQFunctionUser setupgeo, setuprhs, apply, error;
171cb32e2e7SValeria Barra   const char *setupgeofname, *setuprhsfname, *applyfname, *errorfname;
172cb32e2e7SValeria Barra   CeedEvalMode inmode, outmode;
173cb32e2e7SValeria Barra   CeedQuadMode qmode;
174cb32e2e7SValeria Barra } bpData;
175cb32e2e7SValeria Barra 
176cb32e2e7SValeria Barra bpData bpOptions[6] = {
177cb32e2e7SValeria Barra   [CEED_BP1] = {
178cb32e2e7SValeria Barra     .ncompu = 1,
179cb32e2e7SValeria Barra     .qdatasize = 1,
180cb32e2e7SValeria Barra     .qextra = 1,
181cb32e2e7SValeria Barra     .setupgeo = SetupMassGeo,
182cb32e2e7SValeria Barra     .setuprhs = SetupMassRhs,
183cb32e2e7SValeria Barra     .apply = Mass,
184cb32e2e7SValeria Barra     .error = Error,
185cb32e2e7SValeria Barra     .setupgeofname = SetupMassGeo_loc,
186cb32e2e7SValeria Barra     .setuprhsfname = SetupMassRhs_loc,
187cb32e2e7SValeria Barra     .applyfname = Mass_loc,
188cb32e2e7SValeria Barra     .errorfname = Error_loc,
189cb32e2e7SValeria Barra     .inmode = CEED_EVAL_INTERP,
190cb32e2e7SValeria Barra     .outmode = CEED_EVAL_INTERP,
191cb32e2e7SValeria Barra     .qmode = CEED_GAUSS
192cb32e2e7SValeria Barra   },
193cb32e2e7SValeria Barra   [CEED_BP2] = {
194cb32e2e7SValeria Barra     .ncompu = 3,
195cb32e2e7SValeria Barra     .qdatasize = 1,
196cb32e2e7SValeria Barra     .qextra = 1,
197cb32e2e7SValeria Barra     .setupgeo = SetupMassGeo,
198cb32e2e7SValeria Barra     .setuprhs = SetupMassRhs3,
199cb32e2e7SValeria Barra     .apply = Mass3,
200cb32e2e7SValeria Barra     .error = Error3,
201cb32e2e7SValeria Barra     .setupgeofname = SetupMassGeo_loc,
202cb32e2e7SValeria Barra     .setuprhsfname = SetupMassRhs3_loc,
203cb32e2e7SValeria Barra     .applyfname = Mass3_loc,
204cb32e2e7SValeria Barra     .errorfname = Error3_loc,
205cb32e2e7SValeria Barra     .inmode = CEED_EVAL_INTERP,
206cb32e2e7SValeria Barra     .outmode = CEED_EVAL_INTERP,
207cb32e2e7SValeria Barra     .qmode = CEED_GAUSS
208cb32e2e7SValeria Barra   },
209cb32e2e7SValeria Barra   [CEED_BP3] = {
210cb32e2e7SValeria Barra     .ncompu = 1,
211cb32e2e7SValeria Barra     .qdatasize = 6,
212cb32e2e7SValeria Barra     .qextra = 1,
213cb32e2e7SValeria Barra     .setupgeo = SetupDiffGeo,
214cb32e2e7SValeria Barra     .setuprhs = SetupDiffRhs,
215cb32e2e7SValeria Barra     .apply = Diff,
216cb32e2e7SValeria Barra     .error = Error,
217cb32e2e7SValeria Barra     .setupgeofname = SetupDiffGeo_loc,
218cb32e2e7SValeria Barra     .setuprhsfname = SetupDiffRhs_loc,
219cb32e2e7SValeria Barra     .applyfname = Diff_loc,
220cb32e2e7SValeria Barra     .errorfname = Error_loc,
221cb32e2e7SValeria Barra     .inmode = CEED_EVAL_GRAD,
222cb32e2e7SValeria Barra     .outmode = CEED_EVAL_GRAD,
223cb32e2e7SValeria Barra     .qmode = CEED_GAUSS
224cb32e2e7SValeria Barra   },
225cb32e2e7SValeria Barra   [CEED_BP4] = {
226cb32e2e7SValeria Barra     .ncompu = 3,
227cb32e2e7SValeria Barra     .qdatasize = 6,
228cb32e2e7SValeria Barra     .qextra = 1,
229cb32e2e7SValeria Barra     .setupgeo = SetupDiffGeo,
230cb32e2e7SValeria Barra     .setuprhs = SetupDiffRhs3,
231cb32e2e7SValeria Barra     .apply = Diff3,
232cb32e2e7SValeria Barra     .error = Error3,
233cb32e2e7SValeria Barra     .setupgeofname = SetupDiffGeo_loc,
234cb32e2e7SValeria Barra     .setuprhsfname = SetupDiffRhs3_loc,
235cb32e2e7SValeria Barra     .applyfname = Diff_loc,
236cb32e2e7SValeria Barra     .errorfname = Error3_loc,
237cb32e2e7SValeria Barra     .inmode = CEED_EVAL_GRAD,
238cb32e2e7SValeria Barra     .outmode = CEED_EVAL_GRAD,
239cb32e2e7SValeria Barra     .qmode = CEED_GAUSS
240cb32e2e7SValeria Barra   },
241cb32e2e7SValeria Barra   [CEED_BP5] = {
242cb32e2e7SValeria Barra     .ncompu = 1,
243cb32e2e7SValeria Barra     .qdatasize = 6,
244cb32e2e7SValeria Barra     .qextra = 0,
245cb32e2e7SValeria Barra     .setupgeo = SetupDiffGeo,
246cb32e2e7SValeria Barra     .setuprhs = SetupDiffRhs,
247cb32e2e7SValeria Barra     .apply = Diff,
248cb32e2e7SValeria Barra     .error = Error,
249cb32e2e7SValeria Barra     .setupgeofname = SetupDiffGeo_loc,
250cb32e2e7SValeria Barra     .setuprhsfname = SetupDiffRhs_loc,
251cb32e2e7SValeria Barra     .applyfname = Diff_loc,
252cb32e2e7SValeria Barra     .errorfname = Error_loc,
253cb32e2e7SValeria Barra     .inmode = CEED_EVAL_GRAD,
254cb32e2e7SValeria Barra     .outmode = CEED_EVAL_GRAD,
255cb32e2e7SValeria Barra     .qmode = CEED_GAUSS_LOBATTO
256cb32e2e7SValeria Barra   },
257cb32e2e7SValeria Barra   [CEED_BP6] = {
258cb32e2e7SValeria Barra     .ncompu = 3,
259cb32e2e7SValeria Barra     .qdatasize = 6,
260cb32e2e7SValeria Barra     .qextra = 0,
261cb32e2e7SValeria Barra     .setupgeo = SetupDiffGeo,
262cb32e2e7SValeria Barra     .setuprhs = SetupDiffRhs3,
263cb32e2e7SValeria Barra     .apply = Diff3,
264cb32e2e7SValeria Barra     .error = Error3,
265cb32e2e7SValeria Barra     .setupgeofname = SetupDiffGeo_loc,
266cb32e2e7SValeria Barra     .setuprhsfname = SetupDiffRhs3_loc,
267cb32e2e7SValeria Barra     .applyfname = Diff_loc,
268cb32e2e7SValeria Barra     .errorfname = Error3_loc,
269cb32e2e7SValeria Barra     .inmode = CEED_EVAL_GRAD,
270cb32e2e7SValeria Barra     .outmode = CEED_EVAL_GRAD,
271cb32e2e7SValeria Barra     .qmode = CEED_GAUSS_LOBATTO
272cb32e2e7SValeria Barra   }
273cb32e2e7SValeria Barra };
274cb32e2e7SValeria Barra 
275cb32e2e7SValeria Barra // This function uses libCEED to compute the action of the mass matrix
276cb32e2e7SValeria Barra static PetscErrorCode MatMult_Mass(Mat A, Vec X, Vec Y) {
277cb32e2e7SValeria Barra   PetscErrorCode ierr;
278cb32e2e7SValeria Barra   User user;
279cb32e2e7SValeria Barra   PetscScalar *x, *y;
280cb32e2e7SValeria Barra 
281cb32e2e7SValeria Barra   PetscFunctionBeginUser;
2829396343dSjeremylt 
283cb32e2e7SValeria Barra   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
2849396343dSjeremylt 
2859396343dSjeremylt   // Global-to-local
286cb32e2e7SValeria Barra   ierr = VecScatterBegin(user->ltog, X, user->Xloc, INSERT_VALUES,
287cb32e2e7SValeria Barra                          SCATTER_REVERSE); CHKERRQ(ierr);
2889396343dSjeremylt   ierr = VecScatterEnd(user->ltog, X, user->Xloc, INSERT_VALUES,
2899396343dSjeremylt                        SCATTER_REVERSE); CHKERRQ(ierr);
290cb32e2e7SValeria Barra 
2919396343dSjeremylt   // Setup libCEED vectors
2929396343dSjeremylt   ierr = user->VecGetArrayRead(user->Xloc, (const PetscScalar **)&x);
2939396343dSjeremylt   CHKERRQ(ierr);
2949396343dSjeremylt   ierr = user->VecGetArray(user->Yloc, &y); CHKERRQ(ierr);
2959396343dSjeremylt   CeedVectorSetArray(user->xceed, user->memtype, CEED_USE_POINTER, x);
2969396343dSjeremylt   CeedVectorSetArray(user->yceed, user->memtype, CEED_USE_POINTER, y);
297cb32e2e7SValeria Barra 
2989396343dSjeremylt   // Apply libCEED operator
299cb32e2e7SValeria Barra   CeedOperatorApply(user->op, user->xceed, user->yceed,
300cb32e2e7SValeria Barra                     CEED_REQUEST_IMMEDIATE);
3019396343dSjeremylt   CeedVectorSyncArray(user->yceed, user->memtype);
302cb32e2e7SValeria Barra 
3039396343dSjeremylt   // Restore PETSc vectors
3049396343dSjeremylt   ierr = user->VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x);
3059396343dSjeremylt   CHKERRQ(ierr);
3069396343dSjeremylt   ierr = user->VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr);
307cb32e2e7SValeria Barra 
3089396343dSjeremylt   // Local-to-global
309cb32e2e7SValeria Barra   if (Y) {
310cb32e2e7SValeria Barra     ierr = VecZeroEntries(Y); CHKERRQ(ierr);
3119396343dSjeremylt     ierr = VecScatterBegin(user->ltog, user->Yloc, Y, ADD_VALUES,
3129396343dSjeremylt                            SCATTER_FORWARD); CHKERRQ(ierr);
3139396343dSjeremylt     ierr = VecScatterEnd(user->ltog, user->Yloc, Y, ADD_VALUES,
3149396343dSjeremylt                          SCATTER_FORWARD); CHKERRQ(ierr);
315cb32e2e7SValeria Barra   }
316cb32e2e7SValeria Barra   PetscFunctionReturn(0);
317cb32e2e7SValeria Barra }
318cb32e2e7SValeria Barra 
319cb32e2e7SValeria Barra // This function uses libCEED to compute the action of the Laplacian with
320cb32e2e7SValeria Barra // Dirichlet boundary conditions
321cb32e2e7SValeria Barra static PetscErrorCode MatMult_Diff(Mat A, Vec X, Vec Y) {
322cb32e2e7SValeria Barra   PetscErrorCode ierr;
323cb32e2e7SValeria Barra   User user;
324cb32e2e7SValeria Barra   PetscScalar *x, *y;
325cb32e2e7SValeria Barra 
326cb32e2e7SValeria Barra   PetscFunctionBeginUser;
3279396343dSjeremylt 
328cb32e2e7SValeria Barra   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
329cb32e2e7SValeria Barra 
330cb32e2e7SValeria Barra   // Global-to-local
331cb32e2e7SValeria Barra   ierr = VecScatterBegin(user->ltog0, X, user->Xloc, INSERT_VALUES,
332cb32e2e7SValeria Barra                          SCATTER_REVERSE); CHKERRQ(ierr);
333cb32e2e7SValeria Barra   ierr = VecScatterEnd(user->ltog0, X, user->Xloc, INSERT_VALUES,
334cb32e2e7SValeria Barra                        SCATTER_REVERSE);
335cb32e2e7SValeria Barra   CHKERRQ(ierr);
336cb32e2e7SValeria Barra 
3379396343dSjeremylt   // Setup libCEED vectors
3389396343dSjeremylt   ierr = user->VecGetArrayRead(user->Xloc, (const PetscScalar **)&x);
3399396343dSjeremylt   CHKERRQ(ierr);
3409396343dSjeremylt   ierr = user->VecGetArray(user->Yloc, &y); CHKERRQ(ierr);
3419396343dSjeremylt   CeedVectorSetArray(user->xceed, user->memtype, CEED_USE_POINTER, x);
3429396343dSjeremylt   CeedVectorSetArray(user->yceed, user->memtype, CEED_USE_POINTER, y);
343cb32e2e7SValeria Barra 
3449396343dSjeremylt   // Apply libCEED operator
345cb32e2e7SValeria Barra   CeedOperatorApply(user->op, user->xceed, user->yceed,
346cb32e2e7SValeria Barra                     CEED_REQUEST_IMMEDIATE);
3479396343dSjeremylt   CeedVectorSyncArray(user->yceed, user->memtype);
348cb32e2e7SValeria Barra 
349cb32e2e7SValeria Barra   // Restore PETSc vectors
3509396343dSjeremylt   ierr = user->VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x);
3519396343dSjeremylt   CHKERRQ(ierr);
3529396343dSjeremylt   ierr = user->VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr);
353cb32e2e7SValeria Barra 
354cb32e2e7SValeria Barra   // Local-to-global
355cb32e2e7SValeria Barra   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
356cb32e2e7SValeria Barra   ierr = VecScatterBegin(user->gtogD, X, Y, INSERT_VALUES, SCATTER_FORWARD);
357cb32e2e7SValeria Barra   CHKERRQ(ierr);
358cb32e2e7SValeria Barra   ierr = VecScatterEnd(user->gtogD, X, Y, INSERT_VALUES, SCATTER_FORWARD);
359cb32e2e7SValeria Barra   CHKERRQ(ierr);
3609396343dSjeremylt   ierr = VecScatterBegin(user->ltog0, user->Yloc, Y, ADD_VALUES,
3619396343dSjeremylt                          SCATTER_FORWARD); CHKERRQ(ierr);
362cb32e2e7SValeria Barra   ierr = VecScatterEnd(user->ltog0, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD);
363cb32e2e7SValeria Barra   CHKERRQ(ierr);
364cb32e2e7SValeria Barra 
365cb32e2e7SValeria Barra   PetscFunctionReturn(0);
366cb32e2e7SValeria Barra }
367cb32e2e7SValeria Barra 
368cb32e2e7SValeria Barra // This function calculates the error in the final solution
369c70bd2a0Sjeremylt static PetscErrorCode ComputeErrorMax(User user, CeedOperator operror, Vec X,
370cb32e2e7SValeria Barra                                       CeedVector target, PetscReal *maxerror) {
371cb32e2e7SValeria Barra   PetscErrorCode ierr;
372cb32e2e7SValeria Barra   PetscScalar *x;
373cb32e2e7SValeria Barra   CeedVector collocated_error;
374cb32e2e7SValeria Barra   CeedInt length;
375cb32e2e7SValeria Barra 
376cb32e2e7SValeria Barra   PetscFunctionBeginUser;
3779396343dSjeremylt 
378cb32e2e7SValeria Barra   CeedVectorGetLength(target, &length);
379cb32e2e7SValeria Barra   CeedVectorCreate(user->ceed, length, &collocated_error);
380cb32e2e7SValeria Barra 
381cb32e2e7SValeria Barra   // Global-to-local
382cb32e2e7SValeria Barra   ierr = VecScatterBegin(user->ltog, X, user->Xloc, INSERT_VALUES,
383cb32e2e7SValeria Barra                          SCATTER_REVERSE); CHKERRQ(ierr);
3849396343dSjeremylt   ierr = VecScatterEnd(user->ltog, X, user->Xloc, INSERT_VALUES,
3859396343dSjeremylt                        SCATTER_REVERSE); CHKERRQ(ierr);
3869396343dSjeremylt 
3879396343dSjeremylt   // Setup libCEED vector
3889396343dSjeremylt   ierr = user->VecGetArrayRead(user->Xloc, (const PetscScalar **)&x);
389cb32e2e7SValeria Barra   CHKERRQ(ierr);
3909396343dSjeremylt   CeedVectorSetArray(user->xceed, user->memtype, CEED_USE_POINTER, x);
391cb32e2e7SValeria Barra 
3929396343dSjeremylt   // Apply libCEED operator
393c70bd2a0Sjeremylt   CeedOperatorApply(operror, user->xceed, collocated_error,
394cb32e2e7SValeria Barra                     CEED_REQUEST_IMMEDIATE);
395cb32e2e7SValeria Barra 
396cb32e2e7SValeria Barra   // Restore PETSc vector
3979396343dSjeremylt   ierr = user->VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x);
3989396343dSjeremylt   CHKERRQ(ierr);
399cb32e2e7SValeria Barra 
400cb32e2e7SValeria Barra   // Reduce max error
401cb32e2e7SValeria Barra   *maxerror = 0;
402cb32e2e7SValeria Barra   const CeedScalar *e;
403cb32e2e7SValeria Barra   CeedVectorGetArrayRead(collocated_error, CEED_MEM_HOST, &e);
404cb32e2e7SValeria Barra   for (CeedInt i=0; i<length; i++) {
405cb32e2e7SValeria Barra     *maxerror = PetscMax(*maxerror, PetscAbsScalar(e[i]));
406cb32e2e7SValeria Barra   }
407cb32e2e7SValeria Barra   CeedVectorRestoreArrayRead(collocated_error, &e);
4089396343dSjeremylt   ierr = MPI_Allreduce(MPI_IN_PLACE, maxerror, 1, MPIU_REAL, MPIU_MAX,
4099396343dSjeremylt                        user->comm); CHKERRQ(ierr);
410cb32e2e7SValeria Barra 
411cb32e2e7SValeria Barra   // Cleanup
412cb32e2e7SValeria Barra   CeedVectorDestroy(&collocated_error);
413cb32e2e7SValeria Barra 
414cb32e2e7SValeria Barra   PetscFunctionReturn(0);
415cb32e2e7SValeria Barra }
416cb32e2e7SValeria Barra 
417cb32e2e7SValeria Barra int main(int argc, char **argv) {
418cb32e2e7SValeria Barra   PetscInt ierr;
419cb32e2e7SValeria Barra   MPI_Comm comm;
420cb32e2e7SValeria Barra   char ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self";
421cb32e2e7SValeria Barra   double my_rt_start, my_rt, rt_min, rt_max;
422cb32e2e7SValeria Barra   PetscInt degree, qextra, localnodes, localelem, melem[3], mnodes[3], p[3],
423cb32e2e7SValeria Barra            irank[3], lnodes[3], lsize, ncompu = 1;
424cb32e2e7SValeria Barra   PetscScalar *r;
425cb32e2e7SValeria Barra   PetscBool test_mode, benchmark_mode, write_solution;
426cb32e2e7SValeria Barra   PetscMPIInt size, rank;
42709a940d7Sjeremylt   PetscLogStage solvestage;
428cb32e2e7SValeria Barra   Vec X, Xloc, rhs, rhsloc;
429cb32e2e7SValeria Barra   Mat mat;
430cb32e2e7SValeria Barra   KSP ksp;
431cb32e2e7SValeria Barra   VecScatter ltog, ltog0, gtogD;
432cb32e2e7SValeria Barra   User user;
433cb32e2e7SValeria Barra   Ceed ceed;
434cb32e2e7SValeria Barra   CeedBasis basisx, basisu;
43515910d16Sjeremylt   CeedElemRestriction Erestrictx, Erestrictu, Erestrictui, Erestrictqdi;
436c70bd2a0Sjeremylt   CeedQFunction qfsetupgeo, qfsetuprhs, qfapply, qferror;
437c70bd2a0Sjeremylt   CeedOperator opsetupgeo, opsetuprhs, opapply, operror;
438cb32e2e7SValeria Barra   CeedVector xcoord, qdata, rhsceed, target;
4399396343dSjeremylt   CeedMemType memtyperequested;
440cb32e2e7SValeria Barra   CeedInt P, Q;
441cb32e2e7SValeria Barra   const CeedInt dim = 3, ncompx = 3;
442199551f5Sjeremylt   bpType bpchoice;
443cb32e2e7SValeria Barra 
4449396343dSjeremylt   // Check for PETSc CUDA availability
4459396343dSjeremylt   PetscBool petschavecuda, setmemtyperequest = PETSC_FALSE;
4469396343dSjeremylt   // *INDENT-OFF*
4479396343dSjeremylt   #ifdef PETSC_HAVE_CUDA
4489396343dSjeremylt   petschavecuda = PETSC_TRUE;
4499396343dSjeremylt   #else
4509396343dSjeremylt   petschavecuda = PETSC_FALSE;
4519396343dSjeremylt   #endif
4529396343dSjeremylt   // *INDENT-ON*
4539396343dSjeremylt 
454cb32e2e7SValeria Barra   ierr = PetscInitialize(&argc, &argv, NULL, help);
455cb32e2e7SValeria Barra   if (ierr) return ierr;
456cb32e2e7SValeria Barra   comm = PETSC_COMM_WORLD;
45732d2ee49SValeria Barra 
45832d2ee49SValeria Barra   // Read command line options
459cb32e2e7SValeria Barra   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr);
460199551f5Sjeremylt   bpchoice = CEED_BP1;
461cb32e2e7SValeria Barra   ierr = PetscOptionsEnum("-problem",
462cb32e2e7SValeria Barra                           "CEED benchmark problem to solve", NULL,
463199551f5Sjeremylt                           bpTypes, (PetscEnum)bpchoice, (PetscEnum *)&bpchoice,
464cb32e2e7SValeria Barra                           NULL); CHKERRQ(ierr);
465199551f5Sjeremylt   ncompu = bpOptions[bpchoice].ncompu;
466cb32e2e7SValeria Barra   test_mode = PETSC_FALSE;
467cb32e2e7SValeria Barra   ierr = PetscOptionsBool("-test",
468cb32e2e7SValeria Barra                           "Testing mode (do not print unless error is large)",
469cb32e2e7SValeria Barra                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
470cb32e2e7SValeria Barra   benchmark_mode = PETSC_FALSE;
471cb32e2e7SValeria Barra   ierr = PetscOptionsBool("-benchmark",
472cb32e2e7SValeria Barra                           "Benchmarking mode (prints benchmark statistics)",
473cb32e2e7SValeria Barra                           NULL, benchmark_mode, &benchmark_mode, NULL);
474cb32e2e7SValeria Barra   CHKERRQ(ierr);
475cb32e2e7SValeria Barra   write_solution = PETSC_FALSE;
476cb32e2e7SValeria Barra   ierr = PetscOptionsBool("-write_solution",
477cb32e2e7SValeria Barra                           "Write solution for visualization",
478cb32e2e7SValeria Barra                           NULL, write_solution, &write_solution, NULL);
479cb32e2e7SValeria Barra   CHKERRQ(ierr);
480cb32e2e7SValeria Barra   degree = test_mode ? 3 : 1;
481cb32e2e7SValeria Barra   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
482cb32e2e7SValeria Barra                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
483199551f5Sjeremylt   qextra = bpOptions[bpchoice].qextra;
484cb32e2e7SValeria Barra   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
485cb32e2e7SValeria Barra                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
486cb32e2e7SValeria Barra   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
487cb32e2e7SValeria Barra                             NULL, ceedresource, ceedresource,
488cb32e2e7SValeria Barra                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
489cb32e2e7SValeria Barra   localnodes = 1000;
490cb32e2e7SValeria Barra   ierr = PetscOptionsInt("-local",
491cb32e2e7SValeria Barra                          "Target number of locally owned nodes per process",
492cb32e2e7SValeria Barra                          NULL, localnodes, &localnodes, NULL); CHKERRQ(ierr);
4939396343dSjeremylt   memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST;
4949396343dSjeremylt   ierr = PetscOptionsEnum("-memtype",
4959396343dSjeremylt                           "CEED MemType requested", NULL,
4969396343dSjeremylt                           memTypes, (PetscEnum)memtyperequested,
4979396343dSjeremylt                           (PetscEnum *)&memtyperequested, &setmemtyperequest);
4989396343dSjeremylt   CHKERRQ(ierr);
499cb32e2e7SValeria Barra   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
500cb32e2e7SValeria Barra   P = degree + 1;
501cb32e2e7SValeria Barra   Q = P + qextra;
502cb32e2e7SValeria Barra 
5039396343dSjeremylt   // Set up libCEED
5049396343dSjeremylt   CeedInit(ceedresource, &ceed);
5059396343dSjeremylt   CeedMemType memtypebackend;
5069396343dSjeremylt   CeedGetPreferredMemType(ceed, &memtypebackend);
5079396343dSjeremylt 
5089396343dSjeremylt   // Check memtype compatibility
5099396343dSjeremylt   if (!setmemtyperequest)
5109396343dSjeremylt     memtyperequested = memtypebackend;
5119396343dSjeremylt   else if (!petschavecuda && memtyperequested == CEED_MEM_DEVICE)
5129396343dSjeremylt     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS,
5139396343dSjeremylt              "PETSc was not built with CUDA. "
5149396343dSjeremylt              "Requested MemType CEED_MEM_DEVICE is not supported.", NULL);
5159396343dSjeremylt 
516cb32e2e7SValeria Barra   // Determine size of process grid
517cb32e2e7SValeria Barra   ierr = MPI_Comm_size(comm, &size); CHKERRQ(ierr);
518cb32e2e7SValeria Barra   Split3(size, p, false);
519cb32e2e7SValeria Barra 
520cb32e2e7SValeria Barra   // Find a nicely composite number of elements no less than localnodes
521cb32e2e7SValeria Barra   for (localelem = PetscMax(1, localnodes / (degree*degree*degree)); ;
522cb32e2e7SValeria Barra        localelem++) {
523cb32e2e7SValeria Barra     Split3(localelem, melem, true);
524cb32e2e7SValeria Barra     if (Max3(melem) / Min3(melem) <= 2) break;
525cb32e2e7SValeria Barra   }
526cb32e2e7SValeria Barra 
527cb32e2e7SValeria Barra   // Find my location in the process grid
528cb32e2e7SValeria Barra   ierr = MPI_Comm_rank(comm, &rank); CHKERRQ(ierr);
529cb32e2e7SValeria Barra   for (int d=0,rankleft=rank; d<dim; d++) {
530cb32e2e7SValeria Barra     const int pstride[3] = {p[1] *p[2], p[2], 1};
531cb32e2e7SValeria Barra     irank[d] = rankleft / pstride[d];
532cb32e2e7SValeria Barra     rankleft -= irank[d] * pstride[d];
533cb32e2e7SValeria Barra   }
534cb32e2e7SValeria Barra 
535cb32e2e7SValeria Barra   GlobalNodes(p, irank, degree, melem, mnodes);
536cb32e2e7SValeria Barra 
537cb32e2e7SValeria Barra   // Setup global vector
5384a5da23aSJed Brown   ierr = VecCreate(comm, &X); CHKERRQ(ierr);
5399396343dSjeremylt   if (memtyperequested == CEED_MEM_DEVICE) {
5409396343dSjeremylt     ierr = VecSetType(X, VECCUDA); CHKERRQ(ierr);
5419396343dSjeremylt   }
542cb32e2e7SValeria Barra   ierr = VecSetSizes(X, mnodes[0]*mnodes[1]*mnodes[2]*ncompu, PETSC_DECIDE);
543cb32e2e7SValeria Barra   CHKERRQ(ierr);
544cb32e2e7SValeria Barra   ierr = VecSetUp(X); CHKERRQ(ierr);
545cb32e2e7SValeria Barra 
546cb32e2e7SValeria Barra   // Set up libCEED
547cb32e2e7SValeria Barra   CeedInit(ceedresource, &ceed);
548cb32e2e7SValeria Barra 
549cb32e2e7SValeria Barra   // Print summary
550cb32e2e7SValeria Barra   CeedInt gsize;
551cb32e2e7SValeria Barra   ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
552dc7d240cSValeria Barra   if (!test_mode) {
553cb32e2e7SValeria Barra     const char *usedresource;
554cb32e2e7SValeria Barra     CeedGetResource(ceed, &usedresource);
5559396343dSjeremylt 
5569396343dSjeremylt     VecType vectype;
5579396343dSjeremylt     ierr = VecGetType(X, &vectype); CHKERRQ(ierr);
5589396343dSjeremylt 
559cb32e2e7SValeria Barra     ierr = PetscPrintf(comm,
560cb32e2e7SValeria Barra                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n"
5619396343dSjeremylt                        "  PETSc:\n"
5629396343dSjeremylt                        "    PETSc Vec Type                     : %s\n"
563cb32e2e7SValeria Barra                        "  libCEED:\n"
564cb32e2e7SValeria Barra                        "    libCEED Backend                    : %s\n"
5659396343dSjeremylt                        "    libCEED Backend MemType            : %s\n"
5669396343dSjeremylt                        "    libCEED User Requested MemType     : %s\n"
567cb32e2e7SValeria Barra                        "  Mesh:\n"
568*8d0bb2bbSvaleriabarra                        "    Number of 1D Basis Nodes (P)       : %d\n"
569*8d0bb2bbSvaleriabarra                        "    Number of 1D Quadrature Points (Q) : %d\n"
570cb32e2e7SValeria Barra                        "    Global nodes                       : %D\n"
571cb32e2e7SValeria Barra                        "    Process Decomposition              : %D %D %D\n"
572cb32e2e7SValeria Barra                        "    Local Elements                     : %D = %D %D %D\n"
573db419314Sjeremylt                        "    Owned nodes                        : %D = %D %D %D\n"
574db419314Sjeremylt                        "    DoF per node                       : %D\n",
5759396343dSjeremylt                        bpchoice+1, vectype, usedresource,
5769396343dSjeremylt                        CeedMemTypes[memtypebackend],
5779396343dSjeremylt                        (setmemtyperequest) ?
5789396343dSjeremylt                        CeedMemTypes[memtyperequested] : "none",
5799396343dSjeremylt                        P, Q,  gsize/ncompu, p[0], p[1], p[2], localelem,
5809396343dSjeremylt                        melem[0], melem[1], melem[2],
581cb32e2e7SValeria Barra                        mnodes[0]*mnodes[1]*mnodes[2], mnodes[0], mnodes[1],
582db419314Sjeremylt                        mnodes[2], ncompu); CHKERRQ(ierr);
583cb32e2e7SValeria Barra   }
584cb32e2e7SValeria Barra 
585cb32e2e7SValeria Barra   {
586cb32e2e7SValeria Barra     lsize = 1;
587cb32e2e7SValeria Barra     for (int d=0; d<dim; d++) {
588cb32e2e7SValeria Barra       lnodes[d] = melem[d]*degree + 1;
589cb32e2e7SValeria Barra       lsize *= lnodes[d];
590cb32e2e7SValeria Barra     }
5914a5da23aSJed Brown     ierr = VecCreate(PETSC_COMM_SELF, &Xloc); CHKERRQ(ierr);
5929396343dSjeremylt     if (memtyperequested == CEED_MEM_DEVICE) {
5939396343dSjeremylt       ierr = VecSetType(Xloc, VECCUDA); CHKERRQ(ierr);
5949396343dSjeremylt     }
595cb32e2e7SValeria Barra     ierr = VecSetSizes(Xloc, lsize*ncompu, PETSC_DECIDE); CHKERRQ(ierr);
596cb32e2e7SValeria Barra     ierr = VecSetUp(Xloc); CHKERRQ(ierr);
597cb32e2e7SValeria Barra 
598cb32e2e7SValeria Barra     // Create local-to-global scatter
599cb32e2e7SValeria Barra     PetscInt *ltogind, *ltogind0, *locind, l0count;
600cb32e2e7SValeria Barra     IS ltogis, ltogis0, locis;
601cb32e2e7SValeria Barra     PetscInt gstart[2][2][2], gmnodes[2][2][2][dim];
602cb32e2e7SValeria Barra 
603cb32e2e7SValeria Barra     for (int i=0; i<2; i++) {
604cb32e2e7SValeria Barra       for (int j=0; j<2; j++) {
605cb32e2e7SValeria Barra         for (int k=0; k<2; k++) {
606cb32e2e7SValeria Barra           PetscInt ijkrank[3] = {irank[0]+i, irank[1]+j, irank[2]+k};
607cb32e2e7SValeria Barra           gstart[i][j][k] = GlobalStart(p, ijkrank, degree, melem);
608cb32e2e7SValeria Barra           GlobalNodes(p, ijkrank, degree, melem, gmnodes[i][j][k]);
609cb32e2e7SValeria Barra         }
610cb32e2e7SValeria Barra       }
611cb32e2e7SValeria Barra     }
612cb32e2e7SValeria Barra 
613cb32e2e7SValeria Barra     ierr = PetscMalloc1(lsize, &ltogind); CHKERRQ(ierr);
614cb32e2e7SValeria Barra     ierr = PetscMalloc1(lsize, &ltogind0); CHKERRQ(ierr);
615cb32e2e7SValeria Barra     ierr = PetscMalloc1(lsize, &locind); CHKERRQ(ierr);
616cb32e2e7SValeria Barra     l0count = 0;
617cb32e2e7SValeria Barra     for (PetscInt i=0,ir,ii; ir=i>=mnodes[0], ii=i-ir*mnodes[0], i<lnodes[0]; i++)
618cb32e2e7SValeria Barra       for (PetscInt j=0,jr,jj; jr=j>=mnodes[1], jj=j-jr*mnodes[1], j<lnodes[1]; j++)
619cb32e2e7SValeria Barra         for (PetscInt k=0,kr,kk; kr=k>=mnodes[2], kk=k-kr*mnodes[2], k<lnodes[2]; k++) {
620cb32e2e7SValeria Barra           PetscInt here = (i*lnodes[1]+j)*lnodes[2]+k;
621cb32e2e7SValeria Barra           ltogind[here] =
622cb32e2e7SValeria Barra             gstart[ir][jr][kr] + (ii*gmnodes[ir][jr][kr][1]+jj)*gmnodes[ir][jr][kr][2]+kk;
623cb32e2e7SValeria Barra           if ((irank[0] == 0 && i == 0)
624cb32e2e7SValeria Barra               || (irank[1] == 0 && j == 0)
625cb32e2e7SValeria Barra               || (irank[2] == 0 && k == 0)
626cb32e2e7SValeria Barra               || (irank[0]+1 == p[0] && i+1 == lnodes[0])
627cb32e2e7SValeria Barra               || (irank[1]+1 == p[1] && j+1 == lnodes[1])
628cb32e2e7SValeria Barra               || (irank[2]+1 == p[2] && k+1 == lnodes[2]))
629cb32e2e7SValeria Barra             continue;
630cb32e2e7SValeria Barra           ltogind0[l0count] = ltogind[here];
631cb32e2e7SValeria Barra           locind[l0count++] = here;
632cb32e2e7SValeria Barra         }
633cb32e2e7SValeria Barra     ierr = ISCreateBlock(comm, ncompu, lsize, ltogind, PETSC_OWN_POINTER,
634cb32e2e7SValeria Barra                          &ltogis); CHKERRQ(ierr);
635cb32e2e7SValeria Barra     ierr = VecScatterCreate(Xloc, NULL, X, ltogis, &ltog); CHKERRQ(ierr);
636cb32e2e7SValeria Barra     CHKERRQ(ierr);
637cb32e2e7SValeria Barra     ierr = ISCreateBlock(comm, ncompu, l0count, ltogind0, PETSC_OWN_POINTER,
638cb32e2e7SValeria Barra                          &ltogis0); CHKERRQ(ierr);
639cb32e2e7SValeria Barra     ierr = ISCreateBlock(comm, ncompu, l0count, locind, PETSC_OWN_POINTER,
640cb32e2e7SValeria Barra                          &locis); CHKERRQ(ierr);
641cb32e2e7SValeria Barra     ierr = VecScatterCreate(Xloc, locis, X, ltogis0, &ltog0); CHKERRQ(ierr);
642cb32e2e7SValeria Barra     {
643cb32e2e7SValeria Barra       // Create global-to-global scatter for Dirichlet values (everything not in
644cb32e2e7SValeria Barra       // ltogis0, which is the range of ltog0)
645cb32e2e7SValeria Barra       PetscInt xstart, xend, *indD, countD = 0;
646cb32e2e7SValeria Barra       IS isD;
647cb32e2e7SValeria Barra       const PetscScalar *x;
648cb32e2e7SValeria Barra       ierr = VecZeroEntries(Xloc); CHKERRQ(ierr);
649cb32e2e7SValeria Barra       ierr = VecSet(X, 1.0); CHKERRQ(ierr);
650cb32e2e7SValeria Barra       ierr = VecScatterBegin(ltog0, Xloc, X, INSERT_VALUES, SCATTER_FORWARD);
651cb32e2e7SValeria Barra       CHKERRQ(ierr);
652cb32e2e7SValeria Barra       ierr = VecScatterEnd(ltog0, Xloc, X, INSERT_VALUES, SCATTER_FORWARD);
653cb32e2e7SValeria Barra       CHKERRQ(ierr);
654cb32e2e7SValeria Barra       ierr = VecGetOwnershipRange(X, &xstart, &xend); CHKERRQ(ierr);
655cb32e2e7SValeria Barra       ierr = PetscMalloc1(xend-xstart, &indD); CHKERRQ(ierr);
656cb32e2e7SValeria Barra       ierr = VecGetArrayRead(X, &x); CHKERRQ(ierr);
657cb32e2e7SValeria Barra       for (PetscInt i=0; i<xend-xstart; i++) {
658cb32e2e7SValeria Barra         if (x[i] == 1.) indD[countD++] = xstart + i;
659cb32e2e7SValeria Barra       }
660cb32e2e7SValeria Barra       ierr = VecRestoreArrayRead(X, &x); CHKERRQ(ierr);
661cb32e2e7SValeria Barra       ierr = ISCreateGeneral(comm, countD, indD, PETSC_COPY_VALUES, &isD);
662cb32e2e7SValeria Barra       CHKERRQ(ierr);
663cb32e2e7SValeria Barra       ierr = PetscFree(indD); CHKERRQ(ierr);
664cb32e2e7SValeria Barra       ierr = VecScatterCreate(X, isD, X, isD, &gtogD); CHKERRQ(ierr);
665cb32e2e7SValeria Barra       ierr = ISDestroy(&isD); CHKERRQ(ierr);
666cb32e2e7SValeria Barra     }
667cb32e2e7SValeria Barra     ierr = ISDestroy(&ltogis); CHKERRQ(ierr);
668cb32e2e7SValeria Barra     ierr = ISDestroy(&ltogis0); CHKERRQ(ierr);
669cb32e2e7SValeria Barra     ierr = ISDestroy(&locis); CHKERRQ(ierr);
670cb32e2e7SValeria Barra   }
671cb32e2e7SValeria Barra 
672cb32e2e7SValeria Barra   // CEED bases
673cb32e2e7SValeria Barra   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompu, P, Q,
674199551f5Sjeremylt                                   bpOptions[bpchoice].qmode, &basisu);
675cb32e2e7SValeria Barra   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, Q,
676199551f5Sjeremylt                                   bpOptions[bpchoice].qmode, &basisx);
677cb32e2e7SValeria Barra 
678cb32e2e7SValeria Barra   // CEED restrictions
679d979a051Sjeremylt   CreateRestriction(ceed, melem, P, ncompu, &Erestrictu);
680d979a051Sjeremylt   CreateRestriction(ceed, melem, 2, dim, &Erestrictx);
681cb32e2e7SValeria Barra   CeedInt nelem = melem[0]*melem[1]*melem[2];
682d979a051Sjeremylt   CeedElemRestrictionCreateStrided(ceed, nelem, Q*Q*Q, ncompu,
683d979a051Sjeremylt                                    ncompu*nelem*Q*Q*Q,
684523b8ea0Sjeremylt                                    CEED_STRIDES_BACKEND, &Erestrictui);
685d979a051Sjeremylt   CeedElemRestrictionCreateStrided(ceed, nelem, Q*Q*Q,
686199551f5Sjeremylt                                    bpOptions[bpchoice].qdatasize,
687d979a051Sjeremylt                                    bpOptions[bpchoice].qdatasize*nelem*Q*Q*Q,
688523b8ea0Sjeremylt                                    CEED_STRIDES_BACKEND, &Erestrictqdi);
689cb32e2e7SValeria Barra   {
690cb32e2e7SValeria Barra     CeedScalar *xloc;
691cb32e2e7SValeria Barra     CeedInt shape[3] = {melem[0]+1, melem[1]+1, melem[2]+1}, len =
692cb32e2e7SValeria Barra                          shape[0]*shape[1]*shape[2];
693cb32e2e7SValeria Barra     xloc = malloc(len*ncompx*sizeof xloc[0]);
694cb32e2e7SValeria Barra     for (CeedInt i=0; i<shape[0]; i++) {
695cb32e2e7SValeria Barra       for (CeedInt j=0; j<shape[1]; j++) {
696cb32e2e7SValeria Barra         for (CeedInt k=0; k<shape[2]; k++) {
697d979a051Sjeremylt           xloc[dim*((i*shape[1]+j)*shape[2]+k) + 0] = 1.*(irank[0]*melem[0]+i) /
698cb32e2e7SValeria Barra               (p[0]*melem[0]);
699d979a051Sjeremylt           xloc[dim*((i*shape[1]+j)*shape[2]+k) + 1] = 1.*(irank[1]*melem[1]+j) /
700cb32e2e7SValeria Barra               (p[1]*melem[1]);
701d979a051Sjeremylt           xloc[dim*((i*shape[1]+j)*shape[2]+k) + 2] = 1.*(irank[2]*melem[2]+k) /
702cb32e2e7SValeria Barra               (p[2]*melem[2]);
703cb32e2e7SValeria Barra         }
704cb32e2e7SValeria Barra       }
705cb32e2e7SValeria Barra     }
706cb32e2e7SValeria Barra     CeedVectorCreate(ceed, len*ncompx, &xcoord);
707cb32e2e7SValeria Barra     CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_OWN_POINTER, xloc);
708cb32e2e7SValeria Barra   }
709cb32e2e7SValeria Barra 
710cb32e2e7SValeria Barra   // Create the Qfunction that builds the operator quadrature data
711199551f5Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].setupgeo,
712199551f5Sjeremylt                               bpOptions[bpchoice].setupgeofname, &qfsetupgeo);
713c70bd2a0Sjeremylt   CeedQFunctionAddInput(qfsetupgeo, "dx", ncompx*dim, CEED_EVAL_GRAD);
714c70bd2a0Sjeremylt   CeedQFunctionAddInput(qfsetupgeo, "weight", 1, CEED_EVAL_WEIGHT);
715199551f5Sjeremylt   CeedQFunctionAddOutput(qfsetupgeo, "qdata", bpOptions[bpchoice].qdatasize,
716cb32e2e7SValeria Barra                          CEED_EVAL_NONE);
717cb32e2e7SValeria Barra 
718cb32e2e7SValeria Barra   // Create the Qfunction that sets up the RHS and true solution
719199551f5Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].setuprhs,
720199551f5Sjeremylt                               bpOptions[bpchoice].setuprhsfname, &qfsetuprhs);
721c70bd2a0Sjeremylt   CeedQFunctionAddInput(qfsetuprhs, "x", ncompx, CEED_EVAL_INTERP);
722c70bd2a0Sjeremylt   CeedQFunctionAddInput(qfsetuprhs, "dx", ncompx*dim, CEED_EVAL_GRAD);
723c70bd2a0Sjeremylt   CeedQFunctionAddInput(qfsetuprhs, "weight", 1, CEED_EVAL_WEIGHT);
724c70bd2a0Sjeremylt   CeedQFunctionAddOutput(qfsetuprhs, "true_soln", ncompu, CEED_EVAL_NONE);
725c70bd2a0Sjeremylt   CeedQFunctionAddOutput(qfsetuprhs, "rhs", ncompu, CEED_EVAL_INTERP);
726cb32e2e7SValeria Barra 
727cb32e2e7SValeria Barra   // Set up PDE operator
728199551f5Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].apply,
729199551f5Sjeremylt                               bpOptions[bpchoice].applyfname, &qfapply);
730cb32e2e7SValeria Barra   // Add inputs and outputs
731199551f5Sjeremylt   CeedInt inscale = bpOptions[bpchoice].inmode==CEED_EVAL_GRAD ? 3 : 1;
732199551f5Sjeremylt   CeedInt outscale = bpOptions[bpchoice].outmode==CEED_EVAL_GRAD ? 3 : 1;
733c70bd2a0Sjeremylt   CeedQFunctionAddInput(qfapply, "u", ncompu*inscale,
734199551f5Sjeremylt                         bpOptions[bpchoice].inmode);
735199551f5Sjeremylt   CeedQFunctionAddInput(qfapply, "qdata", bpOptions[bpchoice].qdatasize,
736cb32e2e7SValeria Barra                         CEED_EVAL_NONE);
737c70bd2a0Sjeremylt   CeedQFunctionAddOutput(qfapply, "v", ncompu*outscale,
738199551f5Sjeremylt                          bpOptions[bpchoice].outmode);
739cb32e2e7SValeria Barra 
740cb32e2e7SValeria Barra   // Create the error qfunction
741199551f5Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].error,
742199551f5Sjeremylt                               bpOptions[bpchoice].errorfname, &qferror);
743c70bd2a0Sjeremylt   CeedQFunctionAddInput(qferror, "u", ncompu, CEED_EVAL_INTERP);
744c70bd2a0Sjeremylt   CeedQFunctionAddInput(qferror, "true_soln", ncompu, CEED_EVAL_NONE);
745c70bd2a0Sjeremylt   CeedQFunctionAddOutput(qferror, "error", ncompu, CEED_EVAL_NONE);
746cb32e2e7SValeria Barra 
747cb32e2e7SValeria Barra   // Create the persistent vectors that will be needed in setup
748cb32e2e7SValeria Barra   CeedInt nqpts;
749cb32e2e7SValeria Barra   CeedBasisGetNumQuadraturePoints(basisu, &nqpts);
750199551f5Sjeremylt   CeedVectorCreate(ceed, bpOptions[bpchoice].qdatasize*nelem*nqpts, &qdata);
751cb32e2e7SValeria Barra   CeedVectorCreate(ceed, nelem*nqpts*ncompu, &target);
752cb32e2e7SValeria Barra   CeedVectorCreate(ceed, lsize*ncompu, &rhsceed);
753cb32e2e7SValeria Barra 
754cb32e2e7SValeria Barra   // Create the operator that builds the quadrature data for the ceed operator
755c70bd2a0Sjeremylt   CeedOperatorCreate(ceed, qfsetupgeo, CEED_QFUNCTION_NONE,
756c70bd2a0Sjeremylt                      CEED_QFUNCTION_NONE, &opsetupgeo);
757c70bd2a0Sjeremylt   CeedOperatorSetField(opsetupgeo, "dx", Erestrictx, basisx,
758a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
759c70bd2a0Sjeremylt   CeedOperatorSetField(opsetupgeo, "weight", CEED_ELEMRESTRICTION_NONE, basisx,
760a8d32208Sjeremylt                        CEED_VECTOR_NONE);
761c70bd2a0Sjeremylt   CeedOperatorSetField(opsetupgeo, "qdata", Erestrictqdi,
762cb32e2e7SValeria Barra                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
763cb32e2e7SValeria Barra 
764cb32e2e7SValeria Barra   // Create the operator that builds the RHS and true solution
765c70bd2a0Sjeremylt   CeedOperatorCreate(ceed, qfsetuprhs, CEED_QFUNCTION_NONE,
766c70bd2a0Sjeremylt                      CEED_QFUNCTION_NONE, &opsetuprhs);
767c70bd2a0Sjeremylt   CeedOperatorSetField(opsetuprhs, "x", Erestrictx, basisx,
768a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
769c70bd2a0Sjeremylt   CeedOperatorSetField(opsetuprhs, "dx", Erestrictx, basisx,
770a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
771c70bd2a0Sjeremylt   CeedOperatorSetField(opsetuprhs, "weight", CEED_ELEMRESTRICTION_NONE, basisx,
772a8d32208Sjeremylt                        CEED_VECTOR_NONE);
773c70bd2a0Sjeremylt   CeedOperatorSetField(opsetuprhs, "true_soln", Erestrictui,
774cb32e2e7SValeria Barra                        CEED_BASIS_COLLOCATED, target);
775c70bd2a0Sjeremylt   CeedOperatorSetField(opsetuprhs, "rhs", Erestrictu, basisu,
776a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
777cb32e2e7SValeria Barra 
778cb32e2e7SValeria Barra   // Create the mass or diff operator
779c70bd2a0Sjeremylt   CeedOperatorCreate(ceed, qfapply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
780c70bd2a0Sjeremylt                      &opapply);
781c70bd2a0Sjeremylt   CeedOperatorSetField(opapply, "u", Erestrictu, basisu, CEED_VECTOR_ACTIVE);
782c70bd2a0Sjeremylt   CeedOperatorSetField(opapply, "qdata", Erestrictqdi, CEED_BASIS_COLLOCATED,
783a8d32208Sjeremylt                        qdata);
784c70bd2a0Sjeremylt   CeedOperatorSetField(opapply, "v", Erestrictu, basisu, CEED_VECTOR_ACTIVE);
785cb32e2e7SValeria Barra 
786cb32e2e7SValeria Barra   // Create the error operator
787c70bd2a0Sjeremylt   CeedOperatorCreate(ceed, qferror, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
788c70bd2a0Sjeremylt                      &operror);
789c70bd2a0Sjeremylt   CeedOperatorSetField(operror, "u", Erestrictu, basisu, CEED_VECTOR_ACTIVE);
790c70bd2a0Sjeremylt   CeedOperatorSetField(operror, "true_soln", Erestrictui,
791cb32e2e7SValeria Barra                        CEED_BASIS_COLLOCATED, target);
792c70bd2a0Sjeremylt   CeedOperatorSetField(operror, "error", Erestrictui, CEED_BASIS_COLLOCATED,
793a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
794cb32e2e7SValeria Barra 
795cb32e2e7SValeria Barra   // Set up Mat
796cb32e2e7SValeria Barra   ierr = PetscMalloc1(1, &user); CHKERRQ(ierr);
797cb32e2e7SValeria Barra   user->comm = comm;
798cb32e2e7SValeria Barra   user->ltog = ltog;
799199551f5Sjeremylt   if (bpchoice != CEED_BP1 && bpchoice != CEED_BP2) {
800cb32e2e7SValeria Barra     user->ltog0 = ltog0;
801cb32e2e7SValeria Barra     user->gtogD = gtogD;
802cb32e2e7SValeria Barra   }
803cb32e2e7SValeria Barra   user->Xloc = Xloc;
804cb32e2e7SValeria Barra   ierr = VecDuplicate(Xloc, &user->Yloc); CHKERRQ(ierr);
805cb32e2e7SValeria Barra   CeedVectorCreate(ceed, lsize*ncompu, &user->xceed);
806cb32e2e7SValeria Barra   CeedVectorCreate(ceed, lsize*ncompu, &user->yceed);
807c70bd2a0Sjeremylt   user->op = opapply;
808cb32e2e7SValeria Barra   user->qdata = qdata;
809cb32e2e7SValeria Barra   user->ceed = ceed;
8109396343dSjeremylt   user->memtype = memtyperequested;
8119396343dSjeremylt   if (memtyperequested == CEED_MEM_HOST) {
8129396343dSjeremylt     user->VecGetArray = VecGetArray;
8139396343dSjeremylt     user->VecGetArrayRead = VecGetArrayRead;
8149396343dSjeremylt     user->VecRestoreArray = VecRestoreArray;
8159396343dSjeremylt     user->VecRestoreArrayRead = VecRestoreArrayRead;
8169396343dSjeremylt   } else {
8179396343dSjeremylt     user->VecGetArray = VecCUDAGetArray;
8189396343dSjeremylt     user->VecGetArrayRead = VecCUDAGetArrayRead;
8199396343dSjeremylt     user->VecRestoreArray = VecCUDARestoreArray;
8209396343dSjeremylt     user->VecRestoreArrayRead = VecCUDARestoreArrayRead;
8219396343dSjeremylt   }
822cb32e2e7SValeria Barra 
823cb32e2e7SValeria Barra   ierr = MatCreateShell(comm, mnodes[0]*mnodes[1]*mnodes[2]*ncompu,
824cb32e2e7SValeria Barra                         mnodes[0]*mnodes[1]*mnodes[2]*ncompu,
825cb32e2e7SValeria Barra                         PETSC_DECIDE, PETSC_DECIDE, user, &mat); CHKERRQ(ierr);
826199551f5Sjeremylt   if (bpchoice == CEED_BP1 || bpchoice == CEED_BP2) {
827cb32e2e7SValeria Barra     ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Mass);
828cb32e2e7SValeria Barra     CHKERRQ(ierr);
829cb32e2e7SValeria Barra   } else {
830cb32e2e7SValeria Barra     ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Diff);
831cb32e2e7SValeria Barra     CHKERRQ(ierr);
832cb32e2e7SValeria Barra   }
8339396343dSjeremylt   if (user->memtype == CEED_MEM_DEVICE) {
8349396343dSjeremylt     ierr = MatShellSetVecType(mat, VECCUDA); CHKERRQ(ierr);
8359396343dSjeremylt   }
836cb32e2e7SValeria Barra 
837cb32e2e7SValeria Barra   // Get RHS vector
8389396343dSjeremylt   ierr = VecDuplicate(X, &rhs); CHKERRQ(ierr);
839cb32e2e7SValeria Barra   ierr = VecDuplicate(Xloc, &rhsloc); CHKERRQ(ierr);
840cb32e2e7SValeria Barra   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
8419396343dSjeremylt   ierr = user->VecGetArray(rhsloc, &r); CHKERRQ(ierr);
8429396343dSjeremylt   CeedVectorSetArray(rhsceed, user->memtype, CEED_USE_POINTER, r);
843cb32e2e7SValeria Barra 
844cb32e2e7SValeria Barra   // Setup qdata, rhs, and target
845c70bd2a0Sjeremylt   CeedOperatorApply(opsetupgeo, xcoord, qdata, CEED_REQUEST_IMMEDIATE);
846c70bd2a0Sjeremylt   CeedOperatorApply(opsetuprhs, xcoord, rhsceed, CEED_REQUEST_IMMEDIATE);
8479396343dSjeremylt   ierr = CeedVectorSyncArray(rhsceed, user->memtype); CHKERRQ(ierr);
848cb32e2e7SValeria Barra   CeedVectorDestroy(&xcoord);
849cb32e2e7SValeria Barra 
850cb32e2e7SValeria Barra   // Gather RHS
8519396343dSjeremylt   ierr = user->VecRestoreArray(rhsloc, &r); CHKERRQ(ierr);
852cb32e2e7SValeria Barra   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
853cb32e2e7SValeria Barra   ierr = VecScatterBegin(ltog, rhsloc, rhs, ADD_VALUES, SCATTER_FORWARD);
854cb32e2e7SValeria Barra   CHKERRQ(ierr);
855cb32e2e7SValeria Barra   ierr = VecScatterEnd(ltog, rhsloc, rhs, ADD_VALUES, SCATTER_FORWARD);
856cb32e2e7SValeria Barra   CHKERRQ(ierr);
857cb32e2e7SValeria Barra   CeedVectorDestroy(&rhsceed);
858cb32e2e7SValeria Barra 
859cb32e2e7SValeria Barra   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
860cb32e2e7SValeria Barra   {
861cb32e2e7SValeria Barra     PC pc;
862cb32e2e7SValeria Barra     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
863199551f5Sjeremylt     if (bpchoice == CEED_BP1 || bpchoice == CEED_BP2) {
864cb32e2e7SValeria Barra       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
865cb32e2e7SValeria Barra       ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr);
866cb32e2e7SValeria Barra     } else {
867cb32e2e7SValeria Barra       ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr);
868cb32e2e7SValeria Barra     }
869cb32e2e7SValeria Barra     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
870cb32e2e7SValeria Barra     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
871cb32e2e7SValeria Barra     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
872cb32e2e7SValeria Barra                             PETSC_DEFAULT); CHKERRQ(ierr);
873cb32e2e7SValeria Barra   }
874cb32e2e7SValeria Barra   ierr = KSPSetOperators(ksp, mat, mat); CHKERRQ(ierr);
875cb32e2e7SValeria Barra   // First run, if benchmarking
876cb32e2e7SValeria Barra   if (benchmark_mode) {
877cb32e2e7SValeria Barra     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
878cb32e2e7SValeria Barra     CHKERRQ(ierr);
879cb32e2e7SValeria Barra     my_rt_start = MPI_Wtime();
880cb32e2e7SValeria Barra     ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
881cb32e2e7SValeria Barra     my_rt = MPI_Wtime() - my_rt_start;
882cb32e2e7SValeria Barra     ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
883cb32e2e7SValeria Barra     CHKERRQ(ierr);
884cb32e2e7SValeria Barra     // Set maxits based on first iteration timing
885cb32e2e7SValeria Barra     if (my_rt > 0.02) {
886cb32e2e7SValeria Barra       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
887cb32e2e7SValeria Barra       CHKERRQ(ierr);
888cb32e2e7SValeria Barra     } else {
889cb32e2e7SValeria Barra       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
890cb32e2e7SValeria Barra       CHKERRQ(ierr);
891cb32e2e7SValeria Barra     }
892cb32e2e7SValeria Barra   }
893debcf919SJed Brown   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
89409a940d7Sjeremylt 
895cb32e2e7SValeria Barra   // Timed solve
89609a940d7Sjeremylt   ierr = VecZeroEntries(X); CHKERRQ(ierr);
897cb32e2e7SValeria Barra   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
89809a940d7Sjeremylt 
89909a940d7Sjeremylt   // -- Performance logging
90009a940d7Sjeremylt   ierr = PetscLogStageRegister("Solve Stage", &solvestage); CHKERRQ(ierr);
90109a940d7Sjeremylt   ierr = PetscLogStagePush(solvestage); CHKERRQ(ierr);
90209a940d7Sjeremylt 
90309a940d7Sjeremylt   // -- Solve
904cb32e2e7SValeria Barra   my_rt_start = MPI_Wtime();
905cb32e2e7SValeria Barra   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
906cb32e2e7SValeria Barra   my_rt = MPI_Wtime() - my_rt_start;
90709a940d7Sjeremylt 
90809a940d7Sjeremylt   // -- Performance logging
90909a940d7Sjeremylt   ierr = PetscLogStagePop();
91009a940d7Sjeremylt 
9118e87e98bSJed Brown   // Output results
912cb32e2e7SValeria Barra   {
913cb32e2e7SValeria Barra     KSPType ksptype;
914cb32e2e7SValeria Barra     KSPConvergedReason reason;
915cb32e2e7SValeria Barra     PetscReal rnorm;
916cb32e2e7SValeria Barra     PetscInt its;
917cb32e2e7SValeria Barra     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
918cb32e2e7SValeria Barra     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
919cb32e2e7SValeria Barra     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
920cb32e2e7SValeria Barra     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
921cb32e2e7SValeria Barra     if (!test_mode || reason < 0 || rnorm > 1e-8) {
922cb32e2e7SValeria Barra       ierr = PetscPrintf(comm,
923cb32e2e7SValeria Barra                          "  KSP:\n"
924cb32e2e7SValeria Barra                          "    KSP Type                           : %s\n"
925cb32e2e7SValeria Barra                          "    KSP Convergence                    : %s\n"
926cb32e2e7SValeria Barra                          "    Total KSP Iterations               : %D\n"
927cb32e2e7SValeria Barra                          "    Final rnorm                        : %e\n",
928cb32e2e7SValeria Barra                          ksptype, KSPConvergedReasons[reason], its,
929cb32e2e7SValeria Barra                          (double)rnorm); CHKERRQ(ierr);
930cb32e2e7SValeria Barra     }
9318e87e98bSJed Brown     if (!test_mode) {
9328e87e98bSJed Brown       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
9338e87e98bSJed Brown     }
9348e87e98bSJed Brown     {
9358e87e98bSJed Brown       PetscReal maxerror;
9368e87e98bSJed Brown       ierr = ComputeErrorMax(user, operror, X, target, &maxerror);
9378e87e98bSJed Brown       CHKERRQ(ierr);
9388e87e98bSJed Brown       PetscReal tol = 5e-2;
9398e87e98bSJed Brown       if (!test_mode || maxerror > tol) {
940cb32e2e7SValeria Barra         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
941cb32e2e7SValeria Barra         CHKERRQ(ierr);
942cb32e2e7SValeria Barra         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
943cb32e2e7SValeria Barra         CHKERRQ(ierr);
944cb32e2e7SValeria Barra         ierr = PetscPrintf(comm,
9458e87e98bSJed Brown                            "    Pointwise Error (max)              : %e\n"
9468e87e98bSJed Brown                            "    CG Solve Time                      : %g (%g) sec\n",
9478e87e98bSJed Brown                            (double)maxerror, rt_max, rt_min); CHKERRQ(ierr);
948cb32e2e7SValeria Barra       }
949cb32e2e7SValeria Barra     }
950*8d0bb2bbSvaleriabarra     if (!test_mode) {
951cb32e2e7SValeria Barra       ierr = PetscPrintf(comm,
9528e87e98bSJed Brown                          "    DoFs/Sec in CG                     : %g (%g) million\n",
9538e87e98bSJed Brown                          1e-6*gsize*its/rt_max,
9548e87e98bSJed Brown                          1e-6*gsize*its/rt_min); CHKERRQ(ierr);
955cb32e2e7SValeria Barra     }
956cb32e2e7SValeria Barra   }
957cb32e2e7SValeria Barra 
958cb32e2e7SValeria Barra   if (write_solution) {
959cb32e2e7SValeria Barra     PetscViewer vtkviewersoln;
960cb32e2e7SValeria Barra 
961cb32e2e7SValeria Barra     ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr);
962cb32e2e7SValeria Barra     ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr);
963cb32e2e7SValeria Barra     ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr);
964cb32e2e7SValeria Barra     ierr = VecView(X, vtkviewersoln); CHKERRQ(ierr);
965cb32e2e7SValeria Barra     ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr);
966cb32e2e7SValeria Barra   }
967cb32e2e7SValeria Barra 
968cb32e2e7SValeria Barra   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
969cb32e2e7SValeria Barra   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
970cb32e2e7SValeria Barra   ierr = VecDestroy(&X); CHKERRQ(ierr);
971cb32e2e7SValeria Barra   ierr = VecDestroy(&user->Xloc); CHKERRQ(ierr);
972cb32e2e7SValeria Barra   ierr = VecDestroy(&user->Yloc); CHKERRQ(ierr);
973cb32e2e7SValeria Barra   ierr = VecScatterDestroy(&ltog); CHKERRQ(ierr);
974cb32e2e7SValeria Barra   ierr = VecScatterDestroy(&ltog0); CHKERRQ(ierr);
975cb32e2e7SValeria Barra   ierr = VecScatterDestroy(&gtogD); CHKERRQ(ierr);
976cb32e2e7SValeria Barra   ierr = MatDestroy(&mat); CHKERRQ(ierr);
977cb32e2e7SValeria Barra   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
978cb32e2e7SValeria Barra 
979cb32e2e7SValeria Barra   CeedVectorDestroy(&user->xceed);
980cb32e2e7SValeria Barra   CeedVectorDestroy(&user->yceed);
981cb32e2e7SValeria Barra   CeedVectorDestroy(&user->qdata);
982cb32e2e7SValeria Barra   CeedVectorDestroy(&target);
983c70bd2a0Sjeremylt   CeedOperatorDestroy(&opsetupgeo);
984c70bd2a0Sjeremylt   CeedOperatorDestroy(&opsetuprhs);
985c70bd2a0Sjeremylt   CeedOperatorDestroy(&opapply);
986c70bd2a0Sjeremylt   CeedOperatorDestroy(&operror);
987cb32e2e7SValeria Barra   CeedElemRestrictionDestroy(&Erestrictu);
988cb32e2e7SValeria Barra   CeedElemRestrictionDestroy(&Erestrictx);
989cb32e2e7SValeria Barra   CeedElemRestrictionDestroy(&Erestrictui);
990cb32e2e7SValeria Barra   CeedElemRestrictionDestroy(&Erestrictqdi);
991c70bd2a0Sjeremylt   CeedQFunctionDestroy(&qfsetupgeo);
992c70bd2a0Sjeremylt   CeedQFunctionDestroy(&qfsetuprhs);
993c70bd2a0Sjeremylt   CeedQFunctionDestroy(&qfapply);
994c70bd2a0Sjeremylt   CeedQFunctionDestroy(&qferror);
995cb32e2e7SValeria Barra   CeedBasisDestroy(&basisu);
996cb32e2e7SValeria Barra   CeedBasisDestroy(&basisx);
997cb32e2e7SValeria Barra   CeedDestroy(&ceed);
998cb32e2e7SValeria Barra   ierr = PetscFree(user); CHKERRQ(ierr);
999cb32e2e7SValeria Barra   return PetscFinalize();
1000cb32e2e7SValeria Barra }
1001