xref: /libCEED/examples/petsc/bpsraw.c (revision edf0491998c1d524f2f70fdd683669b8904cb3b6)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 //                        libCEED + PETSc Example: CEED BPs
9 //
10 // This example demonstrates a simple usage of libCEED with PETSc to solve the CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps.
11 //
12 // The code is intentionally "raw", using only low-level communication primitives.
13 //
14 // Build with:
15 //
16 //     make bpsraw [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
17 //
18 // Sample runs:
19 //
20 //     ./bpsraw -problem bp1
21 //     ./bpsraw -problem bp2
22 //     ./bpsraw -problem bp3
23 //     ./bpsraw -problem bp4
24 //     ./bpsraw -problem bp5 -ceed /cpu/self
25 //     ./bpsraw -problem bp6 -ceed /gpu/cuda
26 //
27 //TESTARGS -ceed {ceed_resource} -test -problem bp2 -degree 5 -q_extra 1 -ksp_max_it_clip 15,15
28 
29 /// @file
30 /// CEED BPs example using PETSc
31 /// See bps.c for an implementation using DMPlex unstructured grids.
32 const char help[] = "Solve CEED BPs using PETSc\n";
33 
34 #include <ceed.h>
35 #include <petscdm.h>
36 #include <petscksp.h>
37 #include <stdbool.h>
38 #include <string.h>
39 
40 #include "qfunctions/bps/bp1.h"
41 #include "qfunctions/bps/bp2.h"
42 #include "qfunctions/bps/bp3.h"
43 #include "qfunctions/bps/bp4.h"
44 #include "qfunctions/bps/common.h"
45 
46 static CeedMemType MemTypeP2C(PetscMemType mem_type) { return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST; }
47 
48 static void Split3(PetscInt size, PetscInt m[3], bool reverse) {
49   for (PetscInt d = 0, size_left = size; d < 3; d++) {
50     PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(size_left, 1. / (3 - d)));
51     while (try * (size_left / try) != size_left) try++;
52     m[reverse ? 2 - d : d] = try;
53     size_left /= try;
54   }
55 }
56 
57 static PetscInt Max3(const PetscInt a[3]) { return PetscMax(a[0], PetscMax(a[1], a[2])); }
58 static PetscInt Min3(const PetscInt a[3]) { return PetscMin(a[0], PetscMin(a[1], a[2])); }
59 static void     GlobalNodes(const PetscInt p[3], const PetscInt i_rank[3], PetscInt degree, const PetscInt mesh_elem[3], PetscInt m_nodes[3]) {
60   for (int d = 0; d < 3; d++) m_nodes[d] = degree * mesh_elem[d] + (i_rank[d] == p[d] - 1);
61 }
62 static PetscInt GlobalStart(const PetscInt p[3], const PetscInt i_rank[3], PetscInt degree, const PetscInt mesh_elem[3]) {
63   PetscInt start = 0;
64   // Dumb brute-force is easier to read
65   for (PetscInt i = 0; i < p[0]; i++) {
66     for (PetscInt j = 0; j < p[1]; j++) {
67       for (PetscInt k = 0; k < p[2]; k++) {
68         PetscInt m_nodes[3], ijk_rank[] = {i, j, k};
69         if (i == i_rank[0] && j == i_rank[1] && k == i_rank[2]) return start;
70         GlobalNodes(p, ijk_rank, degree, mesh_elem, m_nodes);
71         start += m_nodes[0] * m_nodes[1] * m_nodes[2];
72       }
73     }
74   }
75   return -1;
76 }
77 static PetscErrorCode CreateRestriction(Ceed ceed, const CeedInt mesh_elem[3], CeedInt P, CeedInt num_comp, CeedElemRestriction *elem_restr) {
78   const PetscInt num_elem = mesh_elem[0] * mesh_elem[1] * mesh_elem[2];
79   PetscInt       m_nodes[3], *idx, *idx_p;
80 
81   PetscFunctionBeginUser;
82   // Get indicies
83   for (int d = 0; d < 3; d++) m_nodes[d] = mesh_elem[d] * (P - 1) + 1;
84   idx_p = idx = malloc(num_elem * P * P * P * sizeof idx[0]);
85   for (CeedInt i = 0; i < mesh_elem[0]; i++) {
86     for (CeedInt j = 0; j < mesh_elem[1]; j++) {
87       for (CeedInt k = 0; k < mesh_elem[2]; k++, idx_p += P * P * P) {
88         for (CeedInt ii = 0; ii < P; ii++) {
89           for (CeedInt jj = 0; jj < P; jj++) {
90             for (CeedInt kk = 0; kk < P; kk++) {
91               if (0) {  // This is the C-style (i,j,k) ordering that I prefer
92                 idx_p[(ii * P + jj) * P + kk] = num_comp * (((i * (P - 1) + ii) * m_nodes[1] + (j * (P - 1) + jj)) * m_nodes[2] + (k * (P - 1) + kk));
93               } else {  // (k,j,i) ordering for consistency with MFEM example
94                 idx_p[ii + P * (jj + P * kk)] = num_comp * (((i * (P - 1) + ii) * m_nodes[1] + (j * (P - 1) + jj)) * m_nodes[2] + (k * (P - 1) + kk));
95               }
96             }
97           }
98         }
99       }
100     }
101   }
102 
103   // Setup CEED restriction
104   CeedElemRestrictionCreate(ceed, num_elem, P * P * P, num_comp, 1, m_nodes[0] * m_nodes[1] * m_nodes[2] * num_comp, CEED_MEM_HOST, CEED_OWN_POINTER,
105                             idx, elem_restr);
106 
107   PetscFunctionReturn(PETSC_SUCCESS);
108 }
109 
110 // Data for PETSc
111 typedef struct OperatorApplyContext_ *OperatorApplyContext;
112 struct OperatorApplyContext_ {
113   MPI_Comm     comm;
114   VecScatter   l_to_g;    // Scatter for all entries
115   VecScatter   l_to_g_0;  // Skip Dirichlet values
116   VecScatter   g_to_g_D;  // global-to-global; only Dirichlet values
117   Vec          X_loc, Y_loc;
118   CeedVector   x_ceed, y_ceed;
119   CeedOperator op;
120   CeedVector   q_data;
121   Ceed         ceed;
122 };
123 
124 // BP Options
125 typedef enum { CEED_BP1 = 0, CEED_BP2 = 1, CEED_BP3 = 2, CEED_BP4 = 3, CEED_BP5 = 4, CEED_BP6 = 5 } BPType;
126 static const char *const bp_types[] = {"bp1", "bp2", "bp3", "bp4", "bp5", "bp6", "BPType", "CEED_BP", 0};
127 
128 // BP specific data
129 typedef struct {
130   CeedInt           num_comp_u, q_data_size, q_extra;
131   CeedQFunctionUser setup_geo, setup_rhs, apply, error;
132   const char       *setup_geo_loc, *setup_rhs_loc, *apply_loc, *error_loc;
133   CeedEvalMode      in_mode, out_mode;
134   CeedQuadMode      q_mode;
135 } BPData;
136 
137 BPData bp_options[6] = {
138     [CEED_BP1] = {.num_comp_u    = 1,
139                   .q_data_size   = 1,
140                   .q_extra       = 1,
141                   .setup_geo     = SetupMassGeo,
142                   .setup_rhs     = SetupMassRhs,
143                   .apply         = Mass,
144                   .error         = Error,
145                   .setup_geo_loc = SetupMassGeo_loc,
146                   .setup_rhs_loc = SetupMassRhs_loc,
147                   .apply_loc     = Mass_loc,
148                   .error_loc     = Error_loc,
149                   .in_mode       = CEED_EVAL_INTERP,
150                   .out_mode      = CEED_EVAL_INTERP,
151                   .q_mode        = CEED_GAUSS        },
152     [CEED_BP2] = {.num_comp_u    = 3,
153                   .q_data_size   = 1,
154                   .q_extra       = 1,
155                   .setup_geo     = SetupMassGeo,
156                   .setup_rhs     = SetupMassRhs3,
157                   .apply         = Mass3,
158                   .error         = Error3,
159                   .setup_geo_loc = SetupMassGeo_loc,
160                   .setup_rhs_loc = SetupMassRhs3_loc,
161                   .apply_loc     = Mass3_loc,
162                   .error_loc     = Error3_loc,
163                   .in_mode       = CEED_EVAL_INTERP,
164                   .out_mode      = CEED_EVAL_INTERP,
165                   .q_mode        = CEED_GAUSS        },
166     [CEED_BP3] = {.num_comp_u    = 1,
167                   .q_data_size   = 7,
168                   .q_extra       = 1,
169                   .setup_geo     = SetupDiffGeo,
170                   .setup_rhs     = SetupDiffRhs,
171                   .apply         = Diff,
172                   .error         = Error,
173                   .setup_geo_loc = SetupDiffGeo_loc,
174                   .setup_rhs_loc = SetupDiffRhs_loc,
175                   .apply_loc     = Diff_loc,
176                   .error_loc     = Error_loc,
177                   .in_mode       = CEED_EVAL_GRAD,
178                   .out_mode      = CEED_EVAL_GRAD,
179                   .q_mode        = CEED_GAUSS        },
180     [CEED_BP4] = {.num_comp_u    = 3,
181                   .q_data_size   = 7,
182                   .q_extra       = 1,
183                   .setup_geo     = SetupDiffGeo,
184                   .setup_rhs     = SetupDiffRhs3,
185                   .apply         = Diff3,
186                   .error         = Error3,
187                   .setup_geo_loc = SetupDiffGeo_loc,
188                   .setup_rhs_loc = SetupDiffRhs3_loc,
189                   .apply_loc     = Diff3_loc,
190                   .error_loc     = Error3_loc,
191                   .in_mode       = CEED_EVAL_GRAD,
192                   .out_mode      = CEED_EVAL_GRAD,
193                   .q_mode        = CEED_GAUSS        },
194     [CEED_BP5] = {.num_comp_u    = 1,
195                   .q_data_size   = 7,
196                   .q_extra       = 0,
197                   .setup_geo     = SetupDiffGeo,
198                   .setup_rhs     = SetupDiffRhs,
199                   .apply         = Diff,
200                   .error         = Error,
201                   .setup_geo_loc = SetupDiffGeo_loc,
202                   .setup_rhs_loc = SetupDiffRhs_loc,
203                   .apply_loc     = Diff_loc,
204                   .error_loc     = Error_loc,
205                   .in_mode       = CEED_EVAL_GRAD,
206                   .out_mode      = CEED_EVAL_GRAD,
207                   .q_mode        = CEED_GAUSS_LOBATTO},
208     [CEED_BP6] = {.num_comp_u    = 3,
209                   .q_data_size   = 7,
210                   .q_extra       = 0,
211                   .setup_geo     = SetupDiffGeo,
212                   .setup_rhs     = SetupDiffRhs3,
213                   .apply         = Diff3,
214                   .error         = Error3,
215                   .setup_geo_loc = SetupDiffGeo_loc,
216                   .setup_rhs_loc = SetupDiffRhs3_loc,
217                   .apply_loc     = Diff3_loc,
218                   .error_loc     = Error3_loc,
219                   .in_mode       = CEED_EVAL_GRAD,
220                   .out_mode      = CEED_EVAL_GRAD,
221                   .q_mode        = CEED_GAUSS_LOBATTO}
222 };
223 
224 // This function uses libCEED to compute the action of the mass matrix
225 static PetscErrorCode MatMult_Mass(Mat A, Vec X, Vec Y) {
226   OperatorApplyContext op_apply_ctx;
227   PetscScalar         *x, *y;
228   PetscMemType         x_mem_type, y_mem_type;
229 
230   PetscFunctionBeginUser;
231 
232   PetscCall(MatShellGetContext(A, &op_apply_ctx));
233 
234   // Global-to-local
235   PetscCall(VecScatterBegin(op_apply_ctx->l_to_g, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REVERSE));
236   PetscCall(VecScatterEnd(op_apply_ctx->l_to_g, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REVERSE));
237 
238   // Setup libCEED vectors
239   PetscCall(VecGetArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x, &x_mem_type));
240   PetscCall(VecGetArrayAndMemType(op_apply_ctx->Y_loc, &y, &y_mem_type));
241   CeedVectorSetArray(op_apply_ctx->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x);
242   CeedVectorSetArray(op_apply_ctx->y_ceed, MemTypeP2C(y_mem_type), CEED_USE_POINTER, y);
243 
244   // Apply libCEED operator
245   CeedOperatorApply(op_apply_ctx->op, op_apply_ctx->x_ceed, op_apply_ctx->y_ceed, CEED_REQUEST_IMMEDIATE);
246 
247   // Restore PETSc vectors
248   CeedVectorTakeArray(op_apply_ctx->x_ceed, MemTypeP2C(x_mem_type), NULL);
249   CeedVectorTakeArray(op_apply_ctx->y_ceed, MemTypeP2C(y_mem_type), NULL);
250   PetscCall(VecRestoreArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x));
251   PetscCall(VecRestoreArrayAndMemType(op_apply_ctx->Y_loc, &y));
252 
253   // Local-to-global
254   if (Y) {
255     PetscCall(VecZeroEntries(Y));
256     PetscCall(VecScatterBegin(op_apply_ctx->l_to_g, op_apply_ctx->Y_loc, Y, ADD_VALUES, SCATTER_FORWARD));
257     PetscCall(VecScatterEnd(op_apply_ctx->l_to_g, op_apply_ctx->Y_loc, Y, ADD_VALUES, SCATTER_FORWARD));
258   }
259   PetscFunctionReturn(PETSC_SUCCESS);
260 }
261 
262 // This function uses libCEED to compute the action of the Laplacian with Dirichlet boundary conditions
263 static PetscErrorCode MatMult_Diff(Mat A, Vec X, Vec Y) {
264   OperatorApplyContext op_apply_ctx;
265   PetscScalar         *x, *y;
266   PetscMemType         x_mem_type, y_mem_type;
267 
268   PetscFunctionBeginUser;
269 
270   PetscCall(MatShellGetContext(A, &op_apply_ctx));
271 
272   // Global-to-local
273   PetscCall(VecScatterBegin(op_apply_ctx->l_to_g_0, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REVERSE));
274   PetscCall(VecScatterEnd(op_apply_ctx->l_to_g_0, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REVERSE));
275 
276   // Setup libCEED vectors
277   PetscCall(VecGetArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x, &x_mem_type));
278   PetscCall(VecGetArrayAndMemType(op_apply_ctx->Y_loc, &y, &y_mem_type));
279   CeedVectorSetArray(op_apply_ctx->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x);
280   CeedVectorSetArray(op_apply_ctx->y_ceed, MemTypeP2C(y_mem_type), CEED_USE_POINTER, y);
281 
282   // Apply libCEED operator
283   CeedOperatorApply(op_apply_ctx->op, op_apply_ctx->x_ceed, op_apply_ctx->y_ceed, CEED_REQUEST_IMMEDIATE);
284 
285   // Restore PETSc vectors
286   CeedVectorTakeArray(op_apply_ctx->x_ceed, MemTypeP2C(x_mem_type), NULL);
287   CeedVectorTakeArray(op_apply_ctx->y_ceed, MemTypeP2C(y_mem_type), NULL);
288   PetscCall(VecRestoreArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x));
289   PetscCall(VecRestoreArrayAndMemType(op_apply_ctx->Y_loc, &y));
290 
291   // Local-to-global
292   PetscCall(VecZeroEntries(Y));
293   PetscCall(VecScatterBegin(op_apply_ctx->g_to_g_D, X, Y, INSERT_VALUES, SCATTER_FORWARD));
294   PetscCall(VecScatterEnd(op_apply_ctx->g_to_g_D, X, Y, INSERT_VALUES, SCATTER_FORWARD));
295   PetscCall(VecScatterBegin(op_apply_ctx->l_to_g_0, op_apply_ctx->Y_loc, Y, ADD_VALUES, SCATTER_FORWARD));
296   PetscCall(VecScatterEnd(op_apply_ctx->l_to_g_0, op_apply_ctx->Y_loc, Y, ADD_VALUES, SCATTER_FORWARD));
297 
298   PetscFunctionReturn(PETSC_SUCCESS);
299 }
300 
301 // This function calculates the error in the final solution
302 static PetscErrorCode ComputeErrorMax(OperatorApplyContext op_apply_ctx, CeedOperator op_error, Vec X, CeedVector target, PetscReal *max_error) {
303   PetscScalar *x;
304   PetscMemType mem_type;
305   CeedVector   collocated_error;
306   CeedSize     length;
307 
308   PetscFunctionBeginUser;
309 
310   CeedVectorGetLength(target, &length);
311   CeedVectorCreate(op_apply_ctx->ceed, length, &collocated_error);
312 
313   // Global-to-local
314   PetscCall(VecScatterBegin(op_apply_ctx->l_to_g, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REVERSE));
315   PetscCall(VecScatterEnd(op_apply_ctx->l_to_g, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REVERSE));
316 
317   // Setup libCEED vector
318   PetscCall(VecGetArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x, &mem_type));
319   CeedVectorSetArray(op_apply_ctx->x_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, x);
320 
321   // Apply libCEED operator
322   CeedOperatorApply(op_error, op_apply_ctx->x_ceed, collocated_error, CEED_REQUEST_IMMEDIATE);
323 
324   // Restore PETSc vector
325   CeedVectorTakeArray(op_apply_ctx->x_ceed, MemTypeP2C(mem_type), NULL);
326   PetscCall(VecRestoreArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x));
327 
328   // Reduce max error
329   *max_error = 0;
330   const CeedScalar *e;
331   CeedVectorGetArrayRead(collocated_error, CEED_MEM_HOST, &e);
332   for (CeedInt i = 0; i < length; i++) {
333     *max_error = PetscMax(*max_error, PetscAbsScalar(e[i]));
334   }
335   CeedVectorRestoreArrayRead(collocated_error, &e);
336   PetscCall(MPI_Allreduce(MPI_IN_PLACE, max_error, 1, MPIU_REAL, MPIU_MAX, op_apply_ctx->comm));
337 
338   // Cleanup
339   CeedVectorDestroy(&collocated_error);
340 
341   PetscFunctionReturn(PETSC_SUCCESS);
342 }
343 
344 int main(int argc, char **argv) {
345   MPI_Comm comm;
346   char     ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self";
347   double   my_rt_start, my_rt, rt_min, rt_max;
348   PetscInt degree, q_extra, local_nodes, local_elem, mesh_elem[3], m_nodes[3], p[3], i_rank[3], l_nodes[3], l_size, num_comp_u = 1,
349                                                                                                                     ksp_max_it_clip[2];
350   PetscScalar         *r;
351   PetscBool            test_mode, benchmark_mode, write_solution;
352   PetscMPIInt          size, rank;
353   PetscLogStage        solve_stage;
354   Vec                  X, X_loc, rhs, rhs_loc;
355   Mat                  mat;
356   KSP                  ksp;
357   VecScatter           l_to_g, l_to_g_0, g_to_g_D;
358   PetscMemType         mem_type;
359   OperatorApplyContext op_apply_ctx;
360   Ceed                 ceed;
361   CeedBasis            basis_x, basis_u;
362   CeedElemRestriction  elem_restr_x, elem_restr_u, elem_restr_u_i, elem_restr_qd_i;
363   CeedQFunction        qf_setup_geo, qf_setup_rhs, qf_apply, qf_error;
364   CeedOperator         op_setup_geo, op_setup_rhs, op_apply, op_error;
365   CeedVector           x_coord, q_data, rhs_ceed, target;
366   CeedInt              P, Q;
367   const CeedInt        dim = 3, num_comp_x = 3;
368   BPType               bp_choice;
369 
370   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
371   comm = PETSC_COMM_WORLD;
372 
373   // Read command line options
374   PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
375   bp_choice = CEED_BP1;
376   PetscCall(PetscOptionsEnum("-problem", "CEED benchmark problem to solve", NULL, bp_types, (PetscEnum)bp_choice, (PetscEnum *)&bp_choice, NULL));
377   num_comp_u = bp_options[bp_choice].num_comp_u;
378   test_mode  = PETSC_FALSE;
379   PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, test_mode, &test_mode, NULL));
380   benchmark_mode = PETSC_FALSE;
381   PetscCall(PetscOptionsBool("-benchmark", "Benchmarking mode (prints benchmark statistics)", NULL, benchmark_mode, &benchmark_mode, NULL));
382   write_solution = PETSC_FALSE;
383   PetscCall(PetscOptionsBool("-write_solution", "Write solution for visualization", NULL, write_solution, &write_solution, NULL));
384   degree = test_mode ? 3 : 1;
385   PetscCall(PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", NULL, degree, &degree, NULL));
386   q_extra = bp_options[bp_choice].q_extra;
387   PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, q_extra, &q_extra, NULL));
388   PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, ceed_resource, ceed_resource, sizeof(ceed_resource), NULL));
389   local_nodes = 1000;
390   PetscCall(PetscOptionsInt("-local", "Target number of locally owned nodes per process", NULL, local_nodes, &local_nodes, NULL));
391   PetscInt two       = 2;
392   ksp_max_it_clip[0] = 5;
393   ksp_max_it_clip[1] = 20;
394   PetscCall(
395       PetscOptionsIntArray("-ksp_max_it_clip", "Min and max number of iterations to use during benchmarking", NULL, ksp_max_it_clip, &two, NULL));
396   PetscOptionsEnd();
397   P = degree + 1;
398   Q = P + q_extra;
399 
400   // Set up libCEED
401   CeedInit(ceed_resource, &ceed);
402   CeedMemType mem_type_backend;
403   CeedGetPreferredMemType(ceed, &mem_type_backend);
404 
405   VecType default_vec_type = NULL, vec_type;
406   switch (mem_type_backend) {
407     case CEED_MEM_HOST:
408       default_vec_type = VECSTANDARD;
409       break;
410     case CEED_MEM_DEVICE: {
411       const char *resolved;
412       CeedGetResource(ceed, &resolved);
413       if (strstr(resolved, "/gpu/cuda")) default_vec_type = VECCUDA;
414       else if (strstr(resolved, "/gpu/hip/occa")) default_vec_type = VECSTANDARD;  // https://github.com/CEED/libCEED/issues/678
415       else if (strstr(resolved, "/gpu/hip")) default_vec_type = VECHIP;
416       else default_vec_type = VECSTANDARD;
417     }
418   }
419 
420   // Determine size of process grid
421   PetscCall(MPI_Comm_size(comm, &size));
422   Split3(size, p, false);
423 
424   // Find a nicely composite number of elements no less than local_nodes
425   for (local_elem = PetscMax(1, local_nodes / (degree * degree * degree));; local_elem++) {
426     Split3(local_elem, mesh_elem, true);
427     if (Max3(mesh_elem) / Min3(mesh_elem) <= 2) break;
428   }
429 
430   // Find my location in the process grid
431   PetscCall(MPI_Comm_rank(comm, &rank));
432   for (int d = 0, rank_left = rank; d < dim; d++) {
433     const int pstride[3] = {p[1] * p[2], p[2], 1};
434     i_rank[d]            = rank_left / pstride[d];
435     rank_left -= i_rank[d] * pstride[d];
436   }
437 
438   GlobalNodes(p, i_rank, degree, mesh_elem, m_nodes);
439 
440   // Setup global vector
441   PetscCall(VecCreate(comm, &X));
442   PetscCall(VecSetType(X, default_vec_type));
443   PetscCall(VecSetSizes(X, m_nodes[0] * m_nodes[1] * m_nodes[2] * num_comp_u, PETSC_DECIDE));
444   PetscCall(VecSetFromOptions(X));
445   PetscCall(VecSetUp(X));
446 
447   // Set up libCEED
448   CeedInit(ceed_resource, &ceed);
449 
450   // Print summary
451   CeedInt gsize;
452   PetscCall(VecGetSize(X, &gsize));
453   if (!test_mode) {
454     const char *used_resource;
455     CeedGetResource(ceed, &used_resource);
456 
457     PetscCall(VecGetType(X, &vec_type));
458 
459     PetscCall(PetscPrintf(comm,
460                           "\n-- CEED Benchmark Problem %" CeedInt_FMT " -- libCEED + PETSc --\n"
461                           "  PETSc:\n"
462                           "    PETSc Vec Type                     : %s\n"
463                           "  libCEED:\n"
464                           "    libCEED Backend                    : %s\n"
465                           "    libCEED Backend MemType            : %s\n"
466                           "  Mesh:\n"
467                           "    Solution Order (P)                 : %" CeedInt_FMT "\n"
468                           "    Quadrature  Order (Q)              : %" CeedInt_FMT "\n"
469                           "    Global nodes                       : %" PetscInt_FMT "\n"
470                           "    Process Decomposition              : %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n"
471                           "    Local Elements                     : %" PetscInt_FMT " = %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n"
472                           "    Owned nodes                        : %" PetscInt_FMT " = %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n"
473                           "    DoF per node                       : %" PetscInt_FMT "\n",
474                           bp_choice + 1, vec_type, used_resource, CeedMemTypes[mem_type_backend], P, Q, gsize / num_comp_u, p[0], p[1], p[2],
475                           local_elem, mesh_elem[0], mesh_elem[1], mesh_elem[2], m_nodes[0] * m_nodes[1] * m_nodes[2], m_nodes[0], m_nodes[1],
476                           m_nodes[2], num_comp_u));
477   }
478 
479   {
480     l_size = 1;
481     for (int d = 0; d < dim; d++) {
482       l_nodes[d] = mesh_elem[d] * degree + 1;
483       l_size *= l_nodes[d];
484     }
485     PetscCall(VecCreate(PETSC_COMM_SELF, &X_loc));
486     PetscCall(VecSetType(X_loc, default_vec_type));
487     PetscCall(VecSetSizes(X_loc, l_size * num_comp_u, PETSC_DECIDE));
488     PetscCall(VecSetFromOptions(X_loc));
489     PetscCall(VecSetUp(X_loc));
490 
491     // Create local-to-global scatter
492     PetscInt *l_to_g_ind, *l_to_g_ind_0, *loc_ind, l_0_count;
493     IS        l_to_g_is, l_to_g_is_0, loc_is;
494     PetscInt  g_start[2][2][2], g_m_nodes[2][2][2][dim];
495 
496     for (int i = 0; i < 2; i++) {
497       for (int j = 0; j < 2; j++) {
498         for (int k = 0; k < 2; k++) {
499           PetscInt ijk_rank[3] = {i_rank[0] + i, i_rank[1] + j, i_rank[2] + k};
500           g_start[i][j][k]     = GlobalStart(p, ijk_rank, degree, mesh_elem);
501           GlobalNodes(p, ijk_rank, degree, mesh_elem, g_m_nodes[i][j][k]);
502         }
503       }
504     }
505 
506     PetscCall(PetscMalloc1(l_size, &l_to_g_ind));
507     PetscCall(PetscMalloc1(l_size, &l_to_g_ind_0));
508     PetscCall(PetscMalloc1(l_size, &loc_ind));
509     l_0_count = 0;
510     for (PetscInt i = 0, ir, ii; ir = i >= m_nodes[0], ii = i - ir * m_nodes[0], i < l_nodes[0]; i++) {
511       for (PetscInt j = 0, jr, jj; jr = j >= m_nodes[1], jj = j - jr * m_nodes[1], j < l_nodes[1]; j++) {
512         for (PetscInt k = 0, kr, kk; kr = k >= m_nodes[2], kk = k - kr * m_nodes[2], k < l_nodes[2]; k++) {
513           PetscInt here    = (i * l_nodes[1] + j) * l_nodes[2] + k;
514           l_to_g_ind[here] = g_start[ir][jr][kr] + (ii * g_m_nodes[ir][jr][kr][1] + jj) * g_m_nodes[ir][jr][kr][2] + kk;
515           if ((i_rank[0] == 0 && i == 0) || (i_rank[1] == 0 && j == 0) || (i_rank[2] == 0 && k == 0) ||
516               (i_rank[0] + 1 == p[0] && i + 1 == l_nodes[0]) || (i_rank[1] + 1 == p[1] && j + 1 == l_nodes[1]) ||
517               (i_rank[2] + 1 == p[2] && k + 1 == l_nodes[2]))
518             continue;
519           l_to_g_ind_0[l_0_count] = l_to_g_ind[here];
520           loc_ind[l_0_count++]    = here;
521         }
522       }
523     }
524     PetscCall(ISCreateBlock(comm, num_comp_u, l_size, l_to_g_ind, PETSC_OWN_POINTER, &l_to_g_is));
525     PetscCall(VecScatterCreate(X_loc, NULL, X, l_to_g_is, &l_to_g));
526     PetscCall(ISCreateBlock(comm, num_comp_u, l_0_count, l_to_g_ind_0, PETSC_OWN_POINTER, &l_to_g_is_0));
527     PetscCall(ISCreateBlock(comm, num_comp_u, l_0_count, loc_ind, PETSC_OWN_POINTER, &loc_is));
528     PetscCall(VecScatterCreate(X_loc, loc_is, X, l_to_g_is_0, &l_to_g_0));
529     {
530       // Create global-to-global scatter for Dirichlet values (everything not in l_to_g_is_0, which is the range of l_to_g_0)
531       PetscInt           x_start, x_end, *ind_D, count_D = 0;
532       IS                 is_D;
533       const PetscScalar *x;
534       PetscCall(VecZeroEntries(X_loc));
535       PetscCall(VecSet(X, 1.0));
536       PetscCall(VecScatterBegin(l_to_g_0, X_loc, X, INSERT_VALUES, SCATTER_FORWARD));
537       PetscCall(VecScatterEnd(l_to_g_0, X_loc, X, INSERT_VALUES, SCATTER_FORWARD));
538       PetscCall(VecGetOwnershipRange(X, &x_start, &x_end));
539       PetscCall(PetscMalloc1(x_end - x_start, &ind_D));
540       PetscCall(VecGetArrayRead(X, &x));
541       for (PetscInt i = 0; i < x_end - x_start; i++) {
542         if (x[i] == 1.) ind_D[count_D++] = x_start + i;
543       }
544       PetscCall(VecRestoreArrayRead(X, &x));
545       PetscCall(ISCreateGeneral(comm, count_D, ind_D, PETSC_COPY_VALUES, &is_D));
546       PetscCall(PetscFree(ind_D));
547       PetscCall(VecScatterCreate(X, is_D, X, is_D, &g_to_g_D));
548       PetscCall(ISDestroy(&is_D));
549     }
550     PetscCall(ISDestroy(&l_to_g_is));
551     PetscCall(ISDestroy(&l_to_g_is_0));
552     PetscCall(ISDestroy(&loc_is));
553   }
554 
555   // CEED bases
556   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_u, P, Q, bp_options[bp_choice].q_mode, &basis_u);
557   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, bp_options[bp_choice].q_mode, &basis_x);
558 
559   // CEED restrictions
560   CreateRestriction(ceed, mesh_elem, P, num_comp_u, &elem_restr_u);
561   CreateRestriction(ceed, mesh_elem, 2, dim, &elem_restr_x);
562   CeedInt num_elem = mesh_elem[0] * mesh_elem[1] * mesh_elem[2];
563   CeedElemRestrictionCreateStrided(ceed, num_elem, Q * Q * Q, num_comp_u, num_comp_u * num_elem * Q * Q * Q, CEED_STRIDES_BACKEND, &elem_restr_u_i);
564   CeedElemRestrictionCreateStrided(ceed, num_elem, Q * Q * Q, bp_options[bp_choice].q_data_size,
565                                    bp_options[bp_choice].q_data_size * num_elem * Q * Q * Q, CEED_STRIDES_BACKEND, &elem_restr_qd_i);
566   {
567     CeedScalar *x_loc;
568     CeedInt     shape[3] = {mesh_elem[0] + 1, mesh_elem[1] + 1, mesh_elem[2] + 1}, len = shape[0] * shape[1] * shape[2];
569     x_loc = malloc(len * num_comp_x * sizeof x_loc[0]);
570     for (CeedInt i = 0; i < shape[0]; i++) {
571       for (CeedInt j = 0; j < shape[1]; j++) {
572         for (CeedInt k = 0; k < shape[2]; k++) {
573           x_loc[dim * ((i * shape[1] + j) * shape[2] + k) + 0] = 1. * (i_rank[0] * mesh_elem[0] + i) / (p[0] * mesh_elem[0]);
574           x_loc[dim * ((i * shape[1] + j) * shape[2] + k) + 1] = 1. * (i_rank[1] * mesh_elem[1] + j) / (p[1] * mesh_elem[1]);
575           x_loc[dim * ((i * shape[1] + j) * shape[2] + k) + 2] = 1. * (i_rank[2] * mesh_elem[2] + k) / (p[2] * mesh_elem[2]);
576         }
577       }
578     }
579     CeedVectorCreate(ceed, len * num_comp_x, &x_coord);
580     CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_OWN_POINTER, x_loc);
581   }
582 
583   // Create the QFunction that builds the operator quadrature data
584   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].setup_geo, bp_options[bp_choice].setup_geo_loc, &qf_setup_geo);
585   CeedQFunctionAddInput(qf_setup_geo, "x", num_comp_x, CEED_EVAL_INTERP);
586   CeedQFunctionAddInput(qf_setup_geo, "dx", num_comp_x * dim, CEED_EVAL_GRAD);
587   CeedQFunctionAddInput(qf_setup_geo, "weight", 1, CEED_EVAL_WEIGHT);
588   CeedQFunctionAddOutput(qf_setup_geo, "q_data", bp_options[bp_choice].q_data_size, CEED_EVAL_NONE);
589 
590   // Create the QFunction that sets up the RHS and true solution
591   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].setup_rhs, bp_options[bp_choice].setup_rhs_loc, &qf_setup_rhs);
592   CeedQFunctionAddInput(qf_setup_rhs, "x", num_comp_x, CEED_EVAL_INTERP);
593   CeedQFunctionAddInput(qf_setup_rhs, "q_data", bp_options[bp_choice].q_data_size, CEED_EVAL_NONE);
594   CeedQFunctionAddOutput(qf_setup_rhs, "true_soln", num_comp_u, CEED_EVAL_NONE);
595   CeedQFunctionAddOutput(qf_setup_rhs, "rhs", num_comp_u, CEED_EVAL_INTERP);
596 
597   // Set up PDE operator
598   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].apply, bp_options[bp_choice].apply_loc, &qf_apply);
599   // Add inputs and outputs
600   CeedInt in_scale  = bp_options[bp_choice].in_mode == CEED_EVAL_GRAD ? 3 : 1;
601   CeedInt out_scale = bp_options[bp_choice].out_mode == CEED_EVAL_GRAD ? 3 : 1;
602   CeedQFunctionAddInput(qf_apply, "u", num_comp_u * in_scale, bp_options[bp_choice].in_mode);
603   CeedQFunctionAddInput(qf_apply, "q_data", bp_options[bp_choice].q_data_size, CEED_EVAL_NONE);
604   CeedQFunctionAddOutput(qf_apply, "v", num_comp_u * out_scale, bp_options[bp_choice].out_mode);
605 
606   // Create the error qfunction
607   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].error, bp_options[bp_choice].error_loc, &qf_error);
608   CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP);
609   CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE);
610   CeedQFunctionAddInput(qf_error, "qdata", bp_options[bp_choice].q_data_size, CEED_EVAL_NONE);
611   CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_NONE);
612 
613   // Create the persistent vectors that will be needed in setup
614   CeedInt num_qpts;
615   CeedBasisGetNumQuadraturePoints(basis_u, &num_qpts);
616   CeedVectorCreate(ceed, bp_options[bp_choice].q_data_size * num_elem * num_qpts, &q_data);
617   CeedVectorCreate(ceed, num_elem * num_qpts * num_comp_u, &target);
618   CeedVectorCreate(ceed, l_size * num_comp_u, &rhs_ceed);
619 
620   // Create the operator that builds the quadrature data for the ceed operator
621   CeedOperatorCreate(ceed, qf_setup_geo, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setup_geo);
622   CeedOperatorSetField(op_setup_geo, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
623   CeedOperatorSetField(op_setup_geo, "dx", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
624   CeedOperatorSetField(op_setup_geo, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE);
625   CeedOperatorSetField(op_setup_geo, "q_data", elem_restr_qd_i, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
626 
627   // Create the operator that builds the RHS and true solution
628   CeedOperatorCreate(ceed, qf_setup_rhs, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setup_rhs);
629   CeedOperatorSetField(op_setup_rhs, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
630   CeedOperatorSetField(op_setup_rhs, "q_data", elem_restr_qd_i, CEED_BASIS_COLLOCATED, q_data);
631   CeedOperatorSetField(op_setup_rhs, "true_soln", elem_restr_u_i, CEED_BASIS_COLLOCATED, target);
632   CeedOperatorSetField(op_setup_rhs, "rhs", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
633 
634   // Create the mass or diff operator
635   CeedOperatorCreate(ceed, qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_apply);
636   CeedOperatorSetField(op_apply, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
637   CeedOperatorSetField(op_apply, "q_data", elem_restr_qd_i, CEED_BASIS_COLLOCATED, q_data);
638   CeedOperatorSetField(op_apply, "v", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
639 
640   // Create the error operator
641   CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_error);
642   CeedOperatorSetField(op_error, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
643   CeedOperatorSetField(op_error, "true_soln", elem_restr_u_i, CEED_BASIS_COLLOCATED, target);
644   CeedOperatorSetField(op_error, "qdata", elem_restr_qd_i, CEED_BASIS_COLLOCATED, q_data);
645   CeedOperatorSetField(op_error, "error", elem_restr_u_i, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
646 
647   // Set up Mat
648   PetscCall(PetscMalloc1(1, &op_apply_ctx));
649   op_apply_ctx->comm   = comm;
650   op_apply_ctx->l_to_g = l_to_g;
651   if (bp_choice != CEED_BP1 && bp_choice != CEED_BP2) {
652     op_apply_ctx->l_to_g_0 = l_to_g_0;
653     op_apply_ctx->g_to_g_D = g_to_g_D;
654   }
655   op_apply_ctx->X_loc = X_loc;
656   PetscCall(VecDuplicate(X_loc, &op_apply_ctx->Y_loc));
657   CeedVectorCreate(ceed, l_size * num_comp_u, &op_apply_ctx->x_ceed);
658   CeedVectorCreate(ceed, l_size * num_comp_u, &op_apply_ctx->y_ceed);
659   op_apply_ctx->op     = op_apply;
660   op_apply_ctx->q_data = q_data;
661   op_apply_ctx->ceed   = ceed;
662 
663   PetscCall(MatCreateShell(comm, m_nodes[0] * m_nodes[1] * m_nodes[2] * num_comp_u, m_nodes[0] * m_nodes[1] * m_nodes[2] * num_comp_u, PETSC_DECIDE,
664                            PETSC_DECIDE, op_apply_ctx, &mat));
665   if (bp_choice == CEED_BP1 || bp_choice == CEED_BP2) {
666     PetscCall(MatShellSetOperation(mat, MATOP_MULT, (void (*)(void))MatMult_Mass));
667   } else {
668     PetscCall(MatShellSetOperation(mat, MATOP_MULT, (void (*)(void))MatMult_Diff));
669   }
670   PetscCall(VecGetType(X, &vec_type));
671   PetscCall(MatShellSetVecType(mat, vec_type));
672 
673   // Get RHS vector
674   PetscCall(VecDuplicate(X, &rhs));
675   PetscCall(VecDuplicate(X_loc, &rhs_loc));
676   PetscCall(VecZeroEntries(rhs_loc));
677   PetscCall(VecGetArrayAndMemType(rhs_loc, &r, &mem_type));
678   CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r);
679 
680   // Setup q_data, rhs, and target
681   CeedOperatorApply(op_setup_geo, x_coord, q_data, CEED_REQUEST_IMMEDIATE);
682   CeedOperatorApply(op_setup_rhs, x_coord, rhs_ceed, CEED_REQUEST_IMMEDIATE);
683   CeedVectorDestroy(&x_coord);
684 
685   // Gather RHS
686   PetscCall(CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL));
687   PetscCall(VecRestoreArrayAndMemType(rhs_loc, &r));
688   PetscCall(VecZeroEntries(rhs));
689   PetscCall(VecScatterBegin(l_to_g, rhs_loc, rhs, ADD_VALUES, SCATTER_FORWARD));
690   PetscCall(VecScatterEnd(l_to_g, rhs_loc, rhs, ADD_VALUES, SCATTER_FORWARD));
691   CeedVectorDestroy(&rhs_ceed);
692 
693   PetscCall(KSPCreate(comm, &ksp));
694   {
695     PC pc;
696     PetscCall(KSPGetPC(ksp, &pc));
697     if (bp_choice == CEED_BP1 || bp_choice == CEED_BP2) {
698       PetscCall(PCSetType(pc, PCJACOBI));
699       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
700     } else {
701       PetscCall(PCSetType(pc, PCNONE));
702     }
703     PetscCall(KSPSetType(ksp, KSPCG));
704     PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
705     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
706   }
707   PetscCall(KSPSetOperators(ksp, mat, mat));
708   // First run's performance log is not considered for benchmarking purposes
709   PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1));
710   my_rt_start = MPI_Wtime();
711   PetscCall(KSPSolve(ksp, rhs, X));
712   my_rt = MPI_Wtime() - my_rt_start;
713   PetscCall(MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm));
714   // Set maxits based on first iteration timing
715   if (my_rt > 0.02) {
716     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, ksp_max_it_clip[0]));
717   } else {
718     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, ksp_max_it_clip[1]));
719   }
720   PetscCall(KSPSetFromOptions(ksp));
721 
722   // Timed solve
723   PetscCall(VecZeroEntries(X));
724   PetscCall(PetscBarrier((PetscObject)ksp));
725 
726   // -- Performance logging
727   PetscCall(PetscLogStageRegister("Solve Stage", &solve_stage));
728   PetscCall(PetscLogStagePush(solve_stage));
729 
730   // -- Solve
731   my_rt_start = MPI_Wtime();
732   PetscCall(KSPSolve(ksp, rhs, X));
733   my_rt = MPI_Wtime() - my_rt_start;
734 
735   // -- Performance logging
736   PetscCall(PetscLogStagePop());
737 
738   // Output results
739   {
740     KSPType            ksp_type;
741     KSPConvergedReason reason;
742     PetscReal          rnorm;
743     PetscInt           its;
744     PetscCall(KSPGetType(ksp, &ksp_type));
745     PetscCall(KSPGetConvergedReason(ksp, &reason));
746     PetscCall(KSPGetIterationNumber(ksp, &its));
747     PetscCall(KSPGetResidualNorm(ksp, &rnorm));
748     if (!test_mode || reason < 0 || rnorm > 1e-8) {
749       PetscCall(PetscPrintf(comm,
750                             "  KSP:\n"
751                             "    KSP Type                           : %s\n"
752                             "    KSP Convergence                    : %s\n"
753                             "    Total KSP Iterations               : %" PetscInt_FMT "\n"
754                             "    Final rnorm                        : %e\n",
755                             ksp_type, KSPConvergedReasons[reason], its, (double)rnorm));
756     }
757     if (!test_mode) {
758       PetscCall(PetscPrintf(comm, "  Performance:\n"));
759     }
760     {
761       PetscReal max_error;
762       PetscCall(ComputeErrorMax(op_apply_ctx, op_error, X, target, &max_error));
763       PetscReal tol = 5e-2;
764       if (!test_mode || max_error > tol) {
765         PetscCall(MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm));
766         PetscCall(MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm));
767         PetscCall(PetscPrintf(comm,
768                               "    Pointwise Error (max)              : %e\n"
769                               "    CG Solve Time                      : %g (%g) sec\n",
770                               (double)max_error, rt_max, rt_min));
771       }
772     }
773     if (!test_mode) {
774       PetscCall(
775           PetscPrintf(comm, "    DoFs/Sec in CG                     : %g (%g) million\n", 1e-6 * gsize * its / rt_max, 1e-6 * gsize * its / rt_min));
776     }
777   }
778 
779   if (write_solution) {
780     PetscViewer vtk_viewer_soln;
781 
782     PetscCall(PetscViewerCreate(comm, &vtk_viewer_soln));
783     PetscCall(PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK));
784     PetscCall(PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"));
785     PetscCall(VecView(X, vtk_viewer_soln));
786     PetscCall(PetscViewerDestroy(&vtk_viewer_soln));
787   }
788 
789   PetscCall(VecDestroy(&rhs));
790   PetscCall(VecDestroy(&rhs_loc));
791   PetscCall(VecDestroy(&X));
792   PetscCall(VecDestroy(&op_apply_ctx->X_loc));
793   PetscCall(VecDestroy(&op_apply_ctx->Y_loc));
794   PetscCall(VecScatterDestroy(&l_to_g));
795   PetscCall(VecScatterDestroy(&l_to_g_0));
796   PetscCall(VecScatterDestroy(&g_to_g_D));
797   PetscCall(MatDestroy(&mat));
798   PetscCall(KSPDestroy(&ksp));
799 
800   CeedVectorDestroy(&op_apply_ctx->x_ceed);
801   CeedVectorDestroy(&op_apply_ctx->y_ceed);
802   CeedVectorDestroy(&op_apply_ctx->q_data);
803   CeedVectorDestroy(&target);
804   CeedOperatorDestroy(&op_setup_geo);
805   CeedOperatorDestroy(&op_setup_rhs);
806   CeedOperatorDestroy(&op_apply);
807   CeedOperatorDestroy(&op_error);
808   CeedElemRestrictionDestroy(&elem_restr_u);
809   CeedElemRestrictionDestroy(&elem_restr_x);
810   CeedElemRestrictionDestroy(&elem_restr_u_i);
811   CeedElemRestrictionDestroy(&elem_restr_qd_i);
812   CeedQFunctionDestroy(&qf_setup_geo);
813   CeedQFunctionDestroy(&qf_setup_rhs);
814   CeedQFunctionDestroy(&qf_apply);
815   CeedQFunctionDestroy(&qf_error);
816   CeedBasisDestroy(&basis_u);
817   CeedBasisDestroy(&basis_x);
818   CeedDestroy(&ceed);
819   PetscCall(PetscFree(op_apply_ctx));
820   return PetscFinalize();
821 }
822