xref: /libCEED/examples/petsc/bpsraw.c (revision b68a8d799acb1d44569fb95028e25f895284bc04)
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
3228688798Sjeremylt //     ./bpsraw -problem bp2
3328688798Sjeremylt //     ./bpsraw -problem bp3
3428688798Sjeremylt //     ./bpsraw -problem bp4
3528688798Sjeremylt //     ./bpsraw -problem bp5 -ceed /cpu/self
3628688798Sjeremylt //     ./bpsraw -problem bp6 -ceed /gpu/cuda
37cb32e2e7SValeria Barra //
38f9342789SJeremy L Thompson //TESTARGS -ceed {ceed_resource} -test -problem bp2 -degree 5 -qextra 1 -ksp_max_it_clip 15,15
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 
64*b68a8d79SJed Brown static CeedMemType MemTypeP2C(PetscMemType mtype) {
65*b68a8d79SJed Brown   return PetscMemTypeDevice(mtype) ? CEED_MEM_DEVICE : CEED_MEM_HOST;
66*b68a8d79SJed Brown }
67*b68a8d79SJed Brown 
68cb32e2e7SValeria Barra static void Split3(PetscInt size, PetscInt m[3], bool reverse) {
69cb32e2e7SValeria Barra   for (PetscInt d=0,sizeleft=size; d<3; d++) {
70cb32e2e7SValeria Barra     PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(sizeleft, 1./(3 - d)));
71cb32e2e7SValeria Barra     while (try * (sizeleft / try) != sizeleft) try++;
72cb32e2e7SValeria Barra     m[reverse ? 2-d : d] = try;
73cb32e2e7SValeria Barra     sizeleft /= try;
74cb32e2e7SValeria Barra   }
75cb32e2e7SValeria Barra }
76cb32e2e7SValeria Barra 
77cb32e2e7SValeria Barra static PetscInt Max3(const PetscInt a[3]) {
78cb32e2e7SValeria Barra   return PetscMax(a[0], PetscMax(a[1], a[2]));
79cb32e2e7SValeria Barra }
80cb32e2e7SValeria Barra static PetscInt Min3(const PetscInt a[3]) {
81cb32e2e7SValeria Barra   return PetscMin(a[0], PetscMin(a[1], a[2]));
82cb32e2e7SValeria Barra }
83cb32e2e7SValeria Barra static void GlobalNodes(const PetscInt p[3], const PetscInt irank[3],
84cb32e2e7SValeria Barra                         PetscInt degree, const PetscInt melem[3],
85cb32e2e7SValeria Barra                         PetscInt mnodes[3]) {
86cb32e2e7SValeria Barra   for (int d=0; d<3; d++)
87cb32e2e7SValeria Barra     mnodes[d] = degree*melem[d] + (irank[d] == p[d]-1);
88cb32e2e7SValeria Barra }
89cb32e2e7SValeria Barra static PetscInt GlobalStart(const PetscInt p[3], const PetscInt irank[3],
90cb32e2e7SValeria Barra                             PetscInt degree, const PetscInt melem[3]) {
91cb32e2e7SValeria Barra   PetscInt start = 0;
92cb32e2e7SValeria Barra   // Dumb brute-force is easier to read
93cb32e2e7SValeria Barra   for (PetscInt i=0; i<p[0]; i++) {
94cb32e2e7SValeria Barra     for (PetscInt j=0; j<p[1]; j++) {
95cb32e2e7SValeria Barra       for (PetscInt k=0; k<p[2]; k++) {
96cb32e2e7SValeria Barra         PetscInt mnodes[3], ijkrank[] = {i,j,k};
97cb32e2e7SValeria Barra         if (i == irank[0] && j == irank[1] && k == irank[2]) return start;
98cb32e2e7SValeria Barra         GlobalNodes(p, ijkrank, degree, melem, mnodes);
99cb32e2e7SValeria Barra         start += mnodes[0] * mnodes[1] * mnodes[2];
100cb32e2e7SValeria Barra       }
101cb32e2e7SValeria Barra     }
102cb32e2e7SValeria Barra   }
103cb32e2e7SValeria Barra   return -1;
104cb32e2e7SValeria Barra }
105d979a051Sjeremylt static int CreateRestriction(Ceed ceed, const CeedInt melem[3], CeedInt P,
106d979a051Sjeremylt                              CeedInt ncomp, CeedElemRestriction *Erestrict) {
107cb32e2e7SValeria Barra   const PetscInt nelem = melem[0]*melem[1]*melem[2];
108cb32e2e7SValeria Barra   PetscInt mnodes[3], *idx, *idxp;
109cb32e2e7SValeria Barra 
110cb32e2e7SValeria Barra   // Get indicies
111cb32e2e7SValeria Barra   for (int d=0; d<3; d++) mnodes[d] = melem[d]*(P-1) + 1;
112cb32e2e7SValeria Barra   idxp = idx = malloc(nelem*P*P*P*sizeof idx[0]);
113cb32e2e7SValeria Barra   for (CeedInt i=0; i<melem[0]; i++)
114cb32e2e7SValeria Barra     for (CeedInt j=0; j<melem[1]; j++)
115cb32e2e7SValeria Barra       for (CeedInt k=0; k<melem[2]; k++,idxp += P*P*P)
116cb32e2e7SValeria Barra         for (CeedInt ii=0; ii<P; ii++)
117cb32e2e7SValeria Barra           for (CeedInt jj=0; jj<P; jj++)
118cb32e2e7SValeria Barra             for (CeedInt kk=0; kk<P; kk++) {
119cb32e2e7SValeria Barra               if (0) { // This is the C-style (i,j,k) ordering that I prefer
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               } else { // (k,j,i) ordering for consistency with MFEM example
124d979a051Sjeremylt                 idxp[ii+P*(jj+P*kk)] = ncomp*(((i*(P-1)+ii)*mnodes[1]
125cb32e2e7SValeria Barra                                                + (j*(P-1)+jj))*mnodes[2]
126cb32e2e7SValeria Barra                                               + (k*(P-1)+kk));
127cb32e2e7SValeria Barra               }
128cb32e2e7SValeria Barra             }
129cb32e2e7SValeria Barra 
130cb32e2e7SValeria Barra   // Setup CEED restriction
131d979a051Sjeremylt   CeedElemRestrictionCreate(ceed, nelem, P*P*P, ncomp, 1,
132d979a051Sjeremylt                             mnodes[0]*mnodes[1]*mnodes[2]*ncomp,
133a8d32208Sjeremylt                             CEED_MEM_HOST, CEED_OWN_POINTER, idx, Erestrict);
134cb32e2e7SValeria Barra 
135cb32e2e7SValeria Barra   PetscFunctionReturn(0);
136cb32e2e7SValeria Barra }
137cb32e2e7SValeria Barra 
138cb32e2e7SValeria Barra // Data for PETSc
139cb32e2e7SValeria Barra typedef struct User_ *User;
140cb32e2e7SValeria Barra struct User_ {
141cb32e2e7SValeria Barra   MPI_Comm comm;
142cb32e2e7SValeria Barra   VecScatter ltog;              // Scatter for all entries
143cb32e2e7SValeria Barra   VecScatter ltog0;             // Skip Dirichlet values
144cb32e2e7SValeria Barra   VecScatter gtogD;             // global-to-global; only Dirichlet values
145cb32e2e7SValeria Barra   Vec Xloc, Yloc;
146cb32e2e7SValeria Barra   CeedVector xceed, yceed;
147cb32e2e7SValeria Barra   CeedOperator op;
148cb32e2e7SValeria Barra   CeedVector qdata;
149cb32e2e7SValeria Barra   Ceed ceed;
150cb32e2e7SValeria Barra };
151cb32e2e7SValeria Barra 
152cb32e2e7SValeria Barra // BP Options
153cb32e2e7SValeria Barra typedef enum {
154cb32e2e7SValeria Barra   CEED_BP1 = 0, CEED_BP2 = 1, CEED_BP3 = 2,
155cb32e2e7SValeria Barra   CEED_BP4 = 3, CEED_BP5 = 4, CEED_BP6 = 5
156cb32e2e7SValeria Barra } bpType;
157cb32e2e7SValeria Barra static const char *const bpTypes[] = {"bp1","bp2","bp3","bp4","bp5","bp6",
158cb32e2e7SValeria Barra                                       "bpType","CEED_BP",0
159cb32e2e7SValeria Barra                                      };
160cb32e2e7SValeria Barra 
161cb32e2e7SValeria Barra // BP specific data
162cb32e2e7SValeria Barra typedef struct {
163cb32e2e7SValeria Barra   CeedInt ncompu, qdatasize, qextra;
164cb32e2e7SValeria Barra   CeedQFunctionUser setupgeo, setuprhs, apply, error;
165cb32e2e7SValeria Barra   const char *setupgeofname, *setuprhsfname, *applyfname, *errorfname;
166cb32e2e7SValeria Barra   CeedEvalMode inmode, outmode;
167cb32e2e7SValeria Barra   CeedQuadMode qmode;
168cb32e2e7SValeria Barra } bpData;
169cb32e2e7SValeria Barra 
170cb32e2e7SValeria Barra bpData bpOptions[6] = {
171cb32e2e7SValeria Barra   [CEED_BP1] = {
172cb32e2e7SValeria Barra     .ncompu = 1,
173cb32e2e7SValeria Barra     .qdatasize = 1,
174cb32e2e7SValeria Barra     .qextra = 1,
175cb32e2e7SValeria Barra     .setupgeo = SetupMassGeo,
176cb32e2e7SValeria Barra     .setuprhs = SetupMassRhs,
177cb32e2e7SValeria Barra     .apply = Mass,
178cb32e2e7SValeria Barra     .error = Error,
179cb32e2e7SValeria Barra     .setupgeofname = SetupMassGeo_loc,
180cb32e2e7SValeria Barra     .setuprhsfname = SetupMassRhs_loc,
181cb32e2e7SValeria Barra     .applyfname = Mass_loc,
182cb32e2e7SValeria Barra     .errorfname = Error_loc,
183cb32e2e7SValeria Barra     .inmode = CEED_EVAL_INTERP,
184cb32e2e7SValeria Barra     .outmode = CEED_EVAL_INTERP,
185cb32e2e7SValeria Barra     .qmode = CEED_GAUSS
186cb32e2e7SValeria Barra   },
187cb32e2e7SValeria Barra   [CEED_BP2] = {
188cb32e2e7SValeria Barra     .ncompu = 3,
189cb32e2e7SValeria Barra     .qdatasize = 1,
190cb32e2e7SValeria Barra     .qextra = 1,
191cb32e2e7SValeria Barra     .setupgeo = SetupMassGeo,
192cb32e2e7SValeria Barra     .setuprhs = SetupMassRhs3,
193cb32e2e7SValeria Barra     .apply = Mass3,
194cb32e2e7SValeria Barra     .error = Error3,
195cb32e2e7SValeria Barra     .setupgeofname = SetupMassGeo_loc,
196cb32e2e7SValeria Barra     .setuprhsfname = SetupMassRhs3_loc,
197cb32e2e7SValeria Barra     .applyfname = Mass3_loc,
198cb32e2e7SValeria Barra     .errorfname = Error3_loc,
199cb32e2e7SValeria Barra     .inmode = CEED_EVAL_INTERP,
200cb32e2e7SValeria Barra     .outmode = CEED_EVAL_INTERP,
201cb32e2e7SValeria Barra     .qmode = CEED_GAUSS
202cb32e2e7SValeria Barra   },
203cb32e2e7SValeria Barra   [CEED_BP3] = {
204cb32e2e7SValeria Barra     .ncompu = 1,
205cb32e2e7SValeria Barra     .qdatasize = 6,
206cb32e2e7SValeria Barra     .qextra = 1,
207cb32e2e7SValeria Barra     .setupgeo = SetupDiffGeo,
208cb32e2e7SValeria Barra     .setuprhs = SetupDiffRhs,
209cb32e2e7SValeria Barra     .apply = Diff,
210cb32e2e7SValeria Barra     .error = Error,
211cb32e2e7SValeria Barra     .setupgeofname = SetupDiffGeo_loc,
212cb32e2e7SValeria Barra     .setuprhsfname = SetupDiffRhs_loc,
213cb32e2e7SValeria Barra     .applyfname = Diff_loc,
214cb32e2e7SValeria Barra     .errorfname = Error_loc,
215cb32e2e7SValeria Barra     .inmode = CEED_EVAL_GRAD,
216cb32e2e7SValeria Barra     .outmode = CEED_EVAL_GRAD,
217cb32e2e7SValeria Barra     .qmode = CEED_GAUSS
218cb32e2e7SValeria Barra   },
219cb32e2e7SValeria Barra   [CEED_BP4] = {
220cb32e2e7SValeria Barra     .ncompu = 3,
221cb32e2e7SValeria Barra     .qdatasize = 6,
222cb32e2e7SValeria Barra     .qextra = 1,
223cb32e2e7SValeria Barra     .setupgeo = SetupDiffGeo,
224cb32e2e7SValeria Barra     .setuprhs = SetupDiffRhs3,
225cb32e2e7SValeria Barra     .apply = Diff3,
226cb32e2e7SValeria Barra     .error = Error3,
227cb32e2e7SValeria Barra     .setupgeofname = SetupDiffGeo_loc,
228cb32e2e7SValeria Barra     .setuprhsfname = SetupDiffRhs3_loc,
229f9342789SJeremy L Thompson     .applyfname = Diff3_loc,
230cb32e2e7SValeria Barra     .errorfname = Error3_loc,
231cb32e2e7SValeria Barra     .inmode = CEED_EVAL_GRAD,
232cb32e2e7SValeria Barra     .outmode = CEED_EVAL_GRAD,
233cb32e2e7SValeria Barra     .qmode = CEED_GAUSS
234cb32e2e7SValeria Barra   },
235cb32e2e7SValeria Barra   [CEED_BP5] = {
236cb32e2e7SValeria Barra     .ncompu = 1,
237cb32e2e7SValeria Barra     .qdatasize = 6,
238cb32e2e7SValeria Barra     .qextra = 0,
239cb32e2e7SValeria Barra     .setupgeo = SetupDiffGeo,
240cb32e2e7SValeria Barra     .setuprhs = SetupDiffRhs,
241cb32e2e7SValeria Barra     .apply = Diff,
242cb32e2e7SValeria Barra     .error = Error,
243cb32e2e7SValeria Barra     .setupgeofname = SetupDiffGeo_loc,
244cb32e2e7SValeria Barra     .setuprhsfname = SetupDiffRhs_loc,
245cb32e2e7SValeria Barra     .applyfname = Diff_loc,
246cb32e2e7SValeria Barra     .errorfname = Error_loc,
247cb32e2e7SValeria Barra     .inmode = CEED_EVAL_GRAD,
248cb32e2e7SValeria Barra     .outmode = CEED_EVAL_GRAD,
249cb32e2e7SValeria Barra     .qmode = CEED_GAUSS_LOBATTO
250cb32e2e7SValeria Barra   },
251cb32e2e7SValeria Barra   [CEED_BP6] = {
252cb32e2e7SValeria Barra     .ncompu = 3,
253cb32e2e7SValeria Barra     .qdatasize = 6,
254cb32e2e7SValeria Barra     .qextra = 0,
255cb32e2e7SValeria Barra     .setupgeo = SetupDiffGeo,
256cb32e2e7SValeria Barra     .setuprhs = SetupDiffRhs3,
257cb32e2e7SValeria Barra     .apply = Diff3,
258cb32e2e7SValeria Barra     .error = Error3,
259cb32e2e7SValeria Barra     .setupgeofname = SetupDiffGeo_loc,
260cb32e2e7SValeria Barra     .setuprhsfname = SetupDiffRhs3_loc,
261f9342789SJeremy L Thompson     .applyfname = Diff3_loc,
262cb32e2e7SValeria Barra     .errorfname = Error3_loc,
263cb32e2e7SValeria Barra     .inmode = CEED_EVAL_GRAD,
264cb32e2e7SValeria Barra     .outmode = CEED_EVAL_GRAD,
265cb32e2e7SValeria Barra     .qmode = CEED_GAUSS_LOBATTO
266cb32e2e7SValeria Barra   }
267cb32e2e7SValeria Barra };
268cb32e2e7SValeria Barra 
269cb32e2e7SValeria Barra // This function uses libCEED to compute the action of the mass matrix
270cb32e2e7SValeria Barra static PetscErrorCode MatMult_Mass(Mat A, Vec X, Vec Y) {
271cb32e2e7SValeria Barra   PetscErrorCode ierr;
272cb32e2e7SValeria Barra   User user;
273cb32e2e7SValeria Barra   PetscScalar *x, *y;
274*b68a8d79SJed Brown   PetscMemType xmemtype, ymemtype;
275cb32e2e7SValeria Barra 
276cb32e2e7SValeria Barra   PetscFunctionBeginUser;
2779396343dSjeremylt 
278cb32e2e7SValeria Barra   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
2799396343dSjeremylt 
2809396343dSjeremylt   // Global-to-local
281cb32e2e7SValeria Barra   ierr = VecScatterBegin(user->ltog, X, user->Xloc, INSERT_VALUES,
282cb32e2e7SValeria Barra                          SCATTER_REVERSE); CHKERRQ(ierr);
2839396343dSjeremylt   ierr = VecScatterEnd(user->ltog, X, user->Xloc, INSERT_VALUES,
2849396343dSjeremylt                        SCATTER_REVERSE); CHKERRQ(ierr);
285cb32e2e7SValeria Barra 
2869396343dSjeremylt   // Setup libCEED vectors
287*b68a8d79SJed Brown   ierr = VecGetArrayReadAndMemType(user->Xloc, (const PetscScalar **)&x,
288*b68a8d79SJed Brown                                    &xmemtype); CHKERRQ(ierr);
289*b68a8d79SJed Brown   ierr = VecGetArrayAndMemType(user->Yloc, &y, &ymemtype); CHKERRQ(ierr);
290*b68a8d79SJed Brown   CeedVectorSetArray(user->xceed, MemTypeP2C(xmemtype), CEED_USE_POINTER, x);
291*b68a8d79SJed Brown   CeedVectorSetArray(user->yceed, MemTypeP2C(ymemtype), CEED_USE_POINTER, y);
292cb32e2e7SValeria Barra 
2939396343dSjeremylt   // Apply libCEED operator
294cb32e2e7SValeria Barra   CeedOperatorApply(user->op, user->xceed, user->yceed,
295cb32e2e7SValeria Barra                     CEED_REQUEST_IMMEDIATE);
296cb32e2e7SValeria Barra 
2979396343dSjeremylt   // Restore PETSc vectors
298*b68a8d79SJed Brown   CeedVectorTakeArray(user->xceed, MemTypeP2C(xmemtype), NULL);
299*b68a8d79SJed Brown   CeedVectorTakeArray(user->yceed, MemTypeP2C(ymemtype), NULL);
300*b68a8d79SJed Brown   ierr = VecRestoreArrayReadAndMemType(user->Xloc, (const PetscScalar **)&x);
3019396343dSjeremylt   CHKERRQ(ierr);
302*b68a8d79SJed Brown   ierr = VecRestoreArrayAndMemType(user->Yloc, &y); CHKERRQ(ierr);
303cb32e2e7SValeria Barra 
3049396343dSjeremylt   // Local-to-global
305cb32e2e7SValeria Barra   if (Y) {
306cb32e2e7SValeria Barra     ierr = VecZeroEntries(Y); CHKERRQ(ierr);
3079396343dSjeremylt     ierr = VecScatterBegin(user->ltog, user->Yloc, Y, ADD_VALUES,
3089396343dSjeremylt                            SCATTER_FORWARD); CHKERRQ(ierr);
3099396343dSjeremylt     ierr = VecScatterEnd(user->ltog, user->Yloc, Y, ADD_VALUES,
3109396343dSjeremylt                          SCATTER_FORWARD); CHKERRQ(ierr);
311cb32e2e7SValeria Barra   }
312cb32e2e7SValeria Barra   PetscFunctionReturn(0);
313cb32e2e7SValeria Barra }
314cb32e2e7SValeria Barra 
315cb32e2e7SValeria Barra // This function uses libCEED to compute the action of the Laplacian with
316cb32e2e7SValeria Barra // Dirichlet boundary conditions
317cb32e2e7SValeria Barra static PetscErrorCode MatMult_Diff(Mat A, Vec X, Vec Y) {
318cb32e2e7SValeria Barra   PetscErrorCode ierr;
319cb32e2e7SValeria Barra   User user;
320cb32e2e7SValeria Barra   PetscScalar *x, *y;
321*b68a8d79SJed Brown   PetscMemType xmemtype, ymemtype;
322cb32e2e7SValeria Barra 
323cb32e2e7SValeria Barra   PetscFunctionBeginUser;
3249396343dSjeremylt 
325cb32e2e7SValeria Barra   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
326cb32e2e7SValeria Barra 
327cb32e2e7SValeria Barra   // Global-to-local
328cb32e2e7SValeria Barra   ierr = VecScatterBegin(user->ltog0, X, user->Xloc, INSERT_VALUES,
329cb32e2e7SValeria Barra                          SCATTER_REVERSE); CHKERRQ(ierr);
330cb32e2e7SValeria Barra   ierr = VecScatterEnd(user->ltog0, X, user->Xloc, INSERT_VALUES,
331cb32e2e7SValeria Barra                        SCATTER_REVERSE);
332cb32e2e7SValeria Barra   CHKERRQ(ierr);
333cb32e2e7SValeria Barra 
3349396343dSjeremylt   // Setup libCEED vectors
335*b68a8d79SJed Brown   ierr = VecGetArrayReadAndMemType(user->Xloc, (const PetscScalar **)&x,
336*b68a8d79SJed Brown                                    &xmemtype); CHKERRQ(ierr);
337*b68a8d79SJed Brown   ierr = VecGetArrayAndMemType(user->Yloc, &y, &ymemtype); CHKERRQ(ierr);
338*b68a8d79SJed Brown   CeedVectorSetArray(user->xceed, MemTypeP2C(xmemtype), CEED_USE_POINTER, x);
339*b68a8d79SJed Brown   CeedVectorSetArray(user->yceed, MemTypeP2C(ymemtype), CEED_USE_POINTER, y);
340cb32e2e7SValeria Barra 
3419396343dSjeremylt   // Apply libCEED operator
342cb32e2e7SValeria Barra   CeedOperatorApply(user->op, user->xceed, user->yceed,
343cb32e2e7SValeria Barra                     CEED_REQUEST_IMMEDIATE);
344cb32e2e7SValeria Barra 
345cb32e2e7SValeria Barra   // Restore PETSc vectors
346*b68a8d79SJed Brown   CeedVectorTakeArray(user->xceed, MemTypeP2C(xmemtype), NULL);
347*b68a8d79SJed Brown   CeedVectorTakeArray(user->yceed, MemTypeP2C(ymemtype), NULL);
348*b68a8d79SJed Brown   ierr = VecRestoreArrayReadAndMemType(user->Xloc, (const PetscScalar **)&x);
3499396343dSjeremylt   CHKERRQ(ierr);
350*b68a8d79SJed Brown   ierr = VecRestoreArrayAndMemType(user->Yloc, &y); CHKERRQ(ierr);
351cb32e2e7SValeria Barra 
352cb32e2e7SValeria Barra   // Local-to-global
353cb32e2e7SValeria Barra   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
354cb32e2e7SValeria Barra   ierr = VecScatterBegin(user->gtogD, X, Y, INSERT_VALUES, SCATTER_FORWARD);
355cb32e2e7SValeria Barra   CHKERRQ(ierr);
356cb32e2e7SValeria Barra   ierr = VecScatterEnd(user->gtogD, X, Y, INSERT_VALUES, SCATTER_FORWARD);
357cb32e2e7SValeria Barra   CHKERRQ(ierr);
3589396343dSjeremylt   ierr = VecScatterBegin(user->ltog0, user->Yloc, Y, ADD_VALUES,
3599396343dSjeremylt                          SCATTER_FORWARD); CHKERRQ(ierr);
360cb32e2e7SValeria Barra   ierr = VecScatterEnd(user->ltog0, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD);
361cb32e2e7SValeria Barra   CHKERRQ(ierr);
362cb32e2e7SValeria Barra 
363cb32e2e7SValeria Barra   PetscFunctionReturn(0);
364cb32e2e7SValeria Barra }
365cb32e2e7SValeria Barra 
366cb32e2e7SValeria Barra // This function calculates the error in the final solution
367c70bd2a0Sjeremylt static PetscErrorCode ComputeErrorMax(User user, CeedOperator operror, Vec X,
368cb32e2e7SValeria Barra                                       CeedVector target, PetscReal *maxerror) {
369cb32e2e7SValeria Barra   PetscErrorCode ierr;
370cb32e2e7SValeria Barra   PetscScalar *x;
371*b68a8d79SJed Brown   PetscMemType memtype;
372cb32e2e7SValeria Barra   CeedVector collocated_error;
373cb32e2e7SValeria Barra   CeedInt length;
374cb32e2e7SValeria Barra 
375cb32e2e7SValeria Barra   PetscFunctionBeginUser;
3769396343dSjeremylt 
377cb32e2e7SValeria Barra   CeedVectorGetLength(target, &length);
378cb32e2e7SValeria Barra   CeedVectorCreate(user->ceed, length, &collocated_error);
379cb32e2e7SValeria Barra 
380cb32e2e7SValeria Barra   // Global-to-local
381cb32e2e7SValeria Barra   ierr = VecScatterBegin(user->ltog, X, user->Xloc, INSERT_VALUES,
382cb32e2e7SValeria Barra                          SCATTER_REVERSE); CHKERRQ(ierr);
3839396343dSjeremylt   ierr = VecScatterEnd(user->ltog, X, user->Xloc, INSERT_VALUES,
3849396343dSjeremylt                        SCATTER_REVERSE); CHKERRQ(ierr);
3859396343dSjeremylt 
3869396343dSjeremylt   // Setup libCEED vector
387*b68a8d79SJed Brown   ierr = VecGetArrayReadAndMemType(user->Xloc, (const PetscScalar **)&x,
388*b68a8d79SJed Brown                                    &memtype); CHKERRQ(ierr);
389*b68a8d79SJed Brown   CeedVectorSetArray(user->xceed, MemTypeP2C(memtype), CEED_USE_POINTER, x);
390cb32e2e7SValeria Barra 
3919396343dSjeremylt   // Apply libCEED operator
392c70bd2a0Sjeremylt   CeedOperatorApply(operror, user->xceed, collocated_error,
393cb32e2e7SValeria Barra                     CEED_REQUEST_IMMEDIATE);
394cb32e2e7SValeria Barra 
395cb32e2e7SValeria Barra   // Restore PETSc vector
396*b68a8d79SJed Brown   CeedVectorTakeArray(user->xceed, MemTypeP2C(memtype), NULL);
397*b68a8d79SJed Brown   ierr = VecRestoreArrayReadAndMemType(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],
4232fbc6e41SJeremy L Thompson            irank[3], lnodes[3], lsize, ncompu = 1, ksp_max_it_clip[2];
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;
432*b68a8d79SJed Brown   PetscMemType memtype;
433cb32e2e7SValeria Barra   User user;
434cb32e2e7SValeria Barra   Ceed ceed;
435cb32e2e7SValeria Barra   CeedBasis basisx, basisu;
43615910d16Sjeremylt   CeedElemRestriction Erestrictx, Erestrictu, Erestrictui, Erestrictqdi;
437c70bd2a0Sjeremylt   CeedQFunction qfsetupgeo, qfsetuprhs, qfapply, qferror;
438c70bd2a0Sjeremylt   CeedOperator opsetupgeo, opsetuprhs, opapply, operror;
439cb32e2e7SValeria Barra   CeedVector xcoord, qdata, rhsceed, target;
440cb32e2e7SValeria Barra   CeedInt P, Q;
441cb32e2e7SValeria Barra   const CeedInt dim = 3, ncompx = 3;
442199551f5Sjeremylt   bpType bpchoice;
443cb32e2e7SValeria Barra 
444cb32e2e7SValeria Barra   ierr = PetscInitialize(&argc, &argv, NULL, help);
445cb32e2e7SValeria Barra   if (ierr) return ierr;
446cb32e2e7SValeria Barra   comm = PETSC_COMM_WORLD;
44732d2ee49SValeria Barra 
44832d2ee49SValeria Barra   // Read command line options
449cb32e2e7SValeria Barra   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr);
450199551f5Sjeremylt   bpchoice = CEED_BP1;
451cb32e2e7SValeria Barra   ierr = PetscOptionsEnum("-problem",
452cb32e2e7SValeria Barra                           "CEED benchmark problem to solve", NULL,
453199551f5Sjeremylt                           bpTypes, (PetscEnum)bpchoice, (PetscEnum *)&bpchoice,
454cb32e2e7SValeria Barra                           NULL); CHKERRQ(ierr);
455199551f5Sjeremylt   ncompu = bpOptions[bpchoice].ncompu;
456cb32e2e7SValeria Barra   test_mode = PETSC_FALSE;
457cb32e2e7SValeria Barra   ierr = PetscOptionsBool("-test",
458cb32e2e7SValeria Barra                           "Testing mode (do not print unless error is large)",
459cb32e2e7SValeria Barra                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
460cb32e2e7SValeria Barra   benchmark_mode = PETSC_FALSE;
461cb32e2e7SValeria Barra   ierr = PetscOptionsBool("-benchmark",
462cb32e2e7SValeria Barra                           "Benchmarking mode (prints benchmark statistics)",
463cb32e2e7SValeria Barra                           NULL, benchmark_mode, &benchmark_mode, NULL);
464cb32e2e7SValeria Barra   CHKERRQ(ierr);
465cb32e2e7SValeria Barra   write_solution = PETSC_FALSE;
466cb32e2e7SValeria Barra   ierr = PetscOptionsBool("-write_solution",
467cb32e2e7SValeria Barra                           "Write solution for visualization",
468cb32e2e7SValeria Barra                           NULL, write_solution, &write_solution, NULL);
469cb32e2e7SValeria Barra   CHKERRQ(ierr);
470cb32e2e7SValeria Barra   degree = test_mode ? 3 : 1;
471cb32e2e7SValeria Barra   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
472cb32e2e7SValeria Barra                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
473199551f5Sjeremylt   qextra = bpOptions[bpchoice].qextra;
474cb32e2e7SValeria Barra   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
475cb32e2e7SValeria Barra                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
476cb32e2e7SValeria Barra   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
477cb32e2e7SValeria Barra                             NULL, ceedresource, ceedresource,
478cb32e2e7SValeria Barra                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
479cb32e2e7SValeria Barra   localnodes = 1000;
480cb32e2e7SValeria Barra   ierr = PetscOptionsInt("-local",
481cb32e2e7SValeria Barra                          "Target number of locally owned nodes per process",
482cb32e2e7SValeria Barra                          NULL, localnodes, &localnodes, NULL); CHKERRQ(ierr);
4832fbc6e41SJeremy L Thompson   PetscInt two = 2;
4842fbc6e41SJeremy L Thompson   ksp_max_it_clip[0] = 5;
4852fbc6e41SJeremy L Thompson   ksp_max_it_clip[1] = 20;
4862fbc6e41SJeremy L Thompson   ierr = PetscOptionsIntArray("-ksp_max_it_clip",
4872fbc6e41SJeremy L Thompson                               "Min and max number of iterations to use during benchmarking",
4882fbc6e41SJeremy L Thompson                               NULL, ksp_max_it_clip, &two, NULL);
4892fbc6e41SJeremy L Thompson   CHKERRQ(ierr);
490cb32e2e7SValeria Barra   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
491cb32e2e7SValeria Barra   P = degree + 1;
492cb32e2e7SValeria Barra   Q = P + qextra;
493cb32e2e7SValeria Barra 
4949396343dSjeremylt   // Set up libCEED
4959396343dSjeremylt   CeedInit(ceedresource, &ceed);
4969396343dSjeremylt   CeedMemType memtypebackend;
4979396343dSjeremylt   CeedGetPreferredMemType(ceed, &memtypebackend);
4989396343dSjeremylt 
499*b68a8d79SJed Brown   VecType default_vectype = NULL, vectype;
500*b68a8d79SJed Brown   switch (memtypebackend) {
501*b68a8d79SJed Brown   case CEED_MEM_HOST: default_vectype = VECSTANDARD; break;
502*b68a8d79SJed Brown   case CEED_MEM_DEVICE: {
503*b68a8d79SJed Brown     const char *resolved;
504*b68a8d79SJed Brown     CeedGetResource(ceed, &resolved);
505*b68a8d79SJed Brown     if (strstr(resolved, "/gpu/cuda")) default_vectype = VECCUDA;
506*b68a8d79SJed Brown     else if (strstr(resolved, "/gpu/hip/occa"))
507*b68a8d79SJed Brown       default_vectype = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678
508*b68a8d79SJed Brown     else if (strstr(resolved, "/gpu/hip")) default_vectype = VECHIP;
509*b68a8d79SJed Brown     else default_vectype = VECSTANDARD;
510*b68a8d79SJed Brown   }
511*b68a8d79SJed Brown   }
5129396343dSjeremylt 
513cb32e2e7SValeria Barra   // Determine size of process grid
514cb32e2e7SValeria Barra   ierr = MPI_Comm_size(comm, &size); CHKERRQ(ierr);
515cb32e2e7SValeria Barra   Split3(size, p, false);
516cb32e2e7SValeria Barra 
517cb32e2e7SValeria Barra   // Find a nicely composite number of elements no less than localnodes
518cb32e2e7SValeria Barra   for (localelem = PetscMax(1, localnodes / (degree*degree*degree)); ;
519cb32e2e7SValeria Barra        localelem++) {
520cb32e2e7SValeria Barra     Split3(localelem, melem, true);
521cb32e2e7SValeria Barra     if (Max3(melem) / Min3(melem) <= 2) break;
522cb32e2e7SValeria Barra   }
523cb32e2e7SValeria Barra 
524cb32e2e7SValeria Barra   // Find my location in the process grid
525cb32e2e7SValeria Barra   ierr = MPI_Comm_rank(comm, &rank); CHKERRQ(ierr);
526cb32e2e7SValeria Barra   for (int d=0,rankleft=rank; d<dim; d++) {
527cb32e2e7SValeria Barra     const int pstride[3] = {p[1] *p[2], p[2], 1};
528cb32e2e7SValeria Barra     irank[d] = rankleft / pstride[d];
529cb32e2e7SValeria Barra     rankleft -= irank[d] * pstride[d];
530cb32e2e7SValeria Barra   }
531cb32e2e7SValeria Barra 
532cb32e2e7SValeria Barra   GlobalNodes(p, irank, degree, melem, mnodes);
533cb32e2e7SValeria Barra 
534cb32e2e7SValeria Barra   // Setup global vector
5354a5da23aSJed Brown   ierr = VecCreate(comm, &X); CHKERRQ(ierr);
536*b68a8d79SJed Brown   ierr = VecSetType(X, default_vectype); CHKERRQ(ierr);
537cb32e2e7SValeria Barra   ierr = VecSetSizes(X, mnodes[0]*mnodes[1]*mnodes[2]*ncompu, PETSC_DECIDE);
538cb32e2e7SValeria Barra   CHKERRQ(ierr);
539*b68a8d79SJed Brown   ierr = VecSetFromOptions(X); CHKERRQ(ierr);
540cb32e2e7SValeria Barra   ierr = VecSetUp(X); CHKERRQ(ierr);
541cb32e2e7SValeria Barra 
542cb32e2e7SValeria Barra   // Set up libCEED
543cb32e2e7SValeria Barra   CeedInit(ceedresource, &ceed);
544cb32e2e7SValeria Barra 
545cb32e2e7SValeria Barra   // Print summary
546cb32e2e7SValeria Barra   CeedInt gsize;
547cb32e2e7SValeria Barra   ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
548dc7d240cSValeria Barra   if (!test_mode) {
549cb32e2e7SValeria Barra     const char *usedresource;
550cb32e2e7SValeria Barra     CeedGetResource(ceed, &usedresource);
5519396343dSjeremylt 
5529396343dSjeremylt     ierr = VecGetType(X, &vectype); CHKERRQ(ierr);
5539396343dSjeremylt 
554cb32e2e7SValeria Barra     ierr = PetscPrintf(comm,
555cb32e2e7SValeria Barra                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n"
5569396343dSjeremylt                        "  PETSc:\n"
5579396343dSjeremylt                        "    PETSc Vec Type                     : %s\n"
558cb32e2e7SValeria Barra                        "  libCEED:\n"
559cb32e2e7SValeria Barra                        "    libCEED Backend                    : %s\n"
5609396343dSjeremylt                        "    libCEED Backend MemType            : %s\n"
561cb32e2e7SValeria Barra                        "  Mesh:\n"
5628d0bb2bbSvaleriabarra                        "    Number of 1D Basis Nodes (P)       : %d\n"
5638d0bb2bbSvaleriabarra                        "    Number of 1D Quadrature Points (Q) : %d\n"
564cb32e2e7SValeria Barra                        "    Global nodes                       : %D\n"
565cb32e2e7SValeria Barra                        "    Process Decomposition              : %D %D %D\n"
566cb32e2e7SValeria Barra                        "    Local Elements                     : %D = %D %D %D\n"
567db419314Sjeremylt                        "    Owned nodes                        : %D = %D %D %D\n"
568db419314Sjeremylt                        "    DoF per node                       : %D\n",
5699396343dSjeremylt                        bpchoice+1, vectype, usedresource,
5709396343dSjeremylt                        CeedMemTypes[memtypebackend],
5719396343dSjeremylt                        P, Q,  gsize/ncompu, p[0], p[1], p[2], localelem,
5729396343dSjeremylt                        melem[0], melem[1], melem[2],
573cb32e2e7SValeria Barra                        mnodes[0]*mnodes[1]*mnodes[2], mnodes[0], mnodes[1],
574db419314Sjeremylt                        mnodes[2], ncompu); CHKERRQ(ierr);
575cb32e2e7SValeria Barra   }
576cb32e2e7SValeria Barra 
577cb32e2e7SValeria Barra   {
578cb32e2e7SValeria Barra     lsize = 1;
579cb32e2e7SValeria Barra     for (int d=0; d<dim; d++) {
580cb32e2e7SValeria Barra       lnodes[d] = melem[d]*degree + 1;
581cb32e2e7SValeria Barra       lsize *= lnodes[d];
582cb32e2e7SValeria Barra     }
5834a5da23aSJed Brown     ierr = VecCreate(PETSC_COMM_SELF, &Xloc); CHKERRQ(ierr);
584*b68a8d79SJed Brown     ierr = VecSetType(Xloc, default_vectype); CHKERRQ(ierr);
585cb32e2e7SValeria Barra     ierr = VecSetSizes(Xloc, lsize*ncompu, PETSC_DECIDE); CHKERRQ(ierr);
586*b68a8d79SJed Brown     ierr = VecSetFromOptions(Xloc); CHKERRQ(ierr);
587cb32e2e7SValeria Barra     ierr = VecSetUp(Xloc); CHKERRQ(ierr);
588cb32e2e7SValeria Barra 
589cb32e2e7SValeria Barra     // Create local-to-global scatter
590cb32e2e7SValeria Barra     PetscInt *ltogind, *ltogind0, *locind, l0count;
591cb32e2e7SValeria Barra     IS ltogis, ltogis0, locis;
592cb32e2e7SValeria Barra     PetscInt gstart[2][2][2], gmnodes[2][2][2][dim];
593cb32e2e7SValeria Barra 
594cb32e2e7SValeria Barra     for (int i=0; i<2; i++) {
595cb32e2e7SValeria Barra       for (int j=0; j<2; j++) {
596cb32e2e7SValeria Barra         for (int k=0; k<2; k++) {
597cb32e2e7SValeria Barra           PetscInt ijkrank[3] = {irank[0]+i, irank[1]+j, irank[2]+k};
598cb32e2e7SValeria Barra           gstart[i][j][k] = GlobalStart(p, ijkrank, degree, melem);
599cb32e2e7SValeria Barra           GlobalNodes(p, ijkrank, degree, melem, gmnodes[i][j][k]);
600cb32e2e7SValeria Barra         }
601cb32e2e7SValeria Barra       }
602cb32e2e7SValeria Barra     }
603cb32e2e7SValeria Barra 
604cb32e2e7SValeria Barra     ierr = PetscMalloc1(lsize, &ltogind); CHKERRQ(ierr);
605cb32e2e7SValeria Barra     ierr = PetscMalloc1(lsize, &ltogind0); CHKERRQ(ierr);
606cb32e2e7SValeria Barra     ierr = PetscMalloc1(lsize, &locind); CHKERRQ(ierr);
607cb32e2e7SValeria Barra     l0count = 0;
608cb32e2e7SValeria Barra     for (PetscInt i=0,ir,ii; ir=i>=mnodes[0], ii=i-ir*mnodes[0], i<lnodes[0]; i++)
609cb32e2e7SValeria Barra       for (PetscInt j=0,jr,jj; jr=j>=mnodes[1], jj=j-jr*mnodes[1], j<lnodes[1]; j++)
610cb32e2e7SValeria Barra         for (PetscInt k=0,kr,kk; kr=k>=mnodes[2], kk=k-kr*mnodes[2], k<lnodes[2]; k++) {
611cb32e2e7SValeria Barra           PetscInt here = (i*lnodes[1]+j)*lnodes[2]+k;
612cb32e2e7SValeria Barra           ltogind[here] =
613cb32e2e7SValeria Barra             gstart[ir][jr][kr] + (ii*gmnodes[ir][jr][kr][1]+jj)*gmnodes[ir][jr][kr][2]+kk;
614cb32e2e7SValeria Barra           if ((irank[0] == 0 && i == 0)
615cb32e2e7SValeria Barra               || (irank[1] == 0 && j == 0)
616cb32e2e7SValeria Barra               || (irank[2] == 0 && k == 0)
617cb32e2e7SValeria Barra               || (irank[0]+1 == p[0] && i+1 == lnodes[0])
618cb32e2e7SValeria Barra               || (irank[1]+1 == p[1] && j+1 == lnodes[1])
619cb32e2e7SValeria Barra               || (irank[2]+1 == p[2] && k+1 == lnodes[2]))
620cb32e2e7SValeria Barra             continue;
621cb32e2e7SValeria Barra           ltogind0[l0count] = ltogind[here];
622cb32e2e7SValeria Barra           locind[l0count++] = here;
623cb32e2e7SValeria Barra         }
624cb32e2e7SValeria Barra     ierr = ISCreateBlock(comm, ncompu, lsize, ltogind, PETSC_OWN_POINTER,
625cb32e2e7SValeria Barra                          &ltogis); CHKERRQ(ierr);
626cb32e2e7SValeria Barra     ierr = VecScatterCreate(Xloc, NULL, X, ltogis, &ltog); CHKERRQ(ierr);
627cb32e2e7SValeria Barra     CHKERRQ(ierr);
628cb32e2e7SValeria Barra     ierr = ISCreateBlock(comm, ncompu, l0count, ltogind0, PETSC_OWN_POINTER,
629cb32e2e7SValeria Barra                          &ltogis0); CHKERRQ(ierr);
630cb32e2e7SValeria Barra     ierr = ISCreateBlock(comm, ncompu, l0count, locind, PETSC_OWN_POINTER,
631cb32e2e7SValeria Barra                          &locis); CHKERRQ(ierr);
632cb32e2e7SValeria Barra     ierr = VecScatterCreate(Xloc, locis, X, ltogis0, &ltog0); CHKERRQ(ierr);
633cb32e2e7SValeria Barra     {
634cb32e2e7SValeria Barra       // Create global-to-global scatter for Dirichlet values (everything not in
635cb32e2e7SValeria Barra       // ltogis0, which is the range of ltog0)
636cb32e2e7SValeria Barra       PetscInt xstart, xend, *indD, countD = 0;
637cb32e2e7SValeria Barra       IS isD;
638cb32e2e7SValeria Barra       const PetscScalar *x;
639cb32e2e7SValeria Barra       ierr = VecZeroEntries(Xloc); CHKERRQ(ierr);
640cb32e2e7SValeria Barra       ierr = VecSet(X, 1.0); CHKERRQ(ierr);
641cb32e2e7SValeria Barra       ierr = VecScatterBegin(ltog0, Xloc, X, INSERT_VALUES, SCATTER_FORWARD);
642cb32e2e7SValeria Barra       CHKERRQ(ierr);
643cb32e2e7SValeria Barra       ierr = VecScatterEnd(ltog0, Xloc, X, INSERT_VALUES, SCATTER_FORWARD);
644cb32e2e7SValeria Barra       CHKERRQ(ierr);
645cb32e2e7SValeria Barra       ierr = VecGetOwnershipRange(X, &xstart, &xend); CHKERRQ(ierr);
646cb32e2e7SValeria Barra       ierr = PetscMalloc1(xend-xstart, &indD); CHKERRQ(ierr);
647cb32e2e7SValeria Barra       ierr = VecGetArrayRead(X, &x); CHKERRQ(ierr);
648cb32e2e7SValeria Barra       for (PetscInt i=0; i<xend-xstart; i++) {
649cb32e2e7SValeria Barra         if (x[i] == 1.) indD[countD++] = xstart + i;
650cb32e2e7SValeria Barra       }
651cb32e2e7SValeria Barra       ierr = VecRestoreArrayRead(X, &x); CHKERRQ(ierr);
652cb32e2e7SValeria Barra       ierr = ISCreateGeneral(comm, countD, indD, PETSC_COPY_VALUES, &isD);
653cb32e2e7SValeria Barra       CHKERRQ(ierr);
654cb32e2e7SValeria Barra       ierr = PetscFree(indD); CHKERRQ(ierr);
655cb32e2e7SValeria Barra       ierr = VecScatterCreate(X, isD, X, isD, &gtogD); CHKERRQ(ierr);
656cb32e2e7SValeria Barra       ierr = ISDestroy(&isD); CHKERRQ(ierr);
657cb32e2e7SValeria Barra     }
658cb32e2e7SValeria Barra     ierr = ISDestroy(&ltogis); CHKERRQ(ierr);
659cb32e2e7SValeria Barra     ierr = ISDestroy(&ltogis0); CHKERRQ(ierr);
660cb32e2e7SValeria Barra     ierr = ISDestroy(&locis); CHKERRQ(ierr);
661cb32e2e7SValeria Barra   }
662cb32e2e7SValeria Barra 
663cb32e2e7SValeria Barra   // CEED bases
664cb32e2e7SValeria Barra   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompu, P, Q,
665199551f5Sjeremylt                                   bpOptions[bpchoice].qmode, &basisu);
666cb32e2e7SValeria Barra   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, Q,
667199551f5Sjeremylt                                   bpOptions[bpchoice].qmode, &basisx);
668cb32e2e7SValeria Barra 
669cb32e2e7SValeria Barra   // CEED restrictions
670d979a051Sjeremylt   CreateRestriction(ceed, melem, P, ncompu, &Erestrictu);
671d979a051Sjeremylt   CreateRestriction(ceed, melem, 2, dim, &Erestrictx);
672cb32e2e7SValeria Barra   CeedInt nelem = melem[0]*melem[1]*melem[2];
673d979a051Sjeremylt   CeedElemRestrictionCreateStrided(ceed, nelem, Q*Q*Q, ncompu,
674d979a051Sjeremylt                                    ncompu*nelem*Q*Q*Q,
675523b8ea0Sjeremylt                                    CEED_STRIDES_BACKEND, &Erestrictui);
676d979a051Sjeremylt   CeedElemRestrictionCreateStrided(ceed, nelem, Q*Q*Q,
677199551f5Sjeremylt                                    bpOptions[bpchoice].qdatasize,
678d979a051Sjeremylt                                    bpOptions[bpchoice].qdatasize*nelem*Q*Q*Q,
679523b8ea0Sjeremylt                                    CEED_STRIDES_BACKEND, &Erestrictqdi);
680cb32e2e7SValeria Barra   {
681cb32e2e7SValeria Barra     CeedScalar *xloc;
682cb32e2e7SValeria Barra     CeedInt shape[3] = {melem[0]+1, melem[1]+1, melem[2]+1}, len =
683cb32e2e7SValeria Barra                          shape[0]*shape[1]*shape[2];
684cb32e2e7SValeria Barra     xloc = malloc(len*ncompx*sizeof xloc[0]);
685cb32e2e7SValeria Barra     for (CeedInt i=0; i<shape[0]; i++) {
686cb32e2e7SValeria Barra       for (CeedInt j=0; j<shape[1]; j++) {
687cb32e2e7SValeria Barra         for (CeedInt k=0; k<shape[2]; k++) {
688d979a051Sjeremylt           xloc[dim*((i*shape[1]+j)*shape[2]+k) + 0] = 1.*(irank[0]*melem[0]+i) /
689cb32e2e7SValeria Barra               (p[0]*melem[0]);
690d979a051Sjeremylt           xloc[dim*((i*shape[1]+j)*shape[2]+k) + 1] = 1.*(irank[1]*melem[1]+j) /
691cb32e2e7SValeria Barra               (p[1]*melem[1]);
692d979a051Sjeremylt           xloc[dim*((i*shape[1]+j)*shape[2]+k) + 2] = 1.*(irank[2]*melem[2]+k) /
693cb32e2e7SValeria Barra               (p[2]*melem[2]);
694cb32e2e7SValeria Barra         }
695cb32e2e7SValeria Barra       }
696cb32e2e7SValeria Barra     }
697cb32e2e7SValeria Barra     CeedVectorCreate(ceed, len*ncompx, &xcoord);
698cb32e2e7SValeria Barra     CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_OWN_POINTER, xloc);
699cb32e2e7SValeria Barra   }
700cb32e2e7SValeria Barra 
701cb32e2e7SValeria Barra   // Create the Qfunction that builds the operator quadrature data
702199551f5Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].setupgeo,
703199551f5Sjeremylt                               bpOptions[bpchoice].setupgeofname, &qfsetupgeo);
704c70bd2a0Sjeremylt   CeedQFunctionAddInput(qfsetupgeo, "dx", ncompx*dim, CEED_EVAL_GRAD);
705c70bd2a0Sjeremylt   CeedQFunctionAddInput(qfsetupgeo, "weight", 1, CEED_EVAL_WEIGHT);
706199551f5Sjeremylt   CeedQFunctionAddOutput(qfsetupgeo, "qdata", bpOptions[bpchoice].qdatasize,
707cb32e2e7SValeria Barra                          CEED_EVAL_NONE);
708cb32e2e7SValeria Barra 
709cb32e2e7SValeria Barra   // Create the Qfunction that sets up the RHS and true solution
710199551f5Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].setuprhs,
711199551f5Sjeremylt                               bpOptions[bpchoice].setuprhsfname, &qfsetuprhs);
712c70bd2a0Sjeremylt   CeedQFunctionAddInput(qfsetuprhs, "x", ncompx, CEED_EVAL_INTERP);
713c70bd2a0Sjeremylt   CeedQFunctionAddInput(qfsetuprhs, "dx", ncompx*dim, CEED_EVAL_GRAD);
714c70bd2a0Sjeremylt   CeedQFunctionAddInput(qfsetuprhs, "weight", 1, CEED_EVAL_WEIGHT);
715c70bd2a0Sjeremylt   CeedQFunctionAddOutput(qfsetuprhs, "true_soln", ncompu, CEED_EVAL_NONE);
716c70bd2a0Sjeremylt   CeedQFunctionAddOutput(qfsetuprhs, "rhs", ncompu, CEED_EVAL_INTERP);
717cb32e2e7SValeria Barra 
718cb32e2e7SValeria Barra   // Set up PDE operator
719199551f5Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].apply,
720199551f5Sjeremylt                               bpOptions[bpchoice].applyfname, &qfapply);
721cb32e2e7SValeria Barra   // Add inputs and outputs
722199551f5Sjeremylt   CeedInt inscale = bpOptions[bpchoice].inmode==CEED_EVAL_GRAD ? 3 : 1;
723199551f5Sjeremylt   CeedInt outscale = bpOptions[bpchoice].outmode==CEED_EVAL_GRAD ? 3 : 1;
724c70bd2a0Sjeremylt   CeedQFunctionAddInput(qfapply, "u", ncompu*inscale,
725199551f5Sjeremylt                         bpOptions[bpchoice].inmode);
726199551f5Sjeremylt   CeedQFunctionAddInput(qfapply, "qdata", bpOptions[bpchoice].qdatasize,
727cb32e2e7SValeria Barra                         CEED_EVAL_NONE);
728c70bd2a0Sjeremylt   CeedQFunctionAddOutput(qfapply, "v", ncompu*outscale,
729199551f5Sjeremylt                          bpOptions[bpchoice].outmode);
730cb32e2e7SValeria Barra 
731cb32e2e7SValeria Barra   // Create the error qfunction
732199551f5Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].error,
733199551f5Sjeremylt                               bpOptions[bpchoice].errorfname, &qferror);
734c70bd2a0Sjeremylt   CeedQFunctionAddInput(qferror, "u", ncompu, CEED_EVAL_INTERP);
735c70bd2a0Sjeremylt   CeedQFunctionAddInput(qferror, "true_soln", ncompu, CEED_EVAL_NONE);
736c70bd2a0Sjeremylt   CeedQFunctionAddOutput(qferror, "error", ncompu, CEED_EVAL_NONE);
737cb32e2e7SValeria Barra 
738cb32e2e7SValeria Barra   // Create the persistent vectors that will be needed in setup
739cb32e2e7SValeria Barra   CeedInt nqpts;
740cb32e2e7SValeria Barra   CeedBasisGetNumQuadraturePoints(basisu, &nqpts);
741199551f5Sjeremylt   CeedVectorCreate(ceed, bpOptions[bpchoice].qdatasize*nelem*nqpts, &qdata);
742cb32e2e7SValeria Barra   CeedVectorCreate(ceed, nelem*nqpts*ncompu, &target);
743cb32e2e7SValeria Barra   CeedVectorCreate(ceed, lsize*ncompu, &rhsceed);
744cb32e2e7SValeria Barra 
745cb32e2e7SValeria Barra   // Create the operator that builds the quadrature data for the ceed operator
746c70bd2a0Sjeremylt   CeedOperatorCreate(ceed, qfsetupgeo, CEED_QFUNCTION_NONE,
747c70bd2a0Sjeremylt                      CEED_QFUNCTION_NONE, &opsetupgeo);
748c70bd2a0Sjeremylt   CeedOperatorSetField(opsetupgeo, "dx", Erestrictx, basisx,
749a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
750c70bd2a0Sjeremylt   CeedOperatorSetField(opsetupgeo, "weight", CEED_ELEMRESTRICTION_NONE, basisx,
751a8d32208Sjeremylt                        CEED_VECTOR_NONE);
752c70bd2a0Sjeremylt   CeedOperatorSetField(opsetupgeo, "qdata", Erestrictqdi,
753cb32e2e7SValeria Barra                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
754cb32e2e7SValeria Barra 
755cb32e2e7SValeria Barra   // Create the operator that builds the RHS and true solution
756c70bd2a0Sjeremylt   CeedOperatorCreate(ceed, qfsetuprhs, CEED_QFUNCTION_NONE,
757c70bd2a0Sjeremylt                      CEED_QFUNCTION_NONE, &opsetuprhs);
758c70bd2a0Sjeremylt   CeedOperatorSetField(opsetuprhs, "x", Erestrictx, basisx,
759a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
760c70bd2a0Sjeremylt   CeedOperatorSetField(opsetuprhs, "dx", Erestrictx, basisx,
761a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
762c70bd2a0Sjeremylt   CeedOperatorSetField(opsetuprhs, "weight", CEED_ELEMRESTRICTION_NONE, basisx,
763a8d32208Sjeremylt                        CEED_VECTOR_NONE);
764c70bd2a0Sjeremylt   CeedOperatorSetField(opsetuprhs, "true_soln", Erestrictui,
765cb32e2e7SValeria Barra                        CEED_BASIS_COLLOCATED, target);
766c70bd2a0Sjeremylt   CeedOperatorSetField(opsetuprhs, "rhs", Erestrictu, basisu,
767a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
768cb32e2e7SValeria Barra 
769cb32e2e7SValeria Barra   // Create the mass or diff operator
770c70bd2a0Sjeremylt   CeedOperatorCreate(ceed, qfapply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
771c70bd2a0Sjeremylt                      &opapply);
772c70bd2a0Sjeremylt   CeedOperatorSetField(opapply, "u", Erestrictu, basisu, CEED_VECTOR_ACTIVE);
773c70bd2a0Sjeremylt   CeedOperatorSetField(opapply, "qdata", Erestrictqdi, CEED_BASIS_COLLOCATED,
774a8d32208Sjeremylt                        qdata);
775c70bd2a0Sjeremylt   CeedOperatorSetField(opapply, "v", Erestrictu, basisu, CEED_VECTOR_ACTIVE);
776cb32e2e7SValeria Barra 
777cb32e2e7SValeria Barra   // Create the error operator
778c70bd2a0Sjeremylt   CeedOperatorCreate(ceed, qferror, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
779c70bd2a0Sjeremylt                      &operror);
780c70bd2a0Sjeremylt   CeedOperatorSetField(operror, "u", Erestrictu, basisu, CEED_VECTOR_ACTIVE);
781c70bd2a0Sjeremylt   CeedOperatorSetField(operror, "true_soln", Erestrictui,
782cb32e2e7SValeria Barra                        CEED_BASIS_COLLOCATED, target);
783c70bd2a0Sjeremylt   CeedOperatorSetField(operror, "error", Erestrictui, CEED_BASIS_COLLOCATED,
784a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
785cb32e2e7SValeria Barra 
786cb32e2e7SValeria Barra   // Set up Mat
787cb32e2e7SValeria Barra   ierr = PetscMalloc1(1, &user); CHKERRQ(ierr);
788cb32e2e7SValeria Barra   user->comm = comm;
789cb32e2e7SValeria Barra   user->ltog = ltog;
790199551f5Sjeremylt   if (bpchoice != CEED_BP1 && bpchoice != CEED_BP2) {
791cb32e2e7SValeria Barra     user->ltog0 = ltog0;
792cb32e2e7SValeria Barra     user->gtogD = gtogD;
793cb32e2e7SValeria Barra   }
794cb32e2e7SValeria Barra   user->Xloc = Xloc;
795cb32e2e7SValeria Barra   ierr = VecDuplicate(Xloc, &user->Yloc); CHKERRQ(ierr);
796cb32e2e7SValeria Barra   CeedVectorCreate(ceed, lsize*ncompu, &user->xceed);
797cb32e2e7SValeria Barra   CeedVectorCreate(ceed, lsize*ncompu, &user->yceed);
798c70bd2a0Sjeremylt   user->op = opapply;
799cb32e2e7SValeria Barra   user->qdata = qdata;
800cb32e2e7SValeria Barra   user->ceed = ceed;
801cb32e2e7SValeria Barra 
802cb32e2e7SValeria Barra   ierr = MatCreateShell(comm, mnodes[0]*mnodes[1]*mnodes[2]*ncompu,
803cb32e2e7SValeria Barra                         mnodes[0]*mnodes[1]*mnodes[2]*ncompu,
804cb32e2e7SValeria Barra                         PETSC_DECIDE, PETSC_DECIDE, user, &mat); CHKERRQ(ierr);
805199551f5Sjeremylt   if (bpchoice == CEED_BP1 || bpchoice == CEED_BP2) {
806cb32e2e7SValeria Barra     ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Mass);
807cb32e2e7SValeria Barra     CHKERRQ(ierr);
808cb32e2e7SValeria Barra   } else {
809cb32e2e7SValeria Barra     ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Diff);
810cb32e2e7SValeria Barra     CHKERRQ(ierr);
811cb32e2e7SValeria Barra   }
812*b68a8d79SJed Brown   ierr = VecGetType(X, &vectype); CHKERRQ(ierr);
813*b68a8d79SJed Brown   ierr = MatShellSetVecType(mat, vectype); CHKERRQ(ierr);
814cb32e2e7SValeria Barra 
815cb32e2e7SValeria Barra   // Get RHS vector
8169396343dSjeremylt   ierr = VecDuplicate(X, &rhs); CHKERRQ(ierr);
817cb32e2e7SValeria Barra   ierr = VecDuplicate(Xloc, &rhsloc); CHKERRQ(ierr);
818cb32e2e7SValeria Barra   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
819*b68a8d79SJed Brown   ierr = VecGetArrayAndMemType(rhsloc, &r, &memtype); CHKERRQ(ierr);
820*b68a8d79SJed Brown   CeedVectorSetArray(rhsceed, MemTypeP2C(memtype), CEED_USE_POINTER, r);
821cb32e2e7SValeria Barra 
822cb32e2e7SValeria Barra   // Setup qdata, rhs, and target
823c70bd2a0Sjeremylt   CeedOperatorApply(opsetupgeo, xcoord, qdata, CEED_REQUEST_IMMEDIATE);
824c70bd2a0Sjeremylt   CeedOperatorApply(opsetuprhs, xcoord, rhsceed, CEED_REQUEST_IMMEDIATE);
825cb32e2e7SValeria Barra   CeedVectorDestroy(&xcoord);
826cb32e2e7SValeria Barra 
827cb32e2e7SValeria Barra   // Gather RHS
828*b68a8d79SJed Brown   ierr = CeedVectorTakeArray(rhsceed, MemTypeP2C(memtype), NULL); CHKERRQ(ierr);
829*b68a8d79SJed Brown   ierr = VecRestoreArrayAndMemType(rhsloc, &r); CHKERRQ(ierr);
830cb32e2e7SValeria Barra   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
831cb32e2e7SValeria Barra   ierr = VecScatterBegin(ltog, rhsloc, rhs, ADD_VALUES, SCATTER_FORWARD);
832cb32e2e7SValeria Barra   CHKERRQ(ierr);
833cb32e2e7SValeria Barra   ierr = VecScatterEnd(ltog, rhsloc, rhs, ADD_VALUES, SCATTER_FORWARD);
834cb32e2e7SValeria Barra   CHKERRQ(ierr);
835cb32e2e7SValeria Barra   CeedVectorDestroy(&rhsceed);
836cb32e2e7SValeria Barra 
837cb32e2e7SValeria Barra   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
838cb32e2e7SValeria Barra   {
839cb32e2e7SValeria Barra     PC pc;
840cb32e2e7SValeria Barra     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
841199551f5Sjeremylt     if (bpchoice == CEED_BP1 || bpchoice == CEED_BP2) {
842cb32e2e7SValeria Barra       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
843cb32e2e7SValeria Barra       ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr);
844cb32e2e7SValeria Barra     } else {
845cb32e2e7SValeria Barra       ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr);
846cb32e2e7SValeria Barra     }
847cb32e2e7SValeria Barra     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
848cb32e2e7SValeria Barra     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
849cb32e2e7SValeria Barra     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
850cb32e2e7SValeria Barra                             PETSC_DEFAULT); CHKERRQ(ierr);
851cb32e2e7SValeria Barra   }
852cb32e2e7SValeria Barra   ierr = KSPSetOperators(ksp, mat, mat); CHKERRQ(ierr);
853da9108adSvaleriabarra   // First run's performance log is not considered for benchmarking purposes
854cb32e2e7SValeria Barra   ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
855cb32e2e7SValeria Barra   CHKERRQ(ierr);
856cb32e2e7SValeria Barra   my_rt_start = MPI_Wtime();
857cb32e2e7SValeria Barra   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
858cb32e2e7SValeria Barra   my_rt = MPI_Wtime() - my_rt_start;
859cb32e2e7SValeria Barra   ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
860cb32e2e7SValeria Barra   CHKERRQ(ierr);
861cb32e2e7SValeria Barra   // Set maxits based on first iteration timing
862cb32e2e7SValeria Barra   if (my_rt > 0.02) {
8632fbc6e41SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
8642fbc6e41SJeremy L Thompson                             ksp_max_it_clip[0]); CHKERRQ(ierr);
865cb32e2e7SValeria Barra   } else {
8662fbc6e41SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
8672fbc6e41SJeremy L Thompson                             ksp_max_it_clip[1]); CHKERRQ(ierr);
868cb32e2e7SValeria Barra   }
869debcf919SJed Brown   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
87009a940d7Sjeremylt 
871cb32e2e7SValeria Barra   // Timed solve
87209a940d7Sjeremylt   ierr = VecZeroEntries(X); CHKERRQ(ierr);
873cb32e2e7SValeria Barra   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
87409a940d7Sjeremylt 
87509a940d7Sjeremylt   // -- Performance logging
87609a940d7Sjeremylt   ierr = PetscLogStageRegister("Solve Stage", &solvestage); CHKERRQ(ierr);
87709a940d7Sjeremylt   ierr = PetscLogStagePush(solvestage); CHKERRQ(ierr);
87809a940d7Sjeremylt 
87909a940d7Sjeremylt   // -- Solve
880cb32e2e7SValeria Barra   my_rt_start = MPI_Wtime();
881cb32e2e7SValeria Barra   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
882cb32e2e7SValeria Barra   my_rt = MPI_Wtime() - my_rt_start;
88309a940d7Sjeremylt 
88409a940d7Sjeremylt   // -- Performance logging
88509a940d7Sjeremylt   ierr = PetscLogStagePop();
88609a940d7Sjeremylt 
8878e87e98bSJed Brown   // Output results
888cb32e2e7SValeria Barra   {
889cb32e2e7SValeria Barra     KSPType ksptype;
890cb32e2e7SValeria Barra     KSPConvergedReason reason;
891cb32e2e7SValeria Barra     PetscReal rnorm;
892cb32e2e7SValeria Barra     PetscInt its;
893cb32e2e7SValeria Barra     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
894cb32e2e7SValeria Barra     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
895cb32e2e7SValeria Barra     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
896cb32e2e7SValeria Barra     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
8970ec79729SJeremy L Thompson     if (!test_mode || reason < 0 || rnorm > 1e-9) {
898cb32e2e7SValeria Barra       ierr = PetscPrintf(comm,
899cb32e2e7SValeria Barra                          "  KSP:\n"
900cb32e2e7SValeria Barra                          "    KSP Type                           : %s\n"
901cb32e2e7SValeria Barra                          "    KSP Convergence                    : %s\n"
902cb32e2e7SValeria Barra                          "    Total KSP Iterations               : %D\n"
903cb32e2e7SValeria Barra                          "    Final rnorm                        : %e\n",
904cb32e2e7SValeria Barra                          ksptype, KSPConvergedReasons[reason], its,
905cb32e2e7SValeria Barra                          (double)rnorm); CHKERRQ(ierr);
906cb32e2e7SValeria Barra     }
9078e87e98bSJed Brown     if (!test_mode) {
9088e87e98bSJed Brown       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
9098e87e98bSJed Brown     }
9108e87e98bSJed Brown     {
9118e87e98bSJed Brown       PetscReal maxerror;
9128e87e98bSJed Brown       ierr = ComputeErrorMax(user, operror, X, target, &maxerror);
9138e87e98bSJed Brown       CHKERRQ(ierr);
9148e87e98bSJed Brown       PetscReal tol = 5e-2;
9158e87e98bSJed Brown       if (!test_mode || maxerror > tol) {
916cb32e2e7SValeria Barra         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
917cb32e2e7SValeria Barra         CHKERRQ(ierr);
918cb32e2e7SValeria Barra         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
919cb32e2e7SValeria Barra         CHKERRQ(ierr);
920cb32e2e7SValeria Barra         ierr = PetscPrintf(comm,
9218e87e98bSJed Brown                            "    Pointwise Error (max)              : %e\n"
9228e87e98bSJed Brown                            "    CG Solve Time                      : %g (%g) sec\n",
9238e87e98bSJed Brown                            (double)maxerror, rt_max, rt_min); CHKERRQ(ierr);
924cb32e2e7SValeria Barra       }
925cb32e2e7SValeria Barra     }
9268d0bb2bbSvaleriabarra     if (!test_mode) {
927cb32e2e7SValeria Barra       ierr = PetscPrintf(comm,
9288e87e98bSJed Brown                          "    DoFs/Sec in CG                     : %g (%g) million\n",
9298e87e98bSJed Brown                          1e-6*gsize*its/rt_max,
9308e87e98bSJed Brown                          1e-6*gsize*its/rt_min); CHKERRQ(ierr);
931cb32e2e7SValeria Barra     }
932cb32e2e7SValeria Barra   }
933cb32e2e7SValeria Barra 
934cb32e2e7SValeria Barra   if (write_solution) {
935cb32e2e7SValeria Barra     PetscViewer vtkviewersoln;
936cb32e2e7SValeria Barra 
937cb32e2e7SValeria Barra     ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr);
938cb32e2e7SValeria Barra     ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr);
939cb32e2e7SValeria Barra     ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr);
940cb32e2e7SValeria Barra     ierr = VecView(X, vtkviewersoln); CHKERRQ(ierr);
941cb32e2e7SValeria Barra     ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr);
942cb32e2e7SValeria Barra   }
943cb32e2e7SValeria Barra 
944cb32e2e7SValeria Barra   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
945cb32e2e7SValeria Barra   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
946cb32e2e7SValeria Barra   ierr = VecDestroy(&X); CHKERRQ(ierr);
947cb32e2e7SValeria Barra   ierr = VecDestroy(&user->Xloc); CHKERRQ(ierr);
948cb32e2e7SValeria Barra   ierr = VecDestroy(&user->Yloc); CHKERRQ(ierr);
949cb32e2e7SValeria Barra   ierr = VecScatterDestroy(&ltog); CHKERRQ(ierr);
950cb32e2e7SValeria Barra   ierr = VecScatterDestroy(&ltog0); CHKERRQ(ierr);
951cb32e2e7SValeria Barra   ierr = VecScatterDestroy(&gtogD); CHKERRQ(ierr);
952cb32e2e7SValeria Barra   ierr = MatDestroy(&mat); CHKERRQ(ierr);
953cb32e2e7SValeria Barra   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
954cb32e2e7SValeria Barra 
955cb32e2e7SValeria Barra   CeedVectorDestroy(&user->xceed);
956cb32e2e7SValeria Barra   CeedVectorDestroy(&user->yceed);
957cb32e2e7SValeria Barra   CeedVectorDestroy(&user->qdata);
958cb32e2e7SValeria Barra   CeedVectorDestroy(&target);
959c70bd2a0Sjeremylt   CeedOperatorDestroy(&opsetupgeo);
960c70bd2a0Sjeremylt   CeedOperatorDestroy(&opsetuprhs);
961c70bd2a0Sjeremylt   CeedOperatorDestroy(&opapply);
962c70bd2a0Sjeremylt   CeedOperatorDestroy(&operror);
963cb32e2e7SValeria Barra   CeedElemRestrictionDestroy(&Erestrictu);
964cb32e2e7SValeria Barra   CeedElemRestrictionDestroy(&Erestrictx);
965cb32e2e7SValeria Barra   CeedElemRestrictionDestroy(&Erestrictui);
966cb32e2e7SValeria Barra   CeedElemRestrictionDestroy(&Erestrictqdi);
967c70bd2a0Sjeremylt   CeedQFunctionDestroy(&qfsetupgeo);
968c70bd2a0Sjeremylt   CeedQFunctionDestroy(&qfsetuprhs);
969c70bd2a0Sjeremylt   CeedQFunctionDestroy(&qfapply);
970c70bd2a0Sjeremylt   CeedQFunctionDestroy(&qferror);
971cb32e2e7SValeria Barra   CeedBasisDestroy(&basisu);
972cb32e2e7SValeria Barra   CeedBasisDestroy(&basisx);
973cb32e2e7SValeria Barra   CeedDestroy(&ceed);
974cb32e2e7SValeria Barra   ierr = PetscFree(user); CHKERRQ(ierr);
975cb32e2e7SValeria Barra   return PetscFinalize();
976cb32e2e7SValeria Barra }
977