xref: /libCEED/examples/petsc/bpsraw.c (revision 9b072555b57804a6f4e0fc2b1ad83be89838f0e5)
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 //
38*9b072555Sjeremylt //TESTARGS -ceed {ceed_resource} -test -problem bp2 -degree 5 -q_extra 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 
453d576824SJeremy L Thompson #include <ceed.h>
463d576824SJeremy L Thompson #include <petscksp.h>
473d576824SJeremy L Thompson #include <petscsys.h>
48cb32e2e7SValeria Barra #include <stdbool.h>
49cb32e2e7SValeria Barra #include <string.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"
543d576824SJeremy L Thompson #include "qfunctions/bps/common.h"
55cb32e2e7SValeria Barra 
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*9b072555Sjeremylt static CeedMemType MemTypeP2C(PetscMemType mem_type) {
65*9b072555Sjeremylt   return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST;
66b68a8d79SJed Brown }
67b68a8d79SJed Brown 
68cb32e2e7SValeria Barra static void Split3(PetscInt size, PetscInt m[3], bool reverse) {
69*9b072555Sjeremylt   for (PetscInt d=0, size_left=size; d<3; d++) {
70*9b072555Sjeremylt     PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(size_left, 1./(3 - d)));
71*9b072555Sjeremylt     while (try * (size_left / try) != size_left) try++;
72cb32e2e7SValeria Barra     m[reverse ? 2-d : d] = try;
73*9b072555Sjeremylt     size_left /= 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 }
83*9b072555Sjeremylt static void GlobalNodes(const PetscInt p[3], const PetscInt i_rank[3],
84*9b072555Sjeremylt                         PetscInt degree, const PetscInt mesh_elem[3],
85*9b072555Sjeremylt                         PetscInt m_nodes[3]) {
86cb32e2e7SValeria Barra   for (int d=0; d<3; d++)
87*9b072555Sjeremylt     m_nodes[d] = degree*mesh_elem[d] + (i_rank[d] == p[d]-1);
88cb32e2e7SValeria Barra }
89*9b072555Sjeremylt static PetscInt GlobalStart(const PetscInt p[3], const PetscInt i_rank[3],
90*9b072555Sjeremylt                             PetscInt degree, const PetscInt mesh_elem[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++) {
96*9b072555Sjeremylt         PetscInt m_nodes[3], ijk_rank[] = {i,j,k};
97*9b072555Sjeremylt         if (i == i_rank[0] && j == i_rank[1] && k == i_rank[2]) return start;
98*9b072555Sjeremylt         GlobalNodes(p, ijk_rank, degree, mesh_elem, m_nodes);
99*9b072555Sjeremylt         start += m_nodes[0] * m_nodes[1] * m_nodes[2];
100cb32e2e7SValeria Barra       }
101cb32e2e7SValeria Barra     }
102cb32e2e7SValeria Barra   }
103cb32e2e7SValeria Barra   return -1;
104cb32e2e7SValeria Barra }
105*9b072555Sjeremylt static int CreateRestriction(Ceed ceed, const CeedInt mesh_elem[3], CeedInt P,
106*9b072555Sjeremylt                              CeedInt num_comp, CeedElemRestriction *elem_restr) {
107*9b072555Sjeremylt   const PetscInt num_elem = mesh_elem[0]*mesh_elem[1]*mesh_elem[2];
108*9b072555Sjeremylt   PetscInt m_nodes[3], *idx, *idx_p;
109cb32e2e7SValeria Barra 
110cb32e2e7SValeria Barra   // Get indicies
111*9b072555Sjeremylt   for (int d=0; d<3; d++) m_nodes[d] = mesh_elem[d]*(P-1) + 1;
112*9b072555Sjeremylt   idx_p = idx = malloc(num_elem*P*P*P*sizeof idx[0]);
113*9b072555Sjeremylt   for (CeedInt i=0; i<mesh_elem[0]; i++)
114*9b072555Sjeremylt     for (CeedInt j=0; j<mesh_elem[1]; j++)
115*9b072555Sjeremylt       for (CeedInt k=0; k<mesh_elem[2]; k++,idx_p += 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
120*9b072555Sjeremylt                 idx_p[(ii*P+jj)*P+kk] = num_comp*(((i*(P-1)+ii)*m_nodes[1]
121*9b072555Sjeremylt                                                    + (j*(P-1)+jj))*m_nodes[2]
122cb32e2e7SValeria Barra                                                   + (k*(P-1)+kk));
123cb32e2e7SValeria Barra               } else { // (k,j,i) ordering for consistency with MFEM example
124*9b072555Sjeremylt                 idx_p[ii+P*(jj+P*kk)] = num_comp*(((i*(P-1)+ii)*m_nodes[1]
125*9b072555Sjeremylt                                                    + (j*(P-1)+jj))*m_nodes[2]
126cb32e2e7SValeria Barra                                                   + (k*(P-1)+kk));
127cb32e2e7SValeria Barra               }
128cb32e2e7SValeria Barra             }
129cb32e2e7SValeria Barra 
130cb32e2e7SValeria Barra   // Setup CEED restriction
131*9b072555Sjeremylt   CeedElemRestrictionCreate(ceed, num_elem, P*P*P, num_comp, 1,
132*9b072555Sjeremylt                             m_nodes[0]*m_nodes[1]*m_nodes[2]*num_comp,
133*9b072555Sjeremylt                             CEED_MEM_HOST, CEED_OWN_POINTER, idx, elem_restr);
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;
142*9b072555Sjeremylt   VecScatter l_to_g;              // Scatter for all entries
143*9b072555Sjeremylt   VecScatter l_to_g_0;            // Skip Dirichlet values
144*9b072555Sjeremylt   VecScatter g_to_g_D;            // global-to-global; only Dirichlet values
145*9b072555Sjeremylt   Vec X_loc, Y_loc;
146*9b072555Sjeremylt   CeedVector x_ceed, y_ceed;
147cb32e2e7SValeria Barra   CeedOperator op;
148*9b072555Sjeremylt   CeedVector q_data;
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
156*9b072555Sjeremylt } BPType;
157*9b072555Sjeremylt static const char *const bp_types[] = {"bp1","bp2","bp3","bp4","bp5","bp6",
158*9b072555Sjeremylt                                        "BPType","CEED_BP",0
159cb32e2e7SValeria Barra                                       };
160cb32e2e7SValeria Barra 
161cb32e2e7SValeria Barra // BP specific data
162cb32e2e7SValeria Barra typedef struct {
163*9b072555Sjeremylt   CeedInt num_comp_u, q_data_size, q_extra;
164*9b072555Sjeremylt   CeedQFunctionUser setup_geo, setup_rhs, apply, error;
165*9b072555Sjeremylt   const char *setup_geo_loc, *setup_rhs_loc, *apply_loc, *error_loc;
166*9b072555Sjeremylt   CeedEvalMode in_mode, out_mode;
167*9b072555Sjeremylt   CeedQuadMode q_mode;
168*9b072555Sjeremylt } BPData;
169cb32e2e7SValeria Barra 
170*9b072555Sjeremylt BPData bp_options[6] = {
171cb32e2e7SValeria Barra   [CEED_BP1] = {
172*9b072555Sjeremylt     .num_comp_u = 1,
173*9b072555Sjeremylt     .q_data_size = 1,
174*9b072555Sjeremylt     .q_extra = 1,
175*9b072555Sjeremylt     .setup_geo = SetupMassGeo,
176*9b072555Sjeremylt     .setup_rhs = SetupMassRhs,
177cb32e2e7SValeria Barra     .apply = Mass,
178cb32e2e7SValeria Barra     .error = Error,
179*9b072555Sjeremylt     .setup_geo_loc = SetupMassGeo_loc,
180*9b072555Sjeremylt     .setup_rhs_loc = SetupMassRhs_loc,
181*9b072555Sjeremylt     .apply_loc = Mass_loc,
182*9b072555Sjeremylt     .error_loc = Error_loc,
183*9b072555Sjeremylt     .in_mode = CEED_EVAL_INTERP,
184*9b072555Sjeremylt     .out_mode = CEED_EVAL_INTERP,
185*9b072555Sjeremylt     .q_mode = CEED_GAUSS
186cb32e2e7SValeria Barra   },
187cb32e2e7SValeria Barra   [CEED_BP2] = {
188*9b072555Sjeremylt     .num_comp_u = 3,
189*9b072555Sjeremylt     .q_data_size = 1,
190*9b072555Sjeremylt     .q_extra = 1,
191*9b072555Sjeremylt     .setup_geo = SetupMassGeo,
192*9b072555Sjeremylt     .setup_rhs = SetupMassRhs3,
193cb32e2e7SValeria Barra     .apply = Mass3,
194cb32e2e7SValeria Barra     .error = Error3,
195*9b072555Sjeremylt     .setup_geo_loc = SetupMassGeo_loc,
196*9b072555Sjeremylt     .setup_rhs_loc = SetupMassRhs3_loc,
197*9b072555Sjeremylt     .apply_loc = Mass3_loc,
198*9b072555Sjeremylt     .error_loc = Error3_loc,
199*9b072555Sjeremylt     .in_mode = CEED_EVAL_INTERP,
200*9b072555Sjeremylt     .out_mode = CEED_EVAL_INTERP,
201*9b072555Sjeremylt     .q_mode = CEED_GAUSS
202cb32e2e7SValeria Barra   },
203cb32e2e7SValeria Barra   [CEED_BP3] = {
204*9b072555Sjeremylt     .num_comp_u = 1,
205*9b072555Sjeremylt     .q_data_size = 7,
206*9b072555Sjeremylt     .q_extra = 1,
207*9b072555Sjeremylt     .setup_geo = SetupDiffGeo,
208*9b072555Sjeremylt     .setup_rhs = SetupDiffRhs,
209cb32e2e7SValeria Barra     .apply = Diff,
210cb32e2e7SValeria Barra     .error = Error,
211*9b072555Sjeremylt     .setup_geo_loc = SetupDiffGeo_loc,
212*9b072555Sjeremylt     .setup_rhs_loc = SetupDiffRhs_loc,
213*9b072555Sjeremylt     .apply_loc = Diff_loc,
214*9b072555Sjeremylt     .error_loc = Error_loc,
215*9b072555Sjeremylt     .in_mode = CEED_EVAL_GRAD,
216*9b072555Sjeremylt     .out_mode = CEED_EVAL_GRAD,
217*9b072555Sjeremylt     .q_mode = CEED_GAUSS
218cb32e2e7SValeria Barra   },
219cb32e2e7SValeria Barra   [CEED_BP4] = {
220*9b072555Sjeremylt     .num_comp_u = 3,
221*9b072555Sjeremylt     .q_data_size = 7,
222*9b072555Sjeremylt     .q_extra = 1,
223*9b072555Sjeremylt     .setup_geo = SetupDiffGeo,
224*9b072555Sjeremylt     .setup_rhs = SetupDiffRhs3,
225cb32e2e7SValeria Barra     .apply = Diff3,
226cb32e2e7SValeria Barra     .error = Error3,
227*9b072555Sjeremylt     .setup_geo_loc = SetupDiffGeo_loc,
228*9b072555Sjeremylt     .setup_rhs_loc = SetupDiffRhs3_loc,
229*9b072555Sjeremylt     .apply_loc = Diff3_loc,
230*9b072555Sjeremylt     .error_loc = Error3_loc,
231*9b072555Sjeremylt     .in_mode = CEED_EVAL_GRAD,
232*9b072555Sjeremylt     .out_mode = CEED_EVAL_GRAD,
233*9b072555Sjeremylt     .q_mode = CEED_GAUSS
234cb32e2e7SValeria Barra   },
235cb32e2e7SValeria Barra   [CEED_BP5] = {
236*9b072555Sjeremylt     .num_comp_u = 1,
237*9b072555Sjeremylt     .q_data_size = 7,
238*9b072555Sjeremylt     .q_extra = 0,
239*9b072555Sjeremylt     .setup_geo = SetupDiffGeo,
240*9b072555Sjeremylt     .setup_rhs = SetupDiffRhs,
241cb32e2e7SValeria Barra     .apply = Diff,
242cb32e2e7SValeria Barra     .error = Error,
243*9b072555Sjeremylt     .setup_geo_loc = SetupDiffGeo_loc,
244*9b072555Sjeremylt     .setup_rhs_loc = SetupDiffRhs_loc,
245*9b072555Sjeremylt     .apply_loc = Diff_loc,
246*9b072555Sjeremylt     .error_loc = Error_loc,
247*9b072555Sjeremylt     .in_mode = CEED_EVAL_GRAD,
248*9b072555Sjeremylt     .out_mode = CEED_EVAL_GRAD,
249*9b072555Sjeremylt     .q_mode = CEED_GAUSS_LOBATTO
250cb32e2e7SValeria Barra   },
251cb32e2e7SValeria Barra   [CEED_BP6] = {
252*9b072555Sjeremylt     .num_comp_u = 3,
253*9b072555Sjeremylt     .q_data_size = 7,
254*9b072555Sjeremylt     .q_extra = 0,
255*9b072555Sjeremylt     .setup_geo = SetupDiffGeo,
256*9b072555Sjeremylt     .setup_rhs = SetupDiffRhs3,
257cb32e2e7SValeria Barra     .apply = Diff3,
258cb32e2e7SValeria Barra     .error = Error3,
259*9b072555Sjeremylt     .setup_geo_loc = SetupDiffGeo_loc,
260*9b072555Sjeremylt     .setup_rhs_loc = SetupDiffRhs3_loc,
261*9b072555Sjeremylt     .apply_loc = Diff3_loc,
262*9b072555Sjeremylt     .error_loc = Error3_loc,
263*9b072555Sjeremylt     .in_mode = CEED_EVAL_GRAD,
264*9b072555Sjeremylt     .out_mode = CEED_EVAL_GRAD,
265*9b072555Sjeremylt     .q_mode = 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*9b072555Sjeremylt   PetscMemType x_mem_type, y_mem_type;
275cb32e2e7SValeria Barra 
276cb32e2e7SValeria Barra   PetscFunctionBeginUser;
2779396343dSjeremylt 
278cb32e2e7SValeria Barra   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
2799396343dSjeremylt 
2809396343dSjeremylt   // Global-to-local
281*9b072555Sjeremylt   ierr = VecScatterBegin(user->l_to_g, X, user->X_loc, INSERT_VALUES,
282cb32e2e7SValeria Barra                          SCATTER_REVERSE); CHKERRQ(ierr);
283*9b072555Sjeremylt   ierr = VecScatterEnd(user->l_to_g, X, user->X_loc, INSERT_VALUES,
2849396343dSjeremylt                        SCATTER_REVERSE); CHKERRQ(ierr);
285cb32e2e7SValeria Barra 
2869396343dSjeremylt   // Setup libCEED vectors
287*9b072555Sjeremylt   ierr = VecGetArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x,
288*9b072555Sjeremylt                                    &x_mem_type); CHKERRQ(ierr);
289*9b072555Sjeremylt   ierr = VecGetArrayAndMemType(user->Y_loc, &y, &y_mem_type); CHKERRQ(ierr);
290*9b072555Sjeremylt   CeedVectorSetArray(user->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x);
291*9b072555Sjeremylt   CeedVectorSetArray(user->y_ceed, MemTypeP2C(y_mem_type), CEED_USE_POINTER, y);
292cb32e2e7SValeria Barra 
2939396343dSjeremylt   // Apply libCEED operator
294*9b072555Sjeremylt   CeedOperatorApply(user->op, user->x_ceed, user->y_ceed,
295cb32e2e7SValeria Barra                     CEED_REQUEST_IMMEDIATE);
296cb32e2e7SValeria Barra 
2979396343dSjeremylt   // Restore PETSc vectors
298*9b072555Sjeremylt   CeedVectorTakeArray(user->x_ceed, MemTypeP2C(x_mem_type), NULL);
299*9b072555Sjeremylt   CeedVectorTakeArray(user->y_ceed, MemTypeP2C(y_mem_type), NULL);
300*9b072555Sjeremylt   ierr = VecRestoreArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x);
3019396343dSjeremylt   CHKERRQ(ierr);
302*9b072555Sjeremylt   ierr = VecRestoreArrayAndMemType(user->Y_loc, &y); CHKERRQ(ierr);
303cb32e2e7SValeria Barra 
3049396343dSjeremylt   // Local-to-global
305cb32e2e7SValeria Barra   if (Y) {
306cb32e2e7SValeria Barra     ierr = VecZeroEntries(Y); CHKERRQ(ierr);
307*9b072555Sjeremylt     ierr = VecScatterBegin(user->l_to_g, user->Y_loc, Y, ADD_VALUES,
3089396343dSjeremylt                            SCATTER_FORWARD); CHKERRQ(ierr);
309*9b072555Sjeremylt     ierr = VecScatterEnd(user->l_to_g, user->Y_loc, 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*9b072555Sjeremylt   PetscMemType x_mem_type, y_mem_type;
322cb32e2e7SValeria Barra 
323cb32e2e7SValeria Barra   PetscFunctionBeginUser;
3249396343dSjeremylt 
325cb32e2e7SValeria Barra   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
326cb32e2e7SValeria Barra 
327cb32e2e7SValeria Barra   // Global-to-local
328*9b072555Sjeremylt   ierr = VecScatterBegin(user->l_to_g_0, X, user->X_loc, INSERT_VALUES,
329cb32e2e7SValeria Barra                          SCATTER_REVERSE); CHKERRQ(ierr);
330*9b072555Sjeremylt   ierr = VecScatterEnd(user->l_to_g_0, X, user->X_loc, INSERT_VALUES,
331cb32e2e7SValeria Barra                        SCATTER_REVERSE);
332cb32e2e7SValeria Barra   CHKERRQ(ierr);
333cb32e2e7SValeria Barra 
3349396343dSjeremylt   // Setup libCEED vectors
335*9b072555Sjeremylt   ierr = VecGetArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x,
336*9b072555Sjeremylt                                    &x_mem_type); CHKERRQ(ierr);
337*9b072555Sjeremylt   ierr = VecGetArrayAndMemType(user->Y_loc, &y, &y_mem_type); CHKERRQ(ierr);
338*9b072555Sjeremylt   CeedVectorSetArray(user->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x);
339*9b072555Sjeremylt   CeedVectorSetArray(user->y_ceed, MemTypeP2C(y_mem_type), CEED_USE_POINTER, y);
340cb32e2e7SValeria Barra 
3419396343dSjeremylt   // Apply libCEED operator
342*9b072555Sjeremylt   CeedOperatorApply(user->op, user->x_ceed, user->y_ceed,
343cb32e2e7SValeria Barra                     CEED_REQUEST_IMMEDIATE);
344cb32e2e7SValeria Barra 
345cb32e2e7SValeria Barra   // Restore PETSc vectors
346*9b072555Sjeremylt   CeedVectorTakeArray(user->x_ceed, MemTypeP2C(x_mem_type), NULL);
347*9b072555Sjeremylt   CeedVectorTakeArray(user->y_ceed, MemTypeP2C(y_mem_type), NULL);
348*9b072555Sjeremylt   ierr = VecRestoreArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x);
3499396343dSjeremylt   CHKERRQ(ierr);
350*9b072555Sjeremylt   ierr = VecRestoreArrayAndMemType(user->Y_loc, &y); CHKERRQ(ierr);
351cb32e2e7SValeria Barra 
352cb32e2e7SValeria Barra   // Local-to-global
353cb32e2e7SValeria Barra   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
354*9b072555Sjeremylt   ierr = VecScatterBegin(user->g_to_g_D, X, Y, INSERT_VALUES, SCATTER_FORWARD);
355cb32e2e7SValeria Barra   CHKERRQ(ierr);
356*9b072555Sjeremylt   ierr = VecScatterEnd(user->g_to_g_D, X, Y, INSERT_VALUES, SCATTER_FORWARD);
357cb32e2e7SValeria Barra   CHKERRQ(ierr);
358*9b072555Sjeremylt   ierr = VecScatterBegin(user->l_to_g_0, user->Y_loc, Y, ADD_VALUES,
3599396343dSjeremylt                          SCATTER_FORWARD); CHKERRQ(ierr);
360*9b072555Sjeremylt   ierr = VecScatterEnd(user->l_to_g_0, user->Y_loc, Y, ADD_VALUES,
361*9b072555Sjeremylt                        SCATTER_FORWARD);
362cb32e2e7SValeria Barra   CHKERRQ(ierr);
363cb32e2e7SValeria Barra 
364cb32e2e7SValeria Barra   PetscFunctionReturn(0);
365cb32e2e7SValeria Barra }
366cb32e2e7SValeria Barra 
367cb32e2e7SValeria Barra // This function calculates the error in the final solution
368*9b072555Sjeremylt static PetscErrorCode ComputeErrorMax(User user, CeedOperator op_error, Vec X,
369*9b072555Sjeremylt                                       CeedVector target, PetscReal *max_error) {
370cb32e2e7SValeria Barra   PetscErrorCode ierr;
371cb32e2e7SValeria Barra   PetscScalar *x;
372*9b072555Sjeremylt   PetscMemType mem_type;
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
382*9b072555Sjeremylt   ierr = VecScatterBegin(user->l_to_g, X, user->X_loc, INSERT_VALUES,
383cb32e2e7SValeria Barra                          SCATTER_REVERSE); CHKERRQ(ierr);
384*9b072555Sjeremylt   ierr = VecScatterEnd(user->l_to_g, X, user->X_loc, INSERT_VALUES,
3859396343dSjeremylt                        SCATTER_REVERSE); CHKERRQ(ierr);
3869396343dSjeremylt 
3879396343dSjeremylt   // Setup libCEED vector
388*9b072555Sjeremylt   ierr = VecGetArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x,
389*9b072555Sjeremylt                                    &mem_type); CHKERRQ(ierr);
390*9b072555Sjeremylt   CeedVectorSetArray(user->x_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, x);
391cb32e2e7SValeria Barra 
3929396343dSjeremylt   // Apply libCEED operator
393*9b072555Sjeremylt   CeedOperatorApply(op_error, user->x_ceed, collocated_error,
394cb32e2e7SValeria Barra                     CEED_REQUEST_IMMEDIATE);
395cb32e2e7SValeria Barra 
396cb32e2e7SValeria Barra   // Restore PETSc vector
397*9b072555Sjeremylt   CeedVectorTakeArray(user->x_ceed, MemTypeP2C(mem_type), NULL);
398*9b072555Sjeremylt   ierr = VecRestoreArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x);
3999396343dSjeremylt   CHKERRQ(ierr);
400cb32e2e7SValeria Barra 
401cb32e2e7SValeria Barra   // Reduce max error
402*9b072555Sjeremylt   *max_error = 0;
403cb32e2e7SValeria Barra   const CeedScalar *e;
404cb32e2e7SValeria Barra   CeedVectorGetArrayRead(collocated_error, CEED_MEM_HOST, &e);
405cb32e2e7SValeria Barra   for (CeedInt i=0; i<length; i++) {
406*9b072555Sjeremylt     *max_error = PetscMax(*max_error, PetscAbsScalar(e[i]));
407cb32e2e7SValeria Barra   }
408cb32e2e7SValeria Barra   CeedVectorRestoreArrayRead(collocated_error, &e);
409*9b072555Sjeremylt   ierr = MPI_Allreduce(MPI_IN_PLACE, max_error, 1, MPIU_REAL, MPIU_MAX,
4109396343dSjeremylt                        user->comm); CHKERRQ(ierr);
411cb32e2e7SValeria Barra 
412cb32e2e7SValeria Barra   // Cleanup
413cb32e2e7SValeria Barra   CeedVectorDestroy(&collocated_error);
414cb32e2e7SValeria Barra 
415cb32e2e7SValeria Barra   PetscFunctionReturn(0);
416cb32e2e7SValeria Barra }
417cb32e2e7SValeria Barra 
418cb32e2e7SValeria Barra int main(int argc, char **argv) {
419cb32e2e7SValeria Barra   PetscInt ierr;
420cb32e2e7SValeria Barra   MPI_Comm comm;
421*9b072555Sjeremylt   char ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self";
422cb32e2e7SValeria Barra   double my_rt_start, my_rt, rt_min, rt_max;
423*9b072555Sjeremylt   PetscInt degree, q_extra, local_nodes, local_elem, mesh_elem[3], m_nodes[3],
424*9b072555Sjeremylt            p[3],
425*9b072555Sjeremylt            i_rank[3], l_nodes[3], l_size, num_comp_u = 1, ksp_max_it_clip[2];
426cb32e2e7SValeria Barra   PetscScalar *r;
427cb32e2e7SValeria Barra   PetscBool test_mode, benchmark_mode, write_solution;
428cb32e2e7SValeria Barra   PetscMPIInt size, rank;
429*9b072555Sjeremylt   PetscLogStage solve_stage;
430*9b072555Sjeremylt   Vec X, X_loc, rhs, rhs_loc;
431cb32e2e7SValeria Barra   Mat mat;
432cb32e2e7SValeria Barra   KSP ksp;
433*9b072555Sjeremylt   VecScatter l_to_g, l_to_g_0, g_to_g_D;
434*9b072555Sjeremylt   PetscMemType mem_type;
435cb32e2e7SValeria Barra   User user;
436cb32e2e7SValeria Barra   Ceed ceed;
437*9b072555Sjeremylt   CeedBasis basis_x, basis_u;
438*9b072555Sjeremylt   CeedElemRestriction elem_restr_x, elem_restr_u, elem_restr_u_i, elem_restr_qd_i;
439*9b072555Sjeremylt   CeedQFunction qf_setup_geo, qf_setup_rhs, qf_apply, qf_error;
440*9b072555Sjeremylt   CeedOperator op_setup_geo, op_setup_rhs, op_apply, op_error;
441*9b072555Sjeremylt   CeedVector x_coord, q_data, rhs_ceed, target;
442cb32e2e7SValeria Barra   CeedInt P, Q;
443*9b072555Sjeremylt   const CeedInt dim = 3, num_comp_x = 3;
444*9b072555Sjeremylt   BPType bp_choice;
445cb32e2e7SValeria Barra 
446cb32e2e7SValeria Barra   ierr = PetscInitialize(&argc, &argv, NULL, help);
447cb32e2e7SValeria Barra   if (ierr) return ierr;
448cb32e2e7SValeria Barra   comm = PETSC_COMM_WORLD;
44932d2ee49SValeria Barra 
45032d2ee49SValeria Barra   // Read command line options
451cb32e2e7SValeria Barra   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr);
452*9b072555Sjeremylt   bp_choice = CEED_BP1;
453cb32e2e7SValeria Barra   ierr = PetscOptionsEnum("-problem",
454cb32e2e7SValeria Barra                           "CEED benchmark problem to solve", NULL,
455*9b072555Sjeremylt                           bp_types, (PetscEnum)bp_choice, (PetscEnum *)&bp_choice,
456cb32e2e7SValeria Barra                           NULL); CHKERRQ(ierr);
457*9b072555Sjeremylt   num_comp_u = bp_options[bp_choice].num_comp_u;
458cb32e2e7SValeria Barra   test_mode = PETSC_FALSE;
459cb32e2e7SValeria Barra   ierr = PetscOptionsBool("-test",
460cb32e2e7SValeria Barra                           "Testing mode (do not print unless error is large)",
461cb32e2e7SValeria Barra                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
462cb32e2e7SValeria Barra   benchmark_mode = PETSC_FALSE;
463cb32e2e7SValeria Barra   ierr = PetscOptionsBool("-benchmark",
464cb32e2e7SValeria Barra                           "Benchmarking mode (prints benchmark statistics)",
465cb32e2e7SValeria Barra                           NULL, benchmark_mode, &benchmark_mode, NULL);
466cb32e2e7SValeria Barra   CHKERRQ(ierr);
467cb32e2e7SValeria Barra   write_solution = PETSC_FALSE;
468cb32e2e7SValeria Barra   ierr = PetscOptionsBool("-write_solution",
469cb32e2e7SValeria Barra                           "Write solution for visualization",
470cb32e2e7SValeria Barra                           NULL, write_solution, &write_solution, NULL);
471cb32e2e7SValeria Barra   CHKERRQ(ierr);
472cb32e2e7SValeria Barra   degree = test_mode ? 3 : 1;
473cb32e2e7SValeria Barra   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
474cb32e2e7SValeria Barra                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
475*9b072555Sjeremylt   q_extra = bp_options[bp_choice].q_extra;
476*9b072555Sjeremylt   ierr = PetscOptionsInt("-q_extra", "Number of extra quadrature points",
477*9b072555Sjeremylt                          NULL, q_extra, &q_extra, NULL); CHKERRQ(ierr);
478cb32e2e7SValeria Barra   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
479*9b072555Sjeremylt                             NULL, ceed_resource, ceed_resource,
480*9b072555Sjeremylt                             sizeof(ceed_resource), NULL); CHKERRQ(ierr);
481*9b072555Sjeremylt   local_nodes = 1000;
482cb32e2e7SValeria Barra   ierr = PetscOptionsInt("-local",
483cb32e2e7SValeria Barra                          "Target number of locally owned nodes per process",
484*9b072555Sjeremylt                          NULL, local_nodes, &local_nodes, NULL); CHKERRQ(ierr);
4852fbc6e41SJeremy L Thompson   PetscInt two = 2;
4862fbc6e41SJeremy L Thompson   ksp_max_it_clip[0] = 5;
4872fbc6e41SJeremy L Thompson   ksp_max_it_clip[1] = 20;
4882fbc6e41SJeremy L Thompson   ierr = PetscOptionsIntArray("-ksp_max_it_clip",
4892fbc6e41SJeremy L Thompson                               "Min and max number of iterations to use during benchmarking",
4902fbc6e41SJeremy L Thompson                               NULL, ksp_max_it_clip, &two, NULL);
4912fbc6e41SJeremy L Thompson   CHKERRQ(ierr);
492cb32e2e7SValeria Barra   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
493cb32e2e7SValeria Barra   P = degree + 1;
494*9b072555Sjeremylt   Q = P + q_extra;
495cb32e2e7SValeria Barra 
4969396343dSjeremylt   // Set up libCEED
497*9b072555Sjeremylt   CeedInit(ceed_resource, &ceed);
498*9b072555Sjeremylt   CeedMemType mem_type_backend;
499*9b072555Sjeremylt   CeedGetPreferredMemType(ceed, &mem_type_backend);
5009396343dSjeremylt 
501*9b072555Sjeremylt   VecType default_vec_type = NULL, vec_type;
502*9b072555Sjeremylt   switch (mem_type_backend) {
503*9b072555Sjeremylt   case CEED_MEM_HOST: default_vec_type = VECSTANDARD; break;
504b68a8d79SJed Brown   case CEED_MEM_DEVICE: {
505b68a8d79SJed Brown     const char *resolved;
506b68a8d79SJed Brown     CeedGetResource(ceed, &resolved);
507*9b072555Sjeremylt     if (strstr(resolved, "/gpu/cuda")) default_vec_type = VECCUDA;
508b68a8d79SJed Brown     else if (strstr(resolved, "/gpu/hip/occa"))
509*9b072555Sjeremylt       default_vec_type = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678
510*9b072555Sjeremylt     else if (strstr(resolved, "/gpu/hip")) default_vec_type = VECHIP;
511*9b072555Sjeremylt     else default_vec_type = VECSTANDARD;
512b68a8d79SJed Brown   }
513b68a8d79SJed Brown   }
5149396343dSjeremylt 
515cb32e2e7SValeria Barra   // Determine size of process grid
516cb32e2e7SValeria Barra   ierr = MPI_Comm_size(comm, &size); CHKERRQ(ierr);
517cb32e2e7SValeria Barra   Split3(size, p, false);
518cb32e2e7SValeria Barra 
519*9b072555Sjeremylt   // Find a nicely composite number of elements no less than local_nodes
520*9b072555Sjeremylt   for (local_elem = PetscMax(1, local_nodes / (degree*degree*degree)); ;
521*9b072555Sjeremylt        local_elem++) {
522*9b072555Sjeremylt     Split3(local_elem, mesh_elem, true);
523*9b072555Sjeremylt     if (Max3(mesh_elem) / Min3(mesh_elem) <= 2) break;
524cb32e2e7SValeria Barra   }
525cb32e2e7SValeria Barra 
526cb32e2e7SValeria Barra   // Find my location in the process grid
527cb32e2e7SValeria Barra   ierr = MPI_Comm_rank(comm, &rank); CHKERRQ(ierr);
528*9b072555Sjeremylt   for (int d=0, rank_left=rank; d<dim; d++) {
529cb32e2e7SValeria Barra     const int pstride[3] = {p[1] *p[2], p[2], 1};
530*9b072555Sjeremylt     i_rank[d] = rank_left / pstride[d];
531*9b072555Sjeremylt     rank_left -= i_rank[d] * pstride[d];
532cb32e2e7SValeria Barra   }
533cb32e2e7SValeria Barra 
534*9b072555Sjeremylt   GlobalNodes(p, i_rank, degree, mesh_elem, m_nodes);
535cb32e2e7SValeria Barra 
536cb32e2e7SValeria Barra   // Setup global vector
5374a5da23aSJed Brown   ierr = VecCreate(comm, &X); CHKERRQ(ierr);
538*9b072555Sjeremylt   ierr = VecSetType(X, default_vec_type); CHKERRQ(ierr);
539*9b072555Sjeremylt   ierr = VecSetSizes(X, m_nodes[0]*m_nodes[1]*m_nodes[2]*num_comp_u,
540*9b072555Sjeremylt                      PETSC_DECIDE);
541cb32e2e7SValeria Barra   CHKERRQ(ierr);
542b68a8d79SJed Brown   ierr = VecSetFromOptions(X); CHKERRQ(ierr);
543cb32e2e7SValeria Barra   ierr = VecSetUp(X); CHKERRQ(ierr);
544cb32e2e7SValeria Barra 
545cb32e2e7SValeria Barra   // Set up libCEED
546*9b072555Sjeremylt   CeedInit(ceed_resource, &ceed);
547cb32e2e7SValeria Barra 
548cb32e2e7SValeria Barra   // Print summary
549cb32e2e7SValeria Barra   CeedInt gsize;
550cb32e2e7SValeria Barra   ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
551dc7d240cSValeria Barra   if (!test_mode) {
552*9b072555Sjeremylt     const char *used_resource;
553*9b072555Sjeremylt     CeedGetResource(ceed, &used_resource);
5549396343dSjeremylt 
555*9b072555Sjeremylt     ierr = VecGetType(X, &vec_type); CHKERRQ(ierr);
5569396343dSjeremylt 
557cb32e2e7SValeria Barra     ierr = PetscPrintf(comm,
558cb32e2e7SValeria Barra                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n"
5599396343dSjeremylt                        "  PETSc:\n"
5609396343dSjeremylt                        "    PETSc Vec Type                     : %s\n"
561cb32e2e7SValeria Barra                        "  libCEED:\n"
562cb32e2e7SValeria Barra                        "    libCEED Backend                    : %s\n"
5639396343dSjeremylt                        "    libCEED Backend MemType            : %s\n"
564cb32e2e7SValeria Barra                        "  Mesh:\n"
5658d0bb2bbSvaleriabarra                        "    Number of 1D Basis Nodes (P)       : %d\n"
5668d0bb2bbSvaleriabarra                        "    Number of 1D Quadrature Points (Q) : %d\n"
567cb32e2e7SValeria Barra                        "    Global nodes                       : %D\n"
568cb32e2e7SValeria Barra                        "    Process Decomposition              : %D %D %D\n"
569cb32e2e7SValeria Barra                        "    Local Elements                     : %D = %D %D %D\n"
570db419314Sjeremylt                        "    Owned nodes                        : %D = %D %D %D\n"
571db419314Sjeremylt                        "    DoF per node                       : %D\n",
572*9b072555Sjeremylt                        bp_choice+1, vec_type, used_resource,
573*9b072555Sjeremylt                        CeedMemTypes[mem_type_backend],
574*9b072555Sjeremylt                        P, Q,  gsize/num_comp_u, p[0], p[1], p[2], local_elem,
575*9b072555Sjeremylt                        mesh_elem[0], mesh_elem[1], mesh_elem[2],
576*9b072555Sjeremylt                        m_nodes[0]*m_nodes[1]*m_nodes[2], m_nodes[0], m_nodes[1],
577*9b072555Sjeremylt                        m_nodes[2], num_comp_u); CHKERRQ(ierr);
578cb32e2e7SValeria Barra   }
579cb32e2e7SValeria Barra 
580cb32e2e7SValeria Barra   {
581*9b072555Sjeremylt     l_size = 1;
582cb32e2e7SValeria Barra     for (int d=0; d<dim; d++) {
583*9b072555Sjeremylt       l_nodes[d] = mesh_elem[d]*degree + 1;
584*9b072555Sjeremylt       l_size *= l_nodes[d];
585cb32e2e7SValeria Barra     }
586*9b072555Sjeremylt     ierr = VecCreate(PETSC_COMM_SELF, &X_loc); CHKERRQ(ierr);
587*9b072555Sjeremylt     ierr = VecSetType(X_loc, default_vec_type); CHKERRQ(ierr);
588*9b072555Sjeremylt     ierr = VecSetSizes(X_loc, l_size*num_comp_u, PETSC_DECIDE); CHKERRQ(ierr);
589*9b072555Sjeremylt     ierr = VecSetFromOptions(X_loc); CHKERRQ(ierr);
590*9b072555Sjeremylt     ierr = VecSetUp(X_loc); CHKERRQ(ierr);
591cb32e2e7SValeria Barra 
592cb32e2e7SValeria Barra     // Create local-to-global scatter
593*9b072555Sjeremylt     PetscInt *l_to_g_ind, *l_to_g_ind_0, *loc_ind, l_0_count;
594*9b072555Sjeremylt     IS l_to_g_is, l_to_g_is_0, loc_is;
595*9b072555Sjeremylt     PetscInt g_start[2][2][2], g_m_nodes[2][2][2][dim];
596cb32e2e7SValeria Barra 
597cb32e2e7SValeria Barra     for (int i=0; i<2; i++) {
598cb32e2e7SValeria Barra       for (int j=0; j<2; j++) {
599cb32e2e7SValeria Barra         for (int k=0; k<2; k++) {
600*9b072555Sjeremylt           PetscInt ijk_rank[3] = {i_rank[0]+i, i_rank[1]+j, i_rank[2]+k};
601*9b072555Sjeremylt           g_start[i][j][k] = GlobalStart(p, ijk_rank, degree, mesh_elem);
602*9b072555Sjeremylt           GlobalNodes(p, ijk_rank, degree, mesh_elem, g_m_nodes[i][j][k]);
603cb32e2e7SValeria Barra         }
604cb32e2e7SValeria Barra       }
605cb32e2e7SValeria Barra     }
606cb32e2e7SValeria Barra 
607*9b072555Sjeremylt     ierr = PetscMalloc1(l_size, &l_to_g_ind); CHKERRQ(ierr);
608*9b072555Sjeremylt     ierr = PetscMalloc1(l_size, &l_to_g_ind_0); CHKERRQ(ierr);
609*9b072555Sjeremylt     ierr = PetscMalloc1(l_size, &loc_ind); CHKERRQ(ierr);
610*9b072555Sjeremylt     l_0_count = 0;
611*9b072555Sjeremylt     for (PetscInt i=0,ir,ii; ir=i>=m_nodes[0], ii=i-ir*m_nodes[0], i<l_nodes[0];
612*9b072555Sjeremylt          i++)
613*9b072555Sjeremylt       for (PetscInt j=0,jr,jj; jr=j>=m_nodes[1], jj=j-jr*m_nodes[1], j<l_nodes[1];
614*9b072555Sjeremylt            j++)
615*9b072555Sjeremylt         for (PetscInt k=0,kr,kk; kr=k>=m_nodes[2], kk=k-kr*m_nodes[2], k<l_nodes[2];
616*9b072555Sjeremylt              k++) {
617*9b072555Sjeremylt           PetscInt here = (i*l_nodes[1]+j)*l_nodes[2]+k;
618*9b072555Sjeremylt           l_to_g_ind[here] =
619*9b072555Sjeremylt             g_start[ir][jr][kr] + (ii*g_m_nodes[ir][jr][kr][1]+jj)*g_m_nodes[ir][jr][kr][2]
620*9b072555Sjeremylt             +kk;
621*9b072555Sjeremylt           if ((i_rank[0] == 0 && i == 0)
622*9b072555Sjeremylt               || (i_rank[1] == 0 && j == 0)
623*9b072555Sjeremylt               || (i_rank[2] == 0 && k == 0)
624*9b072555Sjeremylt               || (i_rank[0]+1 == p[0] && i+1 == l_nodes[0])
625*9b072555Sjeremylt               || (i_rank[1]+1 == p[1] && j+1 == l_nodes[1])
626*9b072555Sjeremylt               || (i_rank[2]+1 == p[2] && k+1 == l_nodes[2]))
627cb32e2e7SValeria Barra             continue;
628*9b072555Sjeremylt           l_to_g_ind_0[l_0_count] = l_to_g_ind[here];
629*9b072555Sjeremylt           loc_ind[l_0_count++] = here;
630cb32e2e7SValeria Barra         }
631*9b072555Sjeremylt     ierr = ISCreateBlock(comm, num_comp_u, l_size, l_to_g_ind, PETSC_OWN_POINTER,
632*9b072555Sjeremylt                          &l_to_g_is); CHKERRQ(ierr);
633*9b072555Sjeremylt     ierr = VecScatterCreate(X_loc, NULL, X, l_to_g_is, &l_to_g); CHKERRQ(ierr);
634cb32e2e7SValeria Barra     CHKERRQ(ierr);
635*9b072555Sjeremylt     ierr = ISCreateBlock(comm, num_comp_u, l_0_count, l_to_g_ind_0,
636*9b072555Sjeremylt                          PETSC_OWN_POINTER,
637*9b072555Sjeremylt                          &l_to_g_is_0); CHKERRQ(ierr);
638*9b072555Sjeremylt     ierr = ISCreateBlock(comm, num_comp_u, l_0_count, loc_ind, PETSC_OWN_POINTER,
639*9b072555Sjeremylt                          &loc_is); CHKERRQ(ierr);
640*9b072555Sjeremylt     ierr = VecScatterCreate(X_loc, loc_is, X, l_to_g_is_0, &l_to_g_0);
641*9b072555Sjeremylt     CHKERRQ(ierr);
642cb32e2e7SValeria Barra     {
643cb32e2e7SValeria Barra       // Create global-to-global scatter for Dirichlet values (everything not in
644*9b072555Sjeremylt       // l_to_g_is_0, which is the range of l_to_g_0)
645*9b072555Sjeremylt       PetscInt x_start, x_end, *ind_D, count_D = 0;
646*9b072555Sjeremylt       IS is_D;
647cb32e2e7SValeria Barra       const PetscScalar *x;
648*9b072555Sjeremylt       ierr = VecZeroEntries(X_loc); CHKERRQ(ierr);
649cb32e2e7SValeria Barra       ierr = VecSet(X, 1.0); CHKERRQ(ierr);
650*9b072555Sjeremylt       ierr = VecScatterBegin(l_to_g_0, X_loc, X, INSERT_VALUES, SCATTER_FORWARD);
651cb32e2e7SValeria Barra       CHKERRQ(ierr);
652*9b072555Sjeremylt       ierr = VecScatterEnd(l_to_g_0, X_loc, X, INSERT_VALUES, SCATTER_FORWARD);
653cb32e2e7SValeria Barra       CHKERRQ(ierr);
654*9b072555Sjeremylt       ierr = VecGetOwnershipRange(X, &x_start, &x_end); CHKERRQ(ierr);
655*9b072555Sjeremylt       ierr = PetscMalloc1(x_end-x_start, &ind_D); CHKERRQ(ierr);
656cb32e2e7SValeria Barra       ierr = VecGetArrayRead(X, &x); CHKERRQ(ierr);
657*9b072555Sjeremylt       for (PetscInt i=0; i<x_end-x_start; i++) {
658*9b072555Sjeremylt         if (x[i] == 1.) ind_D[count_D++] = x_start + i;
659cb32e2e7SValeria Barra       }
660cb32e2e7SValeria Barra       ierr = VecRestoreArrayRead(X, &x); CHKERRQ(ierr);
661*9b072555Sjeremylt       ierr = ISCreateGeneral(comm, count_D, ind_D, PETSC_COPY_VALUES, &is_D);
662cb32e2e7SValeria Barra       CHKERRQ(ierr);
663*9b072555Sjeremylt       ierr = PetscFree(ind_D); CHKERRQ(ierr);
664*9b072555Sjeremylt       ierr = VecScatterCreate(X, is_D, X, is_D, &g_to_g_D); CHKERRQ(ierr);
665*9b072555Sjeremylt       ierr = ISDestroy(&is_D); CHKERRQ(ierr);
666cb32e2e7SValeria Barra     }
667*9b072555Sjeremylt     ierr = ISDestroy(&l_to_g_is); CHKERRQ(ierr);
668*9b072555Sjeremylt     ierr = ISDestroy(&l_to_g_is_0); CHKERRQ(ierr);
669*9b072555Sjeremylt     ierr = ISDestroy(&loc_is); CHKERRQ(ierr);
670cb32e2e7SValeria Barra   }
671cb32e2e7SValeria Barra 
672cb32e2e7SValeria Barra   // CEED bases
673*9b072555Sjeremylt   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_u, P, Q,
674*9b072555Sjeremylt                                   bp_options[bp_choice].q_mode, &basis_u);
675*9b072555Sjeremylt   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q,
676*9b072555Sjeremylt                                   bp_options[bp_choice].q_mode, &basis_x);
677cb32e2e7SValeria Barra 
678cb32e2e7SValeria Barra   // CEED restrictions
679*9b072555Sjeremylt   CreateRestriction(ceed, mesh_elem, P, num_comp_u, &elem_restr_u);
680*9b072555Sjeremylt   CreateRestriction(ceed, mesh_elem, 2, dim, &elem_restr_x);
681*9b072555Sjeremylt   CeedInt num_elem = mesh_elem[0]*mesh_elem[1]*mesh_elem[2];
682*9b072555Sjeremylt   CeedElemRestrictionCreateStrided(ceed, num_elem, Q*Q*Q, num_comp_u,
683*9b072555Sjeremylt                                    num_comp_u*num_elem*Q*Q*Q,
684*9b072555Sjeremylt                                    CEED_STRIDES_BACKEND, &elem_restr_u_i);
685*9b072555Sjeremylt   CeedElemRestrictionCreateStrided(ceed, num_elem, Q*Q*Q,
686*9b072555Sjeremylt                                    bp_options[bp_choice].q_data_size,
687*9b072555Sjeremylt                                    bp_options[bp_choice].q_data_size*num_elem*Q*Q*Q,
688*9b072555Sjeremylt                                    CEED_STRIDES_BACKEND, &elem_restr_qd_i);
689cb32e2e7SValeria Barra   {
690*9b072555Sjeremylt     CeedScalar *x_loc;
691*9b072555Sjeremylt     CeedInt shape[3] = {mesh_elem[0]+1, mesh_elem[1]+1, mesh_elem[2]+1}, len =
692cb32e2e7SValeria Barra                          shape[0]*shape[1]*shape[2];
693*9b072555Sjeremylt     x_loc = malloc(len*num_comp_x*sizeof x_loc[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++) {
697*9b072555Sjeremylt           x_loc[dim*((i*shape[1]+j)*shape[2]+k) + 0] = 1.*(i_rank[0]*mesh_elem[0]+i) /
698*9b072555Sjeremylt               (p[0]*mesh_elem[0]);
699*9b072555Sjeremylt           x_loc[dim*((i*shape[1]+j)*shape[2]+k) + 1] = 1.*(i_rank[1]*mesh_elem[1]+j) /
700*9b072555Sjeremylt               (p[1]*mesh_elem[1]);
701*9b072555Sjeremylt           x_loc[dim*((i*shape[1]+j)*shape[2]+k) + 2] = 1.*(i_rank[2]*mesh_elem[2]+k) /
702*9b072555Sjeremylt               (p[2]*mesh_elem[2]);
703cb32e2e7SValeria Barra         }
704cb32e2e7SValeria Barra       }
705cb32e2e7SValeria Barra     }
706*9b072555Sjeremylt     CeedVectorCreate(ceed, len*num_comp_x, &x_coord);
707*9b072555Sjeremylt     CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_OWN_POINTER, x_loc);
708cb32e2e7SValeria Barra   }
709cb32e2e7SValeria Barra 
710*9b072555Sjeremylt   // Create the QFunction that builds the operator quadrature data
711*9b072555Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].setup_geo,
712*9b072555Sjeremylt                               bp_options[bp_choice].setup_geo_loc, &qf_setup_geo);
713*9b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_geo, "x", num_comp_x, CEED_EVAL_INTERP);
714*9b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_geo, "dx", num_comp_x*dim, CEED_EVAL_GRAD);
715*9b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_geo, "weight", 1, CEED_EVAL_WEIGHT);
716*9b072555Sjeremylt   CeedQFunctionAddOutput(qf_setup_geo, "q_data",
717*9b072555Sjeremylt                          bp_options[bp_choice].q_data_size,
718cb32e2e7SValeria Barra                          CEED_EVAL_NONE);
719cb32e2e7SValeria Barra 
720*9b072555Sjeremylt   // Create the QFunction that sets up the RHS and true solution
721*9b072555Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].setup_rhs,
722*9b072555Sjeremylt                               bp_options[bp_choice].setup_rhs_loc, &qf_setup_rhs);
723*9b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_rhs, "x", num_comp_x, CEED_EVAL_INTERP);
724*9b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_rhs, "q_data", bp_options[bp_choice].q_data_size,
725e83e87a5Sjeremylt                         CEED_EVAL_NONE);
726*9b072555Sjeremylt   CeedQFunctionAddOutput(qf_setup_rhs, "true_soln", num_comp_u, CEED_EVAL_NONE);
727*9b072555Sjeremylt   CeedQFunctionAddOutput(qf_setup_rhs, "rhs", num_comp_u, CEED_EVAL_INTERP);
728cb32e2e7SValeria Barra 
729cb32e2e7SValeria Barra   // Set up PDE operator
730*9b072555Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].apply,
731*9b072555Sjeremylt                               bp_options[bp_choice].apply_loc, &qf_apply);
732cb32e2e7SValeria Barra   // Add inputs and outputs
733*9b072555Sjeremylt   CeedInt in_scale = bp_options[bp_choice].in_mode==CEED_EVAL_GRAD ? 3 : 1;
734*9b072555Sjeremylt   CeedInt out_scale = bp_options[bp_choice].out_mode==CEED_EVAL_GRAD ? 3 : 1;
735*9b072555Sjeremylt   CeedQFunctionAddInput(qf_apply, "u", num_comp_u*in_scale,
736*9b072555Sjeremylt                         bp_options[bp_choice].in_mode);
737*9b072555Sjeremylt   CeedQFunctionAddInput(qf_apply, "q_data", bp_options[bp_choice].q_data_size,
738cb32e2e7SValeria Barra                         CEED_EVAL_NONE);
739*9b072555Sjeremylt   CeedQFunctionAddOutput(qf_apply, "v", num_comp_u*out_scale,
740*9b072555Sjeremylt                          bp_options[bp_choice].out_mode);
741cb32e2e7SValeria Barra 
742cb32e2e7SValeria Barra   // Create the error qfunction
743*9b072555Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].error,
744*9b072555Sjeremylt                               bp_options[bp_choice].error_loc, &qf_error);
745*9b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP);
746*9b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE);
747*9b072555Sjeremylt   CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_NONE);
748cb32e2e7SValeria Barra 
749cb32e2e7SValeria Barra   // Create the persistent vectors that will be needed in setup
750*9b072555Sjeremylt   CeedInt num_qpts;
751*9b072555Sjeremylt   CeedBasisGetNumQuadraturePoints(basis_u, &num_qpts);
752*9b072555Sjeremylt   CeedVectorCreate(ceed, bp_options[bp_choice].q_data_size*num_elem*num_qpts,
753*9b072555Sjeremylt                    &q_data);
754*9b072555Sjeremylt   CeedVectorCreate(ceed, num_elem*num_qpts*num_comp_u, &target);
755*9b072555Sjeremylt   CeedVectorCreate(ceed, l_size*num_comp_u, &rhs_ceed);
756cb32e2e7SValeria Barra 
757cb32e2e7SValeria Barra   // Create the operator that builds the quadrature data for the ceed operator
758*9b072555Sjeremylt   CeedOperatorCreate(ceed, qf_setup_geo, CEED_QFUNCTION_NONE,
759*9b072555Sjeremylt                      CEED_QFUNCTION_NONE, &op_setup_geo);
760*9b072555Sjeremylt   CeedOperatorSetField(op_setup_geo, "x", elem_restr_x, basis_x,
761e83e87a5Sjeremylt                        CEED_VECTOR_ACTIVE);
762*9b072555Sjeremylt   CeedOperatorSetField(op_setup_geo, "dx", elem_restr_x, basis_x,
763a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
764*9b072555Sjeremylt   CeedOperatorSetField(op_setup_geo, "weight", CEED_ELEMRESTRICTION_NONE, basis_x,
765a8d32208Sjeremylt                        CEED_VECTOR_NONE);
766*9b072555Sjeremylt   CeedOperatorSetField(op_setup_geo, "q_data", elem_restr_qd_i,
767cb32e2e7SValeria Barra                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
768cb32e2e7SValeria Barra 
769cb32e2e7SValeria Barra   // Create the operator that builds the RHS and true solution
770*9b072555Sjeremylt   CeedOperatorCreate(ceed, qf_setup_rhs, CEED_QFUNCTION_NONE,
771*9b072555Sjeremylt                      CEED_QFUNCTION_NONE, &op_setup_rhs);
772*9b072555Sjeremylt   CeedOperatorSetField(op_setup_rhs, "x", elem_restr_x, basis_x,
773a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
774*9b072555Sjeremylt   CeedOperatorSetField(op_setup_rhs, "q_data", elem_restr_qd_i,
775*9b072555Sjeremylt                        CEED_BASIS_COLLOCATED,
776*9b072555Sjeremylt                        q_data);
777*9b072555Sjeremylt   CeedOperatorSetField(op_setup_rhs, "true_soln", elem_restr_u_i,
778cb32e2e7SValeria Barra                        CEED_BASIS_COLLOCATED, target);
779*9b072555Sjeremylt   CeedOperatorSetField(op_setup_rhs, "rhs", elem_restr_u, basis_u,
780a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
781cb32e2e7SValeria Barra 
782cb32e2e7SValeria Barra   // Create the mass or diff operator
783*9b072555Sjeremylt   CeedOperatorCreate(ceed, qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
784*9b072555Sjeremylt                      &op_apply);
785*9b072555Sjeremylt   CeedOperatorSetField(op_apply, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
786*9b072555Sjeremylt   CeedOperatorSetField(op_apply, "q_data", elem_restr_qd_i, CEED_BASIS_COLLOCATED,
787*9b072555Sjeremylt                        q_data);
788*9b072555Sjeremylt   CeedOperatorSetField(op_apply, "v", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
789cb32e2e7SValeria Barra 
790cb32e2e7SValeria Barra   // Create the error operator
791*9b072555Sjeremylt   CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
792*9b072555Sjeremylt                      &op_error);
793*9b072555Sjeremylt   CeedOperatorSetField(op_error, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
794*9b072555Sjeremylt   CeedOperatorSetField(op_error, "true_soln", elem_restr_u_i,
795cb32e2e7SValeria Barra                        CEED_BASIS_COLLOCATED, target);
796*9b072555Sjeremylt   CeedOperatorSetField(op_error, "error", elem_restr_u_i, CEED_BASIS_COLLOCATED,
797a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
798cb32e2e7SValeria Barra 
799cb32e2e7SValeria Barra   // Set up Mat
800cb32e2e7SValeria Barra   ierr = PetscMalloc1(1, &user); CHKERRQ(ierr);
801cb32e2e7SValeria Barra   user->comm = comm;
802*9b072555Sjeremylt   user->l_to_g = l_to_g;
803*9b072555Sjeremylt   if (bp_choice != CEED_BP1 && bp_choice != CEED_BP2) {
804*9b072555Sjeremylt     user->l_to_g_0 = l_to_g_0;
805*9b072555Sjeremylt     user->g_to_g_D = g_to_g_D;
806cb32e2e7SValeria Barra   }
807*9b072555Sjeremylt   user->X_loc = X_loc;
808*9b072555Sjeremylt   ierr = VecDuplicate(X_loc, &user->Y_loc); CHKERRQ(ierr);
809*9b072555Sjeremylt   CeedVectorCreate(ceed, l_size*num_comp_u, &user->x_ceed);
810*9b072555Sjeremylt   CeedVectorCreate(ceed, l_size*num_comp_u, &user->y_ceed);
811*9b072555Sjeremylt   user->op = op_apply;
812*9b072555Sjeremylt   user->q_data = q_data;
813cb32e2e7SValeria Barra   user->ceed = ceed;
814cb32e2e7SValeria Barra 
815*9b072555Sjeremylt   ierr = MatCreateShell(comm, m_nodes[0]*m_nodes[1]*m_nodes[2]*num_comp_u,
816*9b072555Sjeremylt                         m_nodes[0]*m_nodes[1]*m_nodes[2]*num_comp_u,
817cb32e2e7SValeria Barra                         PETSC_DECIDE, PETSC_DECIDE, user, &mat); CHKERRQ(ierr);
818*9b072555Sjeremylt   if (bp_choice == CEED_BP1 || bp_choice == CEED_BP2) {
819cb32e2e7SValeria Barra     ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Mass);
820cb32e2e7SValeria Barra     CHKERRQ(ierr);
821cb32e2e7SValeria Barra   } else {
822cb32e2e7SValeria Barra     ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Diff);
823cb32e2e7SValeria Barra     CHKERRQ(ierr);
824cb32e2e7SValeria Barra   }
825*9b072555Sjeremylt   ierr = VecGetType(X, &vec_type); CHKERRQ(ierr);
826*9b072555Sjeremylt   ierr = MatShellSetVecType(mat, vec_type); CHKERRQ(ierr);
827cb32e2e7SValeria Barra 
828cb32e2e7SValeria Barra   // Get RHS vector
8299396343dSjeremylt   ierr = VecDuplicate(X, &rhs); CHKERRQ(ierr);
830*9b072555Sjeremylt   ierr = VecDuplicate(X_loc, &rhs_loc); CHKERRQ(ierr);
831*9b072555Sjeremylt   ierr = VecZeroEntries(rhs_loc); CHKERRQ(ierr);
832*9b072555Sjeremylt   ierr = VecGetArrayAndMemType(rhs_loc, &r, &mem_type); CHKERRQ(ierr);
833*9b072555Sjeremylt   CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r);
834cb32e2e7SValeria Barra 
835*9b072555Sjeremylt   // Setup q_data, rhs, and target
836*9b072555Sjeremylt   CeedOperatorApply(op_setup_geo, x_coord, q_data, CEED_REQUEST_IMMEDIATE);
837*9b072555Sjeremylt   CeedOperatorApply(op_setup_rhs, x_coord, rhs_ceed, CEED_REQUEST_IMMEDIATE);
838*9b072555Sjeremylt   CeedVectorDestroy(&x_coord);
839cb32e2e7SValeria Barra 
840cb32e2e7SValeria Barra   // Gather RHS
841*9b072555Sjeremylt   ierr = CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL); CHKERRQ(ierr);
842*9b072555Sjeremylt   ierr = VecRestoreArrayAndMemType(rhs_loc, &r); CHKERRQ(ierr);
843cb32e2e7SValeria Barra   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
844*9b072555Sjeremylt   ierr = VecScatterBegin(l_to_g, rhs_loc, rhs, ADD_VALUES, SCATTER_FORWARD);
845cb32e2e7SValeria Barra   CHKERRQ(ierr);
846*9b072555Sjeremylt   ierr = VecScatterEnd(l_to_g, rhs_loc, rhs, ADD_VALUES, SCATTER_FORWARD);
847cb32e2e7SValeria Barra   CHKERRQ(ierr);
848*9b072555Sjeremylt   CeedVectorDestroy(&rhs_ceed);
849cb32e2e7SValeria Barra 
850cb32e2e7SValeria Barra   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
851cb32e2e7SValeria Barra   {
852cb32e2e7SValeria Barra     PC pc;
853cb32e2e7SValeria Barra     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
854*9b072555Sjeremylt     if (bp_choice == CEED_BP1 || bp_choice == CEED_BP2) {
855cb32e2e7SValeria Barra       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
856cb32e2e7SValeria Barra       ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr);
857cb32e2e7SValeria Barra     } else {
858cb32e2e7SValeria Barra       ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr);
859cb32e2e7SValeria Barra     }
860cb32e2e7SValeria Barra     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
861cb32e2e7SValeria Barra     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
862cb32e2e7SValeria Barra     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
863cb32e2e7SValeria Barra                             PETSC_DEFAULT); CHKERRQ(ierr);
864cb32e2e7SValeria Barra   }
865cb32e2e7SValeria Barra   ierr = KSPSetOperators(ksp, mat, mat); CHKERRQ(ierr);
866da9108adSvaleriabarra   // First run's performance log is not considered for benchmarking purposes
867cb32e2e7SValeria Barra   ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
868cb32e2e7SValeria Barra   CHKERRQ(ierr);
869cb32e2e7SValeria Barra   my_rt_start = MPI_Wtime();
870cb32e2e7SValeria Barra   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
871cb32e2e7SValeria Barra   my_rt = MPI_Wtime() - my_rt_start;
872cb32e2e7SValeria Barra   ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
873cb32e2e7SValeria Barra   CHKERRQ(ierr);
874cb32e2e7SValeria Barra   // Set maxits based on first iteration timing
875cb32e2e7SValeria Barra   if (my_rt > 0.02) {
8762fbc6e41SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
8772fbc6e41SJeremy L Thompson                             ksp_max_it_clip[0]); CHKERRQ(ierr);
878cb32e2e7SValeria Barra   } else {
8792fbc6e41SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
8802fbc6e41SJeremy L Thompson                             ksp_max_it_clip[1]); CHKERRQ(ierr);
881cb32e2e7SValeria Barra   }
882debcf919SJed Brown   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
88309a940d7Sjeremylt 
884cb32e2e7SValeria Barra   // Timed solve
88509a940d7Sjeremylt   ierr = VecZeroEntries(X); CHKERRQ(ierr);
886cb32e2e7SValeria Barra   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
88709a940d7Sjeremylt 
88809a940d7Sjeremylt   // -- Performance logging
889*9b072555Sjeremylt   ierr = PetscLogStageRegister("Solve Stage", &solve_stage); CHKERRQ(ierr);
890*9b072555Sjeremylt   ierr = PetscLogStagePush(solve_stage); CHKERRQ(ierr);
89109a940d7Sjeremylt 
89209a940d7Sjeremylt   // -- Solve
893cb32e2e7SValeria Barra   my_rt_start = MPI_Wtime();
894cb32e2e7SValeria Barra   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
895cb32e2e7SValeria Barra   my_rt = MPI_Wtime() - my_rt_start;
89609a940d7Sjeremylt 
89709a940d7Sjeremylt   // -- Performance logging
89809a940d7Sjeremylt   ierr = PetscLogStagePop();
89909a940d7Sjeremylt 
9008e87e98bSJed Brown   // Output results
901cb32e2e7SValeria Barra   {
902*9b072555Sjeremylt     KSPType ksp_type;
903cb32e2e7SValeria Barra     KSPConvergedReason reason;
904cb32e2e7SValeria Barra     PetscReal rnorm;
905cb32e2e7SValeria Barra     PetscInt its;
906*9b072555Sjeremylt     ierr = KSPGetType(ksp, &ksp_type); CHKERRQ(ierr);
907cb32e2e7SValeria Barra     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
908cb32e2e7SValeria Barra     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
909cb32e2e7SValeria Barra     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
9100ec79729SJeremy L Thompson     if (!test_mode || reason < 0 || rnorm > 1e-9) {
911cb32e2e7SValeria Barra       ierr = PetscPrintf(comm,
912cb32e2e7SValeria Barra                          "  KSP:\n"
913cb32e2e7SValeria Barra                          "    KSP Type                           : %s\n"
914cb32e2e7SValeria Barra                          "    KSP Convergence                    : %s\n"
915cb32e2e7SValeria Barra                          "    Total KSP Iterations               : %D\n"
916cb32e2e7SValeria Barra                          "    Final rnorm                        : %e\n",
917*9b072555Sjeremylt                          ksp_type, KSPConvergedReasons[reason], its,
918cb32e2e7SValeria Barra                          (double)rnorm); CHKERRQ(ierr);
919cb32e2e7SValeria Barra     }
9208e87e98bSJed Brown     if (!test_mode) {
9218e87e98bSJed Brown       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
9228e87e98bSJed Brown     }
9238e87e98bSJed Brown     {
924*9b072555Sjeremylt       PetscReal max_error;
925*9b072555Sjeremylt       ierr = ComputeErrorMax(user, op_error, X, target, &max_error);
9268e87e98bSJed Brown       CHKERRQ(ierr);
9278e87e98bSJed Brown       PetscReal tol = 5e-2;
928*9b072555Sjeremylt       if (!test_mode || max_error > tol) {
929cb32e2e7SValeria Barra         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
930cb32e2e7SValeria Barra         CHKERRQ(ierr);
931cb32e2e7SValeria Barra         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
932cb32e2e7SValeria Barra         CHKERRQ(ierr);
933cb32e2e7SValeria Barra         ierr = PetscPrintf(comm,
9348e87e98bSJed Brown                            "    Pointwise Error (max)              : %e\n"
9358e87e98bSJed Brown                            "    CG Solve Time                      : %g (%g) sec\n",
936*9b072555Sjeremylt                            (double)max_error, rt_max, rt_min); CHKERRQ(ierr);
937cb32e2e7SValeria Barra       }
938cb32e2e7SValeria Barra     }
9398d0bb2bbSvaleriabarra     if (!test_mode) {
940cb32e2e7SValeria Barra       ierr = PetscPrintf(comm,
9418e87e98bSJed Brown                          "    DoFs/Sec in CG                     : %g (%g) million\n",
9428e87e98bSJed Brown                          1e-6*gsize*its/rt_max,
9438e87e98bSJed Brown                          1e-6*gsize*its/rt_min); CHKERRQ(ierr);
944cb32e2e7SValeria Barra     }
945cb32e2e7SValeria Barra   }
946cb32e2e7SValeria Barra 
947cb32e2e7SValeria Barra   if (write_solution) {
948*9b072555Sjeremylt     PetscViewer vtk_viewer_soln;
949cb32e2e7SValeria Barra 
950*9b072555Sjeremylt     ierr = PetscViewerCreate(comm, &vtk_viewer_soln); CHKERRQ(ierr);
951*9b072555Sjeremylt     ierr = PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK); CHKERRQ(ierr);
952*9b072555Sjeremylt     ierr = PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"); CHKERRQ(ierr);
953*9b072555Sjeremylt     ierr = VecView(X, vtk_viewer_soln); CHKERRQ(ierr);
954*9b072555Sjeremylt     ierr = PetscViewerDestroy(&vtk_viewer_soln); CHKERRQ(ierr);
955cb32e2e7SValeria Barra   }
956cb32e2e7SValeria Barra 
957cb32e2e7SValeria Barra   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
958*9b072555Sjeremylt   ierr = VecDestroy(&rhs_loc); CHKERRQ(ierr);
959cb32e2e7SValeria Barra   ierr = VecDestroy(&X); CHKERRQ(ierr);
960*9b072555Sjeremylt   ierr = VecDestroy(&user->X_loc); CHKERRQ(ierr);
961*9b072555Sjeremylt   ierr = VecDestroy(&user->Y_loc); CHKERRQ(ierr);
962*9b072555Sjeremylt   ierr = VecScatterDestroy(&l_to_g); CHKERRQ(ierr);
963*9b072555Sjeremylt   ierr = VecScatterDestroy(&l_to_g_0); CHKERRQ(ierr);
964*9b072555Sjeremylt   ierr = VecScatterDestroy(&g_to_g_D); CHKERRQ(ierr);
965cb32e2e7SValeria Barra   ierr = MatDestroy(&mat); CHKERRQ(ierr);
966cb32e2e7SValeria Barra   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
967cb32e2e7SValeria Barra 
968*9b072555Sjeremylt   CeedVectorDestroy(&user->x_ceed);
969*9b072555Sjeremylt   CeedVectorDestroy(&user->y_ceed);
970*9b072555Sjeremylt   CeedVectorDestroy(&user->q_data);
971cb32e2e7SValeria Barra   CeedVectorDestroy(&target);
972*9b072555Sjeremylt   CeedOperatorDestroy(&op_setup_geo);
973*9b072555Sjeremylt   CeedOperatorDestroy(&op_setup_rhs);
974*9b072555Sjeremylt   CeedOperatorDestroy(&op_apply);
975*9b072555Sjeremylt   CeedOperatorDestroy(&op_error);
976*9b072555Sjeremylt   CeedElemRestrictionDestroy(&elem_restr_u);
977*9b072555Sjeremylt   CeedElemRestrictionDestroy(&elem_restr_x);
978*9b072555Sjeremylt   CeedElemRestrictionDestroy(&elem_restr_u_i);
979*9b072555Sjeremylt   CeedElemRestrictionDestroy(&elem_restr_qd_i);
980*9b072555Sjeremylt   CeedQFunctionDestroy(&qf_setup_geo);
981*9b072555Sjeremylt   CeedQFunctionDestroy(&qf_setup_rhs);
982*9b072555Sjeremylt   CeedQFunctionDestroy(&qf_apply);
983*9b072555Sjeremylt   CeedQFunctionDestroy(&qf_error);
984*9b072555Sjeremylt   CeedBasisDestroy(&basis_u);
985*9b072555Sjeremylt   CeedBasisDestroy(&basis_x);
986cb32e2e7SValeria Barra   CeedDestroy(&ceed);
987cb32e2e7SValeria Barra   ierr = PetscFree(user); CHKERRQ(ierr);
988cb32e2e7SValeria Barra   return PetscFinalize();
989cb32e2e7SValeria Barra }
990