xref: /libCEED/examples/petsc/bpsraw.c (revision ea61e9ac44808524e4667c1525a05976f536c19c)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3cb32e2e7SValeria Barra //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5cb32e2e7SValeria Barra //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7cb32e2e7SValeria Barra 
8cb32e2e7SValeria Barra //                        libCEED + PETSc Example: CEED BPs
9cb32e2e7SValeria Barra //
10*ea61e9acSJeremy L Thompson // This example demonstrates a simple usage of libCEED with PETSc to solve the CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps.
11cb32e2e7SValeria Barra //
12*ea61e9acSJeremy L Thompson // The code is intentionally "raw", using only low-level communication primitives.
13cb32e2e7SValeria Barra //
14cb32e2e7SValeria Barra // Build with:
15cb32e2e7SValeria Barra //
16cb32e2e7SValeria Barra //     make bpsraw [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
17cb32e2e7SValeria Barra //
18cb32e2e7SValeria Barra // Sample runs:
19cb32e2e7SValeria Barra //
2032d2ee49SValeria Barra //     ./bpsraw -problem bp1
2128688798Sjeremylt //     ./bpsraw -problem bp2
2228688798Sjeremylt //     ./bpsraw -problem bp3
2328688798Sjeremylt //     ./bpsraw -problem bp4
2428688798Sjeremylt //     ./bpsraw -problem bp5 -ceed /cpu/self
2528688798Sjeremylt //     ./bpsraw -problem bp6 -ceed /gpu/cuda
26cb32e2e7SValeria Barra //
279b072555Sjeremylt //TESTARGS -ceed {ceed_resource} -test -problem bp2 -degree 5 -q_extra 1 -ksp_max_it_clip 15,15
28cb32e2e7SValeria Barra 
29cb32e2e7SValeria Barra /// @file
30cb32e2e7SValeria Barra /// CEED BPs example using PETSc
31cb32e2e7SValeria Barra /// See bps.c for an implementation using DMPlex unstructured grids.
32cb32e2e7SValeria Barra const char help[] = "Solve CEED BPs using PETSc\n";
33cb32e2e7SValeria Barra 
343d576824SJeremy L Thompson #include <ceed.h>
353d576824SJeremy L Thompson #include <petscksp.h>
363d576824SJeremy L Thompson #include <petscsys.h>
37cb32e2e7SValeria Barra #include <stdbool.h>
38cb32e2e7SValeria Barra #include <string.h>
392b730f8bSJeremy L Thompson 
40cb32e2e7SValeria Barra #include "qfunctions/bps/bp1.h"
41cb32e2e7SValeria Barra #include "qfunctions/bps/bp2.h"
42cb32e2e7SValeria Barra #include "qfunctions/bps/bp3.h"
43cb32e2e7SValeria Barra #include "qfunctions/bps/bp4.h"
443d576824SJeremy L Thompson #include "qfunctions/bps/common.h"
45cb32e2e7SValeria Barra 
469396343dSjeremylt #if PETSC_VERSION_LT(3, 12, 0)
479396343dSjeremylt #ifdef PETSC_HAVE_CUDA
489396343dSjeremylt #include <petsccuda.h>
49*ea61e9acSJeremy L Thompson // Note: With PETSc prior to version 3.12.0, providing the source path to include 'cublas_v2.h' will be needed to use 'petsccuda.h'.
509396343dSjeremylt #endif
519396343dSjeremylt #endif
529396343dSjeremylt 
532b730f8bSJeremy L Thompson static CeedMemType MemTypeP2C(PetscMemType mem_type) { return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST; }
54b68a8d79SJed Brown 
55cb32e2e7SValeria Barra static void Split3(PetscInt size, PetscInt m[3], bool reverse) {
569b072555Sjeremylt   for (PetscInt d = 0, size_left = size; d < 3; d++) {
579b072555Sjeremylt     PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(size_left, 1. / (3 - d)));
589b072555Sjeremylt     while (try * (size_left / try) != size_left) try++;
59cb32e2e7SValeria Barra     m[reverse ? 2 - d : d] = try;
609b072555Sjeremylt     size_left /= try;
61cb32e2e7SValeria Barra   }
62cb32e2e7SValeria Barra }
63cb32e2e7SValeria Barra 
642b730f8bSJeremy L Thompson static PetscInt Max3(const PetscInt a[3]) { return PetscMax(a[0], PetscMax(a[1], a[2])); }
652b730f8bSJeremy L Thompson static PetscInt Min3(const PetscInt a[3]) { return PetscMin(a[0], PetscMin(a[1], a[2])); }
662b730f8bSJeremy L Thompson static void     GlobalNodes(const PetscInt p[3], const PetscInt i_rank[3], PetscInt degree, const PetscInt mesh_elem[3], PetscInt m_nodes[3]) {
672b730f8bSJeremy L Thompson       for (int d = 0; d < 3; d++) m_nodes[d] = degree * mesh_elem[d] + (i_rank[d] == p[d] - 1);
68cb32e2e7SValeria Barra }
692b730f8bSJeremy L Thompson static PetscInt GlobalStart(const PetscInt p[3], const PetscInt i_rank[3], PetscInt degree, const PetscInt mesh_elem[3]) {
70cb32e2e7SValeria Barra   PetscInt start = 0;
71cb32e2e7SValeria Barra   // Dumb brute-force is easier to read
72cb32e2e7SValeria Barra   for (PetscInt i = 0; i < p[0]; i++) {
73cb32e2e7SValeria Barra     for (PetscInt j = 0; j < p[1]; j++) {
74cb32e2e7SValeria Barra       for (PetscInt k = 0; k < p[2]; k++) {
759b072555Sjeremylt         PetscInt m_nodes[3], ijk_rank[] = {i, j, k};
769b072555Sjeremylt         if (i == i_rank[0] && j == i_rank[1] && k == i_rank[2]) return start;
779b072555Sjeremylt         GlobalNodes(p, ijk_rank, degree, mesh_elem, m_nodes);
789b072555Sjeremylt         start += m_nodes[0] * m_nodes[1] * m_nodes[2];
79cb32e2e7SValeria Barra       }
80cb32e2e7SValeria Barra     }
81cb32e2e7SValeria Barra   }
82cb32e2e7SValeria Barra   return -1;
83cb32e2e7SValeria Barra }
842b730f8bSJeremy L Thompson static PetscErrorCode CreateRestriction(Ceed ceed, const CeedInt mesh_elem[3], CeedInt P, CeedInt num_comp, CeedElemRestriction *elem_restr) {
859b072555Sjeremylt   const PetscInt num_elem = mesh_elem[0] * mesh_elem[1] * mesh_elem[2];
869b072555Sjeremylt   PetscInt       m_nodes[3], *idx, *idx_p;
87cb32e2e7SValeria Barra 
885dfaedb8SJed Brown   PetscFunctionBeginUser;
89cb32e2e7SValeria Barra   // Get indicies
909b072555Sjeremylt   for (int d = 0; d < 3; d++) m_nodes[d] = mesh_elem[d] * (P - 1) + 1;
919b072555Sjeremylt   idx_p = idx = malloc(num_elem * P * P * P * sizeof idx[0]);
922b730f8bSJeremy L Thompson   for (CeedInt i = 0; i < mesh_elem[0]; i++) {
932b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < mesh_elem[1]; j++) {
942b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < mesh_elem[2]; k++, idx_p += P * P * P) {
952b730f8bSJeremy L Thompson         for (CeedInt ii = 0; ii < P; ii++) {
962b730f8bSJeremy L Thompson           for (CeedInt jj = 0; jj < P; jj++) {
97cb32e2e7SValeria Barra             for (CeedInt kk = 0; kk < P; kk++) {
98cb32e2e7SValeria Barra               if (0) {  // This is the C-style (i,j,k) ordering that I prefer
992b730f8bSJeremy L Thompson                 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));
100cb32e2e7SValeria Barra               } else {  // (k,j,i) ordering for consistency with MFEM example
1012b730f8bSJeremy L Thompson                 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));
1022b730f8bSJeremy L Thompson               }
1032b730f8bSJeremy L Thompson             }
1042b730f8bSJeremy L Thompson           }
1052b730f8bSJeremy L Thompson         }
1062b730f8bSJeremy L Thompson       }
107cb32e2e7SValeria Barra     }
108cb32e2e7SValeria Barra   }
109cb32e2e7SValeria Barra 
110cb32e2e7SValeria Barra   // Setup CEED restriction
1112b730f8bSJeremy L Thompson   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,
1122b730f8bSJeremy L Thompson                             idx, elem_restr);
113cb32e2e7SValeria Barra 
114cb32e2e7SValeria Barra   PetscFunctionReturn(0);
115cb32e2e7SValeria Barra }
116cb32e2e7SValeria Barra 
117cb32e2e7SValeria Barra // Data for PETSc
118d4d45553Srezgarshakeri typedef struct OperatorApplyContext_ *OperatorApplyContext;
119d4d45553Srezgarshakeri struct OperatorApplyContext_ {
120cb32e2e7SValeria Barra   MPI_Comm     comm;
1219b072555Sjeremylt   VecScatter   l_to_g;    // Scatter for all entries
1229b072555Sjeremylt   VecScatter   l_to_g_0;  // Skip Dirichlet values
1239b072555Sjeremylt   VecScatter   g_to_g_D;  // global-to-global; only Dirichlet values
1249b072555Sjeremylt   Vec          X_loc, Y_loc;
1259b072555Sjeremylt   CeedVector   x_ceed, y_ceed;
126cb32e2e7SValeria Barra   CeedOperator op;
1279b072555Sjeremylt   CeedVector   q_data;
128cb32e2e7SValeria Barra   Ceed         ceed;
129cb32e2e7SValeria Barra };
130cb32e2e7SValeria Barra 
131cb32e2e7SValeria Barra // BP Options
1322b730f8bSJeremy L Thompson typedef enum { CEED_BP1 = 0, CEED_BP2 = 1, CEED_BP3 = 2, CEED_BP4 = 3, CEED_BP5 = 4, CEED_BP6 = 5 } BPType;
1332b730f8bSJeremy L Thompson static const char *const bp_types[] = {"bp1", "bp2", "bp3", "bp4", "bp5", "bp6", "BPType", "CEED_BP", 0};
134cb32e2e7SValeria Barra 
135cb32e2e7SValeria Barra // BP specific data
136cb32e2e7SValeria Barra typedef struct {
1379b072555Sjeremylt   CeedInt           num_comp_u, q_data_size, q_extra;
1389b072555Sjeremylt   CeedQFunctionUser setup_geo, setup_rhs, apply, error;
1399b072555Sjeremylt   const char       *setup_geo_loc, *setup_rhs_loc, *apply_loc, *error_loc;
1409b072555Sjeremylt   CeedEvalMode      in_mode, out_mode;
1419b072555Sjeremylt   CeedQuadMode      q_mode;
1429b072555Sjeremylt } BPData;
143cb32e2e7SValeria Barra 
1449b072555Sjeremylt BPData bp_options[6] = {
1452b730f8bSJeremy L Thompson     [CEED_BP1] = {.num_comp_u    = 1,
1469b072555Sjeremylt                   .q_data_size   = 1,
1479b072555Sjeremylt                   .q_extra       = 1,
1489b072555Sjeremylt                   .setup_geo     = SetupMassGeo,
1499b072555Sjeremylt                   .setup_rhs     = SetupMassRhs,
150cb32e2e7SValeria Barra                   .apply         = Mass,
151cb32e2e7SValeria Barra                   .error         = Error,
1529b072555Sjeremylt                   .setup_geo_loc = SetupMassGeo_loc,
1539b072555Sjeremylt                   .setup_rhs_loc = SetupMassRhs_loc,
1549b072555Sjeremylt                   .apply_loc     = Mass_loc,
1559b072555Sjeremylt                   .error_loc     = Error_loc,
1569b072555Sjeremylt                   .in_mode       = CEED_EVAL_INTERP,
1579b072555Sjeremylt                   .out_mode      = CEED_EVAL_INTERP,
1582b730f8bSJeremy L Thompson                   .q_mode        = CEED_GAUSS        },
1592b730f8bSJeremy L Thompson     [CEED_BP2] = {.num_comp_u    = 3,
1609b072555Sjeremylt                   .q_data_size   = 1,
1619b072555Sjeremylt                   .q_extra       = 1,
1629b072555Sjeremylt                   .setup_geo     = SetupMassGeo,
1639b072555Sjeremylt                   .setup_rhs     = SetupMassRhs3,
164cb32e2e7SValeria Barra                   .apply         = Mass3,
165cb32e2e7SValeria Barra                   .error         = Error3,
1669b072555Sjeremylt                   .setup_geo_loc = SetupMassGeo_loc,
1679b072555Sjeremylt                   .setup_rhs_loc = SetupMassRhs3_loc,
1689b072555Sjeremylt                   .apply_loc     = Mass3_loc,
1699b072555Sjeremylt                   .error_loc     = Error3_loc,
1709b072555Sjeremylt                   .in_mode       = CEED_EVAL_INTERP,
1719b072555Sjeremylt                   .out_mode      = CEED_EVAL_INTERP,
1722b730f8bSJeremy L Thompson                   .q_mode        = CEED_GAUSS        },
1732b730f8bSJeremy L Thompson     [CEED_BP3] = {.num_comp_u    = 1,
1749b072555Sjeremylt                   .q_data_size   = 7,
1759b072555Sjeremylt                   .q_extra       = 1,
1769b072555Sjeremylt                   .setup_geo     = SetupDiffGeo,
1779b072555Sjeremylt                   .setup_rhs     = SetupDiffRhs,
178cb32e2e7SValeria Barra                   .apply         = Diff,
179cb32e2e7SValeria Barra                   .error         = Error,
1809b072555Sjeremylt                   .setup_geo_loc = SetupDiffGeo_loc,
1819b072555Sjeremylt                   .setup_rhs_loc = SetupDiffRhs_loc,
1829b072555Sjeremylt                   .apply_loc     = Diff_loc,
1839b072555Sjeremylt                   .error_loc     = Error_loc,
1849b072555Sjeremylt                   .in_mode       = CEED_EVAL_GRAD,
1859b072555Sjeremylt                   .out_mode      = CEED_EVAL_GRAD,
1862b730f8bSJeremy L Thompson                   .q_mode        = CEED_GAUSS        },
1872b730f8bSJeremy L Thompson     [CEED_BP4] = {.num_comp_u    = 3,
1889b072555Sjeremylt                   .q_data_size   = 7,
1899b072555Sjeremylt                   .q_extra       = 1,
1909b072555Sjeremylt                   .setup_geo     = SetupDiffGeo,
1919b072555Sjeremylt                   .setup_rhs     = SetupDiffRhs3,
192cb32e2e7SValeria Barra                   .apply         = Diff3,
193cb32e2e7SValeria Barra                   .error         = Error3,
1949b072555Sjeremylt                   .setup_geo_loc = SetupDiffGeo_loc,
1959b072555Sjeremylt                   .setup_rhs_loc = SetupDiffRhs3_loc,
1969b072555Sjeremylt                   .apply_loc     = Diff3_loc,
1979b072555Sjeremylt                   .error_loc     = Error3_loc,
1989b072555Sjeremylt                   .in_mode       = CEED_EVAL_GRAD,
1999b072555Sjeremylt                   .out_mode      = CEED_EVAL_GRAD,
2002b730f8bSJeremy L Thompson                   .q_mode        = CEED_GAUSS        },
2012b730f8bSJeremy L Thompson     [CEED_BP5] = {.num_comp_u    = 1,
2029b072555Sjeremylt                   .q_data_size   = 7,
2039b072555Sjeremylt                   .q_extra       = 0,
2049b072555Sjeremylt                   .setup_geo     = SetupDiffGeo,
2059b072555Sjeremylt                   .setup_rhs     = SetupDiffRhs,
206cb32e2e7SValeria Barra                   .apply         = Diff,
207cb32e2e7SValeria Barra                   .error         = Error,
2089b072555Sjeremylt                   .setup_geo_loc = SetupDiffGeo_loc,
2099b072555Sjeremylt                   .setup_rhs_loc = SetupDiffRhs_loc,
2109b072555Sjeremylt                   .apply_loc     = Diff_loc,
2119b072555Sjeremylt                   .error_loc     = Error_loc,
2129b072555Sjeremylt                   .in_mode       = CEED_EVAL_GRAD,
2139b072555Sjeremylt                   .out_mode      = CEED_EVAL_GRAD,
2142b730f8bSJeremy L Thompson                   .q_mode        = CEED_GAUSS_LOBATTO},
2152b730f8bSJeremy L Thompson     [CEED_BP6] = {.num_comp_u    = 3,
2169b072555Sjeremylt                   .q_data_size   = 7,
2179b072555Sjeremylt                   .q_extra       = 0,
2189b072555Sjeremylt                   .setup_geo     = SetupDiffGeo,
2199b072555Sjeremylt                   .setup_rhs     = SetupDiffRhs3,
220cb32e2e7SValeria Barra                   .apply         = Diff3,
221cb32e2e7SValeria Barra                   .error         = Error3,
2229b072555Sjeremylt                   .setup_geo_loc = SetupDiffGeo_loc,
2239b072555Sjeremylt                   .setup_rhs_loc = SetupDiffRhs3_loc,
2249b072555Sjeremylt                   .apply_loc     = Diff3_loc,
2259b072555Sjeremylt                   .error_loc     = Error3_loc,
2269b072555Sjeremylt                   .in_mode       = CEED_EVAL_GRAD,
2279b072555Sjeremylt                   .out_mode      = CEED_EVAL_GRAD,
2282b730f8bSJeremy L Thompson                   .q_mode        = CEED_GAUSS_LOBATTO}
229cb32e2e7SValeria Barra };
230cb32e2e7SValeria Barra 
231cb32e2e7SValeria Barra // This function uses libCEED to compute the action of the mass matrix
232cb32e2e7SValeria Barra static PetscErrorCode MatMult_Mass(Mat A, Vec X, Vec Y) {
233d4d45553Srezgarshakeri   OperatorApplyContext op_apply_ctx;
234cb32e2e7SValeria Barra   PetscScalar         *x, *y;
2359b072555Sjeremylt   PetscMemType         x_mem_type, y_mem_type;
236cb32e2e7SValeria Barra 
237cb32e2e7SValeria Barra   PetscFunctionBeginUser;
2389396343dSjeremylt 
2392b730f8bSJeremy L Thompson   PetscCall(MatShellGetContext(A, &op_apply_ctx));
2409396343dSjeremylt 
2419396343dSjeremylt   // Global-to-local
2422b730f8bSJeremy L Thompson   PetscCall(VecScatterBegin(op_apply_ctx->l_to_g, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REVERSE));
2432b730f8bSJeremy L Thompson   PetscCall(VecScatterEnd(op_apply_ctx->l_to_g, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REVERSE));
244cb32e2e7SValeria Barra 
2459396343dSjeremylt   // Setup libCEED vectors
2462b730f8bSJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x, &x_mem_type));
2472b730f8bSJeremy L Thompson   PetscCall(VecGetArrayAndMemType(op_apply_ctx->Y_loc, &y, &y_mem_type));
2482b730f8bSJeremy L Thompson   CeedVectorSetArray(op_apply_ctx->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x);
2492b730f8bSJeremy L Thompson   CeedVectorSetArray(op_apply_ctx->y_ceed, MemTypeP2C(y_mem_type), CEED_USE_POINTER, y);
250cb32e2e7SValeria Barra 
2519396343dSjeremylt   // Apply libCEED operator
2522b730f8bSJeremy L Thompson   CeedOperatorApply(op_apply_ctx->op, op_apply_ctx->x_ceed, op_apply_ctx->y_ceed, CEED_REQUEST_IMMEDIATE);
253cb32e2e7SValeria Barra 
2549396343dSjeremylt   // Restore PETSc vectors
255d4d45553Srezgarshakeri   CeedVectorTakeArray(op_apply_ctx->x_ceed, MemTypeP2C(x_mem_type), NULL);
256d4d45553Srezgarshakeri   CeedVectorTakeArray(op_apply_ctx->y_ceed, MemTypeP2C(y_mem_type), NULL);
2572b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x));
2582b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(op_apply_ctx->Y_loc, &y));
259cb32e2e7SValeria Barra 
2609396343dSjeremylt   // Local-to-global
261cb32e2e7SValeria Barra   if (Y) {
2622b730f8bSJeremy L Thompson     PetscCall(VecZeroEntries(Y));
2632b730f8bSJeremy L Thompson     PetscCall(VecScatterBegin(op_apply_ctx->l_to_g, op_apply_ctx->Y_loc, Y, ADD_VALUES, SCATTER_FORWARD));
2642b730f8bSJeremy L Thompson     PetscCall(VecScatterEnd(op_apply_ctx->l_to_g, op_apply_ctx->Y_loc, Y, ADD_VALUES, SCATTER_FORWARD));
265cb32e2e7SValeria Barra   }
266cb32e2e7SValeria Barra   PetscFunctionReturn(0);
267cb32e2e7SValeria Barra }
268cb32e2e7SValeria Barra 
269*ea61e9acSJeremy L Thompson // This function uses libCEED to compute the action of the Laplacian with Dirichlet boundary conditions
270cb32e2e7SValeria Barra static PetscErrorCode MatMult_Diff(Mat A, Vec X, Vec Y) {
271d4d45553Srezgarshakeri   OperatorApplyContext op_apply_ctx;
272cb32e2e7SValeria Barra   PetscScalar         *x, *y;
2739b072555Sjeremylt   PetscMemType         x_mem_type, y_mem_type;
274cb32e2e7SValeria Barra 
275cb32e2e7SValeria Barra   PetscFunctionBeginUser;
2769396343dSjeremylt 
2772b730f8bSJeremy L Thompson   PetscCall(MatShellGetContext(A, &op_apply_ctx));
278cb32e2e7SValeria Barra 
279cb32e2e7SValeria Barra   // Global-to-local
2802b730f8bSJeremy L Thompson   PetscCall(VecScatterBegin(op_apply_ctx->l_to_g_0, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REVERSE));
2812b730f8bSJeremy L Thompson   PetscCall(VecScatterEnd(op_apply_ctx->l_to_g_0, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REVERSE));
282cb32e2e7SValeria Barra 
2839396343dSjeremylt   // Setup libCEED vectors
2842b730f8bSJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x, &x_mem_type));
2852b730f8bSJeremy L Thompson   PetscCall(VecGetArrayAndMemType(op_apply_ctx->Y_loc, &y, &y_mem_type));
2862b730f8bSJeremy L Thompson   CeedVectorSetArray(op_apply_ctx->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x);
2872b730f8bSJeremy L Thompson   CeedVectorSetArray(op_apply_ctx->y_ceed, MemTypeP2C(y_mem_type), CEED_USE_POINTER, y);
288cb32e2e7SValeria Barra 
2899396343dSjeremylt   // Apply libCEED operator
2902b730f8bSJeremy L Thompson   CeedOperatorApply(op_apply_ctx->op, op_apply_ctx->x_ceed, op_apply_ctx->y_ceed, CEED_REQUEST_IMMEDIATE);
291cb32e2e7SValeria Barra 
292cb32e2e7SValeria Barra   // Restore PETSc vectors
293d4d45553Srezgarshakeri   CeedVectorTakeArray(op_apply_ctx->x_ceed, MemTypeP2C(x_mem_type), NULL);
294d4d45553Srezgarshakeri   CeedVectorTakeArray(op_apply_ctx->y_ceed, MemTypeP2C(y_mem_type), NULL);
2952b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x));
2962b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(op_apply_ctx->Y_loc, &y));
297cb32e2e7SValeria Barra 
298cb32e2e7SValeria Barra   // Local-to-global
2992b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(Y));
3002b730f8bSJeremy L Thompson   PetscCall(VecScatterBegin(op_apply_ctx->g_to_g_D, X, Y, INSERT_VALUES, SCATTER_FORWARD));
3012b730f8bSJeremy L Thompson   PetscCall(VecScatterEnd(op_apply_ctx->g_to_g_D, X, Y, INSERT_VALUES, SCATTER_FORWARD));
3022b730f8bSJeremy L Thompson   PetscCall(VecScatterBegin(op_apply_ctx->l_to_g_0, op_apply_ctx->Y_loc, Y, ADD_VALUES, SCATTER_FORWARD));
3032b730f8bSJeremy L Thompson   PetscCall(VecScatterEnd(op_apply_ctx->l_to_g_0, op_apply_ctx->Y_loc, Y, ADD_VALUES, SCATTER_FORWARD));
304cb32e2e7SValeria Barra 
305cb32e2e7SValeria Barra   PetscFunctionReturn(0);
306cb32e2e7SValeria Barra }
307cb32e2e7SValeria Barra 
308cb32e2e7SValeria Barra // This function calculates the error in the final solution
3092b730f8bSJeremy L Thompson static PetscErrorCode ComputeErrorMax(OperatorApplyContext op_apply_ctx, CeedOperator op_error, Vec X, CeedVector target, PetscReal *max_error) {
310cb32e2e7SValeria Barra   PetscScalar *x;
3119b072555Sjeremylt   PetscMemType mem_type;
312cb32e2e7SValeria Barra   CeedVector   collocated_error;
3131f9221feSJeremy L Thompson   CeedSize     length;
314cb32e2e7SValeria Barra 
315cb32e2e7SValeria Barra   PetscFunctionBeginUser;
3169396343dSjeremylt 
317cb32e2e7SValeria Barra   CeedVectorGetLength(target, &length);
318d4d45553Srezgarshakeri   CeedVectorCreate(op_apply_ctx->ceed, length, &collocated_error);
319cb32e2e7SValeria Barra 
320cb32e2e7SValeria Barra   // Global-to-local
3212b730f8bSJeremy L Thompson   PetscCall(VecScatterBegin(op_apply_ctx->l_to_g, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REVERSE));
3222b730f8bSJeremy L Thompson   PetscCall(VecScatterEnd(op_apply_ctx->l_to_g, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REVERSE));
3239396343dSjeremylt 
3249396343dSjeremylt   // Setup libCEED vector
3252b730f8bSJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x, &mem_type));
3262b730f8bSJeremy L Thompson   CeedVectorSetArray(op_apply_ctx->x_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, x);
327cb32e2e7SValeria Barra 
3289396343dSjeremylt   // Apply libCEED operator
3292b730f8bSJeremy L Thompson   CeedOperatorApply(op_error, op_apply_ctx->x_ceed, collocated_error, CEED_REQUEST_IMMEDIATE);
330cb32e2e7SValeria Barra 
331cb32e2e7SValeria Barra   // Restore PETSc vector
332d4d45553Srezgarshakeri   CeedVectorTakeArray(op_apply_ctx->x_ceed, MemTypeP2C(mem_type), NULL);
3332b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x));
334cb32e2e7SValeria Barra 
335cb32e2e7SValeria Barra   // Reduce max error
3369b072555Sjeremylt   *max_error = 0;
337cb32e2e7SValeria Barra   const CeedScalar *e;
338cb32e2e7SValeria Barra   CeedVectorGetArrayRead(collocated_error, CEED_MEM_HOST, &e);
339cb32e2e7SValeria Barra   for (CeedInt i = 0; i < length; i++) {
3409b072555Sjeremylt     *max_error = PetscMax(*max_error, PetscAbsScalar(e[i]));
341cb32e2e7SValeria Barra   }
342cb32e2e7SValeria Barra   CeedVectorRestoreArrayRead(collocated_error, &e);
3432b730f8bSJeremy L Thompson   PetscCall(MPI_Allreduce(MPI_IN_PLACE, max_error, 1, MPIU_REAL, MPIU_MAX, op_apply_ctx->comm));
344cb32e2e7SValeria Barra 
345cb32e2e7SValeria Barra   // Cleanup
346cb32e2e7SValeria Barra   CeedVectorDestroy(&collocated_error);
347cb32e2e7SValeria Barra 
348cb32e2e7SValeria Barra   PetscFunctionReturn(0);
349cb32e2e7SValeria Barra }
350cb32e2e7SValeria Barra 
351cb32e2e7SValeria Barra int main(int argc, char **argv) {
352cb32e2e7SValeria Barra   MPI_Comm comm;
3539b072555Sjeremylt   char     ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self";
354cb32e2e7SValeria Barra   double   my_rt_start, my_rt, rt_min, rt_max;
3552b730f8bSJeremy L Thompson   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,
3562b730f8bSJeremy L Thompson                                                                                                                     ksp_max_it_clip[2];
357cb32e2e7SValeria Barra   PetscScalar         *r;
358cb32e2e7SValeria Barra   PetscBool            test_mode, benchmark_mode, write_solution;
359cb32e2e7SValeria Barra   PetscMPIInt          size, rank;
3609b072555Sjeremylt   PetscLogStage        solve_stage;
3619b072555Sjeremylt   Vec                  X, X_loc, rhs, rhs_loc;
362cb32e2e7SValeria Barra   Mat                  mat;
363cb32e2e7SValeria Barra   KSP                  ksp;
3649b072555Sjeremylt   VecScatter           l_to_g, l_to_g_0, g_to_g_D;
3659b072555Sjeremylt   PetscMemType         mem_type;
366d4d45553Srezgarshakeri   OperatorApplyContext op_apply_ctx;
367cb32e2e7SValeria Barra   Ceed                 ceed;
3689b072555Sjeremylt   CeedBasis            basis_x, basis_u;
3699b072555Sjeremylt   CeedElemRestriction  elem_restr_x, elem_restr_u, elem_restr_u_i, elem_restr_qd_i;
3709b072555Sjeremylt   CeedQFunction        qf_setup_geo, qf_setup_rhs, qf_apply, qf_error;
3719b072555Sjeremylt   CeedOperator         op_setup_geo, op_setup_rhs, op_apply, op_error;
3729b072555Sjeremylt   CeedVector           x_coord, q_data, rhs_ceed, target;
373cb32e2e7SValeria Barra   CeedInt              P, Q;
3749b072555Sjeremylt   const CeedInt        dim = 3, num_comp_x = 3;
3759b072555Sjeremylt   BPType               bp_choice;
376cb32e2e7SValeria Barra 
3772b730f8bSJeremy L Thompson   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
378cb32e2e7SValeria Barra   comm = PETSC_COMM_WORLD;
37932d2ee49SValeria Barra 
38032d2ee49SValeria Barra   // Read command line options
38167490bc6SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
3829b072555Sjeremylt   bp_choice = CEED_BP1;
3832b730f8bSJeremy L Thompson   PetscCall(PetscOptionsEnum("-problem", "CEED benchmark problem to solve", NULL, bp_types, (PetscEnum)bp_choice, (PetscEnum *)&bp_choice, NULL));
3849b072555Sjeremylt   num_comp_u = bp_options[bp_choice].num_comp_u;
385cb32e2e7SValeria Barra   test_mode  = PETSC_FALSE;
3862b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, test_mode, &test_mode, NULL));
387cb32e2e7SValeria Barra   benchmark_mode = PETSC_FALSE;
3882b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-benchmark", "Benchmarking mode (prints benchmark statistics)", NULL, benchmark_mode, &benchmark_mode, NULL));
389cb32e2e7SValeria Barra   write_solution = PETSC_FALSE;
3902b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-write_solution", "Write solution for visualization", NULL, write_solution, &write_solution, NULL));
391cb32e2e7SValeria Barra   degree = test_mode ? 3 : 1;
3922b730f8bSJeremy L Thompson   PetscCall(PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", NULL, degree, &degree, NULL));
3939b072555Sjeremylt   q_extra = bp_options[bp_choice].q_extra;
3942b730f8bSJeremy L Thompson   PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, q_extra, &q_extra, NULL));
3952b730f8bSJeremy L Thompson   PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, ceed_resource, ceed_resource, sizeof(ceed_resource), NULL));
3969b072555Sjeremylt   local_nodes = 1000;
3972b730f8bSJeremy L Thompson   PetscCall(PetscOptionsInt("-local", "Target number of locally owned nodes per process", NULL, local_nodes, &local_nodes, NULL));
3982fbc6e41SJeremy L Thompson   PetscInt two       = 2;
3992fbc6e41SJeremy L Thompson   ksp_max_it_clip[0] = 5;
4002fbc6e41SJeremy L Thompson   ksp_max_it_clip[1] = 20;
4012b730f8bSJeremy L Thompson   PetscCall(
4022b730f8bSJeremy L Thompson       PetscOptionsIntArray("-ksp_max_it_clip", "Min and max number of iterations to use during benchmarking", NULL, ksp_max_it_clip, &two, NULL));
40367490bc6SJeremy L Thompson   PetscOptionsEnd();
404cb32e2e7SValeria Barra   P = degree + 1;
4059b072555Sjeremylt   Q = P + q_extra;
406cb32e2e7SValeria Barra 
4079396343dSjeremylt   // Set up libCEED
4089b072555Sjeremylt   CeedInit(ceed_resource, &ceed);
4099b072555Sjeremylt   CeedMemType mem_type_backend;
4109b072555Sjeremylt   CeedGetPreferredMemType(ceed, &mem_type_backend);
4119396343dSjeremylt 
4129b072555Sjeremylt   VecType default_vec_type = NULL, vec_type;
4139b072555Sjeremylt   switch (mem_type_backend) {
4142b730f8bSJeremy L Thompson     case CEED_MEM_HOST:
4152b730f8bSJeremy L Thompson       default_vec_type = VECSTANDARD;
4162b730f8bSJeremy L Thompson       break;
417b68a8d79SJed Brown     case CEED_MEM_DEVICE: {
418b68a8d79SJed Brown       const char *resolved;
419b68a8d79SJed Brown       CeedGetResource(ceed, &resolved);
4209b072555Sjeremylt       if (strstr(resolved, "/gpu/cuda")) default_vec_type = VECCUDA;
4212b730f8bSJeremy L Thompson       else if (strstr(resolved, "/gpu/hip/occa")) default_vec_type = VECSTANDARD;  // https://github.com/CEED/libCEED/issues/678
4229b072555Sjeremylt       else if (strstr(resolved, "/gpu/hip")) default_vec_type = VECHIP;
4239b072555Sjeremylt       else default_vec_type = VECSTANDARD;
424b68a8d79SJed Brown     }
425b68a8d79SJed Brown   }
4269396343dSjeremylt 
427cb32e2e7SValeria Barra   // Determine size of process grid
4282b730f8bSJeremy L Thompson   PetscCall(MPI_Comm_size(comm, &size));
429cb32e2e7SValeria Barra   Split3(size, p, false);
430cb32e2e7SValeria Barra 
4319b072555Sjeremylt   // Find a nicely composite number of elements no less than local_nodes
4322b730f8bSJeremy L Thompson   for (local_elem = PetscMax(1, local_nodes / (degree * degree * degree));; local_elem++) {
4339b072555Sjeremylt     Split3(local_elem, mesh_elem, true);
4349b072555Sjeremylt     if (Max3(mesh_elem) / Min3(mesh_elem) <= 2) break;
435cb32e2e7SValeria Barra   }
436cb32e2e7SValeria Barra 
437cb32e2e7SValeria Barra   // Find my location in the process grid
4382b730f8bSJeremy L Thompson   PetscCall(MPI_Comm_rank(comm, &rank));
4399b072555Sjeremylt   for (int d = 0, rank_left = rank; d < dim; d++) {
440cb32e2e7SValeria Barra     const int pstride[3] = {p[1] * p[2], p[2], 1};
4419b072555Sjeremylt     i_rank[d]            = rank_left / pstride[d];
4429b072555Sjeremylt     rank_left -= i_rank[d] * pstride[d];
443cb32e2e7SValeria Barra   }
444cb32e2e7SValeria Barra 
4459b072555Sjeremylt   GlobalNodes(p, i_rank, degree, mesh_elem, m_nodes);
446cb32e2e7SValeria Barra 
447cb32e2e7SValeria Barra   // Setup global vector
4482b730f8bSJeremy L Thompson   PetscCall(VecCreate(comm, &X));
4492b730f8bSJeremy L Thompson   PetscCall(VecSetType(X, default_vec_type));
4502b730f8bSJeremy L Thompson   PetscCall(VecSetSizes(X, m_nodes[0] * m_nodes[1] * m_nodes[2] * num_comp_u, PETSC_DECIDE));
4512b730f8bSJeremy L Thompson   PetscCall(VecSetFromOptions(X));
4522b730f8bSJeremy L Thompson   PetscCall(VecSetUp(X));
453cb32e2e7SValeria Barra 
454cb32e2e7SValeria Barra   // Set up libCEED
4559b072555Sjeremylt   CeedInit(ceed_resource, &ceed);
456cb32e2e7SValeria Barra 
457cb32e2e7SValeria Barra   // Print summary
458cb32e2e7SValeria Barra   CeedInt gsize;
4592b730f8bSJeremy L Thompson   PetscCall(VecGetSize(X, &gsize));
460dc7d240cSValeria Barra   if (!test_mode) {
4619b072555Sjeremylt     const char *used_resource;
4629b072555Sjeremylt     CeedGetResource(ceed, &used_resource);
4639396343dSjeremylt 
4642b730f8bSJeremy L Thompson     PetscCall(VecGetType(X, &vec_type));
4659396343dSjeremylt 
4662b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm,
467990fdeb6SJeremy L Thompson                           "\n-- CEED Benchmark Problem %" CeedInt_FMT " -- libCEED + PETSc --\n"
4689396343dSjeremylt                           "  PETSc:\n"
4699396343dSjeremylt                           "    PETSc Vec Type                     : %s\n"
470cb32e2e7SValeria Barra                           "  libCEED:\n"
471cb32e2e7SValeria Barra                           "    libCEED Backend                    : %s\n"
4729396343dSjeremylt                           "    libCEED Backend MemType            : %s\n"
473cb32e2e7SValeria Barra                           "  Mesh:\n"
474751eb813Srezgarshakeri                           "    Solution Order (P)                 : %" CeedInt_FMT "\n"
475751eb813Srezgarshakeri                           "    Quadrature  Order (Q)              : %" CeedInt_FMT "\n"
47608140895SJed Brown                           "    Global nodes                       : %" PetscInt_FMT "\n"
4772b730f8bSJeremy L Thompson                           "    Process Decomposition              : %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n"
4782b730f8bSJeremy L Thompson                           "    Local Elements                     : %" PetscInt_FMT " = %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n"
4792b730f8bSJeremy L Thompson                           "    Owned nodes                        : %" PetscInt_FMT " = %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n"
48008140895SJed Brown                           "    DoF per node                       : %" PetscInt_FMT "\n",
4812b730f8bSJeremy L Thompson                           bp_choice + 1, vec_type, used_resource, CeedMemTypes[mem_type_backend], P, Q, gsize / num_comp_u, p[0], p[1], p[2],
4822b730f8bSJeremy L Thompson                           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],
4832b730f8bSJeremy L Thompson                           m_nodes[2], num_comp_u));
484cb32e2e7SValeria Barra   }
485cb32e2e7SValeria Barra 
486cb32e2e7SValeria Barra   {
4879b072555Sjeremylt     l_size = 1;
488cb32e2e7SValeria Barra     for (int d = 0; d < dim; d++) {
4899b072555Sjeremylt       l_nodes[d] = mesh_elem[d] * degree + 1;
4909b072555Sjeremylt       l_size *= l_nodes[d];
491cb32e2e7SValeria Barra     }
4922b730f8bSJeremy L Thompson     PetscCall(VecCreate(PETSC_COMM_SELF, &X_loc));
4932b730f8bSJeremy L Thompson     PetscCall(VecSetType(X_loc, default_vec_type));
4942b730f8bSJeremy L Thompson     PetscCall(VecSetSizes(X_loc, l_size * num_comp_u, PETSC_DECIDE));
4952b730f8bSJeremy L Thompson     PetscCall(VecSetFromOptions(X_loc));
4962b730f8bSJeremy L Thompson     PetscCall(VecSetUp(X_loc));
497cb32e2e7SValeria Barra 
498cb32e2e7SValeria Barra     // Create local-to-global scatter
4999b072555Sjeremylt     PetscInt *l_to_g_ind, *l_to_g_ind_0, *loc_ind, l_0_count;
5009b072555Sjeremylt     IS        l_to_g_is, l_to_g_is_0, loc_is;
5019b072555Sjeremylt     PetscInt  g_start[2][2][2], g_m_nodes[2][2][2][dim];
502cb32e2e7SValeria Barra 
503cb32e2e7SValeria Barra     for (int i = 0; i < 2; i++) {
504cb32e2e7SValeria Barra       for (int j = 0; j < 2; j++) {
505cb32e2e7SValeria Barra         for (int k = 0; k < 2; k++) {
5069b072555Sjeremylt           PetscInt ijk_rank[3] = {i_rank[0] + i, i_rank[1] + j, i_rank[2] + k};
5079b072555Sjeremylt           g_start[i][j][k]     = GlobalStart(p, ijk_rank, degree, mesh_elem);
5089b072555Sjeremylt           GlobalNodes(p, ijk_rank, degree, mesh_elem, g_m_nodes[i][j][k]);
509cb32e2e7SValeria Barra         }
510cb32e2e7SValeria Barra       }
511cb32e2e7SValeria Barra     }
512cb32e2e7SValeria Barra 
5132b730f8bSJeremy L Thompson     PetscCall(PetscMalloc1(l_size, &l_to_g_ind));
5142b730f8bSJeremy L Thompson     PetscCall(PetscMalloc1(l_size, &l_to_g_ind_0));
5152b730f8bSJeremy L Thompson     PetscCall(PetscMalloc1(l_size, &loc_ind));
5169b072555Sjeremylt     l_0_count = 0;
5172b730f8bSJeremy L Thompson     for (PetscInt i = 0, ir, ii; ir = i >= m_nodes[0], ii = i - ir * m_nodes[0], i < l_nodes[0]; i++) {
5182b730f8bSJeremy L Thompson       for (PetscInt j = 0, jr, jj; jr = j >= m_nodes[1], jj = j - jr * m_nodes[1], j < l_nodes[1]; j++) {
5192b730f8bSJeremy L Thompson         for (PetscInt k = 0, kr, kk; kr = k >= m_nodes[2], kk = k - kr * m_nodes[2], k < l_nodes[2]; k++) {
5209b072555Sjeremylt           PetscInt here    = (i * l_nodes[1] + j) * l_nodes[2] + k;
5212b730f8bSJeremy L Thompson           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;
5222b730f8bSJeremy L Thompson           if ((i_rank[0] == 0 && i == 0) || (i_rank[1] == 0 && j == 0) || (i_rank[2] == 0 && k == 0) ||
5232b730f8bSJeremy L Thompson               (i_rank[0] + 1 == p[0] && i + 1 == l_nodes[0]) || (i_rank[1] + 1 == p[1] && j + 1 == l_nodes[1]) ||
5242b730f8bSJeremy L Thompson               (i_rank[2] + 1 == p[2] && k + 1 == l_nodes[2]))
525cb32e2e7SValeria Barra             continue;
5269b072555Sjeremylt           l_to_g_ind_0[l_0_count] = l_to_g_ind[here];
5279b072555Sjeremylt           loc_ind[l_0_count++]    = here;
528cb32e2e7SValeria Barra         }
5292b730f8bSJeremy L Thompson       }
5302b730f8bSJeremy L Thompson     }
5312b730f8bSJeremy L Thompson     PetscCall(ISCreateBlock(comm, num_comp_u, l_size, l_to_g_ind, PETSC_OWN_POINTER, &l_to_g_is));
5322b730f8bSJeremy L Thompson     PetscCall(VecScatterCreate(X_loc, NULL, X, l_to_g_is, &l_to_g));
5332b730f8bSJeremy L Thompson     PetscCall(ISCreateBlock(comm, num_comp_u, l_0_count, l_to_g_ind_0, PETSC_OWN_POINTER, &l_to_g_is_0));
5342b730f8bSJeremy L Thompson     PetscCall(ISCreateBlock(comm, num_comp_u, l_0_count, loc_ind, PETSC_OWN_POINTER, &loc_is));
5352b730f8bSJeremy L Thompson     PetscCall(VecScatterCreate(X_loc, loc_is, X, l_to_g_is_0, &l_to_g_0));
536cb32e2e7SValeria Barra     {
537*ea61e9acSJeremy L Thompson       // 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)
5389b072555Sjeremylt       PetscInt           x_start, x_end, *ind_D, count_D = 0;
5399b072555Sjeremylt       IS                 is_D;
540cb32e2e7SValeria Barra       const PetscScalar *x;
5412b730f8bSJeremy L Thompson       PetscCall(VecZeroEntries(X_loc));
5422b730f8bSJeremy L Thompson       PetscCall(VecSet(X, 1.0));
5432b730f8bSJeremy L Thompson       PetscCall(VecScatterBegin(l_to_g_0, X_loc, X, INSERT_VALUES, SCATTER_FORWARD));
5442b730f8bSJeremy L Thompson       PetscCall(VecScatterEnd(l_to_g_0, X_loc, X, INSERT_VALUES, SCATTER_FORWARD));
5452b730f8bSJeremy L Thompson       PetscCall(VecGetOwnershipRange(X, &x_start, &x_end));
5462b730f8bSJeremy L Thompson       PetscCall(PetscMalloc1(x_end - x_start, &ind_D));
5472b730f8bSJeremy L Thompson       PetscCall(VecGetArrayRead(X, &x));
5489b072555Sjeremylt       for (PetscInt i = 0; i < x_end - x_start; i++) {
5499b072555Sjeremylt         if (x[i] == 1.) ind_D[count_D++] = x_start + i;
550cb32e2e7SValeria Barra       }
5512b730f8bSJeremy L Thompson       PetscCall(VecRestoreArrayRead(X, &x));
5522b730f8bSJeremy L Thompson       PetscCall(ISCreateGeneral(comm, count_D, ind_D, PETSC_COPY_VALUES, &is_D));
5532b730f8bSJeremy L Thompson       PetscCall(PetscFree(ind_D));
5542b730f8bSJeremy L Thompson       PetscCall(VecScatterCreate(X, is_D, X, is_D, &g_to_g_D));
5552b730f8bSJeremy L Thompson       PetscCall(ISDestroy(&is_D));
556cb32e2e7SValeria Barra     }
5572b730f8bSJeremy L Thompson     PetscCall(ISDestroy(&l_to_g_is));
5582b730f8bSJeremy L Thompson     PetscCall(ISDestroy(&l_to_g_is_0));
5592b730f8bSJeremy L Thompson     PetscCall(ISDestroy(&loc_is));
560cb32e2e7SValeria Barra   }
561cb32e2e7SValeria Barra 
562cb32e2e7SValeria Barra   // CEED bases
5632b730f8bSJeremy L Thompson   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_u, P, Q, bp_options[bp_choice].q_mode, &basis_u);
5642b730f8bSJeremy L Thompson   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, bp_options[bp_choice].q_mode, &basis_x);
565cb32e2e7SValeria Barra 
566cb32e2e7SValeria Barra   // CEED restrictions
5679b072555Sjeremylt   CreateRestriction(ceed, mesh_elem, P, num_comp_u, &elem_restr_u);
5689b072555Sjeremylt   CreateRestriction(ceed, mesh_elem, 2, dim, &elem_restr_x);
5699b072555Sjeremylt   CeedInt num_elem = mesh_elem[0] * mesh_elem[1] * mesh_elem[2];
5702b730f8bSJeremy L Thompson   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);
5712b730f8bSJeremy L Thompson   CeedElemRestrictionCreateStrided(ceed, num_elem, Q * Q * Q, bp_options[bp_choice].q_data_size,
5722b730f8bSJeremy L Thompson                                    bp_options[bp_choice].q_data_size * num_elem * Q * Q * Q, CEED_STRIDES_BACKEND, &elem_restr_qd_i);
573cb32e2e7SValeria Barra   {
5749b072555Sjeremylt     CeedScalar *x_loc;
5752b730f8bSJeremy L Thompson     CeedInt     shape[3] = {mesh_elem[0] + 1, mesh_elem[1] + 1, mesh_elem[2] + 1}, len = shape[0] * shape[1] * shape[2];
5769b072555Sjeremylt     x_loc = malloc(len * num_comp_x * sizeof x_loc[0]);
577cb32e2e7SValeria Barra     for (CeedInt i = 0; i < shape[0]; i++) {
578cb32e2e7SValeria Barra       for (CeedInt j = 0; j < shape[1]; j++) {
579cb32e2e7SValeria Barra         for (CeedInt k = 0; k < shape[2]; k++) {
5802b730f8bSJeremy L Thompson           x_loc[dim * ((i * shape[1] + j) * shape[2] + k) + 0] = 1. * (i_rank[0] * mesh_elem[0] + i) / (p[0] * mesh_elem[0]);
5812b730f8bSJeremy L Thompson           x_loc[dim * ((i * shape[1] + j) * shape[2] + k) + 1] = 1. * (i_rank[1] * mesh_elem[1] + j) / (p[1] * mesh_elem[1]);
5822b730f8bSJeremy L Thompson           x_loc[dim * ((i * shape[1] + j) * shape[2] + k) + 2] = 1. * (i_rank[2] * mesh_elem[2] + k) / (p[2] * mesh_elem[2]);
583cb32e2e7SValeria Barra         }
584cb32e2e7SValeria Barra       }
585cb32e2e7SValeria Barra     }
5869b072555Sjeremylt     CeedVectorCreate(ceed, len * num_comp_x, &x_coord);
5879b072555Sjeremylt     CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_OWN_POINTER, x_loc);
588cb32e2e7SValeria Barra   }
589cb32e2e7SValeria Barra 
5909b072555Sjeremylt   // Create the QFunction that builds the operator quadrature data
5912b730f8bSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].setup_geo, bp_options[bp_choice].setup_geo_loc, &qf_setup_geo);
5929b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_geo, "x", num_comp_x, CEED_EVAL_INTERP);
5939b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_geo, "dx", num_comp_x * dim, CEED_EVAL_GRAD);
5949b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_geo, "weight", 1, CEED_EVAL_WEIGHT);
5952b730f8bSJeremy L Thompson   CeedQFunctionAddOutput(qf_setup_geo, "q_data", bp_options[bp_choice].q_data_size, CEED_EVAL_NONE);
596cb32e2e7SValeria Barra 
5979b072555Sjeremylt   // Create the QFunction that sets up the RHS and true solution
5982b730f8bSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].setup_rhs, bp_options[bp_choice].setup_rhs_loc, &qf_setup_rhs);
5999b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_rhs, "x", num_comp_x, CEED_EVAL_INTERP);
6002b730f8bSJeremy L Thompson   CeedQFunctionAddInput(qf_setup_rhs, "q_data", bp_options[bp_choice].q_data_size, CEED_EVAL_NONE);
6019b072555Sjeremylt   CeedQFunctionAddOutput(qf_setup_rhs, "true_soln", num_comp_u, CEED_EVAL_NONE);
6029b072555Sjeremylt   CeedQFunctionAddOutput(qf_setup_rhs, "rhs", num_comp_u, CEED_EVAL_INTERP);
603cb32e2e7SValeria Barra 
604cb32e2e7SValeria Barra   // Set up PDE operator
6052b730f8bSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].apply, bp_options[bp_choice].apply_loc, &qf_apply);
606cb32e2e7SValeria Barra   // Add inputs and outputs
6079b072555Sjeremylt   CeedInt in_scale  = bp_options[bp_choice].in_mode == CEED_EVAL_GRAD ? 3 : 1;
6089b072555Sjeremylt   CeedInt out_scale = bp_options[bp_choice].out_mode == CEED_EVAL_GRAD ? 3 : 1;
6092b730f8bSJeremy L Thompson   CeedQFunctionAddInput(qf_apply, "u", num_comp_u * in_scale, bp_options[bp_choice].in_mode);
6102b730f8bSJeremy L Thompson   CeedQFunctionAddInput(qf_apply, "q_data", bp_options[bp_choice].q_data_size, CEED_EVAL_NONE);
6112b730f8bSJeremy L Thompson   CeedQFunctionAddOutput(qf_apply, "v", num_comp_u * out_scale, bp_options[bp_choice].out_mode);
612cb32e2e7SValeria Barra 
613cb32e2e7SValeria Barra   // Create the error qfunction
6142b730f8bSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].error, bp_options[bp_choice].error_loc, &qf_error);
6159b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP);
6169b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE);
6172b730f8bSJeremy L Thompson   CeedQFunctionAddInput(qf_error, "qdata", bp_options[bp_choice].q_data_size, CEED_EVAL_NONE);
6189b072555Sjeremylt   CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_NONE);
619cb32e2e7SValeria Barra 
620cb32e2e7SValeria Barra   // Create the persistent vectors that will be needed in setup
6219b072555Sjeremylt   CeedInt num_qpts;
6229b072555Sjeremylt   CeedBasisGetNumQuadraturePoints(basis_u, &num_qpts);
6232b730f8bSJeremy L Thompson   CeedVectorCreate(ceed, bp_options[bp_choice].q_data_size * num_elem * num_qpts, &q_data);
6249b072555Sjeremylt   CeedVectorCreate(ceed, num_elem * num_qpts * num_comp_u, &target);
6259b072555Sjeremylt   CeedVectorCreate(ceed, l_size * num_comp_u, &rhs_ceed);
626cb32e2e7SValeria Barra 
627cb32e2e7SValeria Barra   // Create the operator that builds the quadrature data for the ceed operator
6282b730f8bSJeremy L Thompson   CeedOperatorCreate(ceed, qf_setup_geo, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setup_geo);
6292b730f8bSJeremy L Thompson   CeedOperatorSetField(op_setup_geo, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
6302b730f8bSJeremy L Thompson   CeedOperatorSetField(op_setup_geo, "dx", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
6312b730f8bSJeremy L Thompson   CeedOperatorSetField(op_setup_geo, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE);
6322b730f8bSJeremy L Thompson   CeedOperatorSetField(op_setup_geo, "q_data", elem_restr_qd_i, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
633cb32e2e7SValeria Barra 
634cb32e2e7SValeria Barra   // Create the operator that builds the RHS and true solution
6352b730f8bSJeremy L Thompson   CeedOperatorCreate(ceed, qf_setup_rhs, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setup_rhs);
6362b730f8bSJeremy L Thompson   CeedOperatorSetField(op_setup_rhs, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
6372b730f8bSJeremy L Thompson   CeedOperatorSetField(op_setup_rhs, "q_data", elem_restr_qd_i, CEED_BASIS_COLLOCATED, q_data);
6382b730f8bSJeremy L Thompson   CeedOperatorSetField(op_setup_rhs, "true_soln", elem_restr_u_i, CEED_BASIS_COLLOCATED, target);
6392b730f8bSJeremy L Thompson   CeedOperatorSetField(op_setup_rhs, "rhs", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
640cb32e2e7SValeria Barra 
641cb32e2e7SValeria Barra   // Create the mass or diff operator
6422b730f8bSJeremy L Thompson   CeedOperatorCreate(ceed, qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_apply);
6439b072555Sjeremylt   CeedOperatorSetField(op_apply, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
6442b730f8bSJeremy L Thompson   CeedOperatorSetField(op_apply, "q_data", elem_restr_qd_i, CEED_BASIS_COLLOCATED, q_data);
6459b072555Sjeremylt   CeedOperatorSetField(op_apply, "v", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
646cb32e2e7SValeria Barra 
647cb32e2e7SValeria Barra   // Create the error operator
6482b730f8bSJeremy L Thompson   CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_error);
6499b072555Sjeremylt   CeedOperatorSetField(op_error, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
6502b730f8bSJeremy L Thompson   CeedOperatorSetField(op_error, "true_soln", elem_restr_u_i, CEED_BASIS_COLLOCATED, target);
6512b730f8bSJeremy L Thompson   CeedOperatorSetField(op_error, "qdata", elem_restr_qd_i, CEED_BASIS_COLLOCATED, q_data);
6522b730f8bSJeremy L Thompson   CeedOperatorSetField(op_error, "error", elem_restr_u_i, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
653cb32e2e7SValeria Barra 
654cb32e2e7SValeria Barra   // Set up Mat
6552b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(1, &op_apply_ctx));
656d4d45553Srezgarshakeri   op_apply_ctx->comm   = comm;
657d4d45553Srezgarshakeri   op_apply_ctx->l_to_g = l_to_g;
6589b072555Sjeremylt   if (bp_choice != CEED_BP1 && bp_choice != CEED_BP2) {
659d4d45553Srezgarshakeri     op_apply_ctx->l_to_g_0 = l_to_g_0;
660d4d45553Srezgarshakeri     op_apply_ctx->g_to_g_D = g_to_g_D;
661cb32e2e7SValeria Barra   }
662d4d45553Srezgarshakeri   op_apply_ctx->X_loc = X_loc;
6632b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(X_loc, &op_apply_ctx->Y_loc));
664d4d45553Srezgarshakeri   CeedVectorCreate(ceed, l_size * num_comp_u, &op_apply_ctx->x_ceed);
665d4d45553Srezgarshakeri   CeedVectorCreate(ceed, l_size * num_comp_u, &op_apply_ctx->y_ceed);
666d4d45553Srezgarshakeri   op_apply_ctx->op     = op_apply;
667d4d45553Srezgarshakeri   op_apply_ctx->q_data = q_data;
668d4d45553Srezgarshakeri   op_apply_ctx->ceed   = ceed;
669cb32e2e7SValeria Barra 
6702b730f8bSJeremy L Thompson   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,
6712b730f8bSJeremy L Thompson                            PETSC_DECIDE, op_apply_ctx, &mat));
6729b072555Sjeremylt   if (bp_choice == CEED_BP1 || bp_choice == CEED_BP2) {
6732b730f8bSJeremy L Thompson     PetscCall(MatShellSetOperation(mat, MATOP_MULT, (void (*)(void))MatMult_Mass));
674cb32e2e7SValeria Barra   } else {
6752b730f8bSJeremy L Thompson     PetscCall(MatShellSetOperation(mat, MATOP_MULT, (void (*)(void))MatMult_Diff));
676cb32e2e7SValeria Barra   }
6772b730f8bSJeremy L Thompson   PetscCall(VecGetType(X, &vec_type));
6782b730f8bSJeremy L Thompson   PetscCall(MatShellSetVecType(mat, vec_type));
679cb32e2e7SValeria Barra 
680cb32e2e7SValeria Barra   // Get RHS vector
6812b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(X, &rhs));
6822b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(X_loc, &rhs_loc));
6832b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(rhs_loc));
6842b730f8bSJeremy L Thompson   PetscCall(VecGetArrayAndMemType(rhs_loc, &r, &mem_type));
6859b072555Sjeremylt   CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r);
686cb32e2e7SValeria Barra 
6879b072555Sjeremylt   // Setup q_data, rhs, and target
6889b072555Sjeremylt   CeedOperatorApply(op_setup_geo, x_coord, q_data, CEED_REQUEST_IMMEDIATE);
6899b072555Sjeremylt   CeedOperatorApply(op_setup_rhs, x_coord, rhs_ceed, CEED_REQUEST_IMMEDIATE);
6909b072555Sjeremylt   CeedVectorDestroy(&x_coord);
691cb32e2e7SValeria Barra 
692cb32e2e7SValeria Barra   // Gather RHS
6932b730f8bSJeremy L Thompson   PetscCall(CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL));
6942b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(rhs_loc, &r));
6952b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(rhs));
6962b730f8bSJeremy L Thompson   PetscCall(VecScatterBegin(l_to_g, rhs_loc, rhs, ADD_VALUES, SCATTER_FORWARD));
6972b730f8bSJeremy L Thompson   PetscCall(VecScatterEnd(l_to_g, rhs_loc, rhs, ADD_VALUES, SCATTER_FORWARD));
6989b072555Sjeremylt   CeedVectorDestroy(&rhs_ceed);
699cb32e2e7SValeria Barra 
7002b730f8bSJeremy L Thompson   PetscCall(KSPCreate(comm, &ksp));
701cb32e2e7SValeria Barra   {
702cb32e2e7SValeria Barra     PC pc;
7032b730f8bSJeremy L Thompson     PetscCall(KSPGetPC(ksp, &pc));
7049b072555Sjeremylt     if (bp_choice == CEED_BP1 || bp_choice == CEED_BP2) {
7052b730f8bSJeremy L Thompson       PetscCall(PCSetType(pc, PCJACOBI));
7062b730f8bSJeremy L Thompson       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
707cb32e2e7SValeria Barra     } else {
7082b730f8bSJeremy L Thompson       PetscCall(PCSetType(pc, PCNONE));
709cb32e2e7SValeria Barra     }
7102b730f8bSJeremy L Thompson     PetscCall(KSPSetType(ksp, KSPCG));
7112b730f8bSJeremy L Thompson     PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
7122b730f8bSJeremy L Thompson     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
713cb32e2e7SValeria Barra   }
7142b730f8bSJeremy L Thompson   PetscCall(KSPSetOperators(ksp, mat, mat));
715da9108adSvaleriabarra   // First run's performance log is not considered for benchmarking purposes
7162b730f8bSJeremy L Thompson   PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1));
717cb32e2e7SValeria Barra   my_rt_start = MPI_Wtime();
7182b730f8bSJeremy L Thompson   PetscCall(KSPSolve(ksp, rhs, X));
719cb32e2e7SValeria Barra   my_rt = MPI_Wtime() - my_rt_start;
7202b730f8bSJeremy L Thompson   PetscCall(MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm));
721cb32e2e7SValeria Barra   // Set maxits based on first iteration timing
722cb32e2e7SValeria Barra   if (my_rt > 0.02) {
7232b730f8bSJeremy L Thompson     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, ksp_max_it_clip[0]));
724cb32e2e7SValeria Barra   } else {
7252b730f8bSJeremy L Thompson     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, ksp_max_it_clip[1]));
726cb32e2e7SValeria Barra   }
7272b730f8bSJeremy L Thompson   PetscCall(KSPSetFromOptions(ksp));
72809a940d7Sjeremylt 
729cb32e2e7SValeria Barra   // Timed solve
7302b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(X));
7312b730f8bSJeremy L Thompson   PetscCall(PetscBarrier((PetscObject)ksp));
73209a940d7Sjeremylt 
73309a940d7Sjeremylt   // -- Performance logging
7342b730f8bSJeremy L Thompson   PetscCall(PetscLogStageRegister("Solve Stage", &solve_stage));
7352b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePush(solve_stage));
73609a940d7Sjeremylt 
73709a940d7Sjeremylt   // -- Solve
738cb32e2e7SValeria Barra   my_rt_start = MPI_Wtime();
7392b730f8bSJeremy L Thompson   PetscCall(KSPSolve(ksp, rhs, X));
740cb32e2e7SValeria Barra   my_rt = MPI_Wtime() - my_rt_start;
74109a940d7Sjeremylt 
74209a940d7Sjeremylt   // -- Performance logging
7432b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePop());
74409a940d7Sjeremylt 
7458e87e98bSJed Brown   // Output results
746cb32e2e7SValeria Barra   {
7479b072555Sjeremylt     KSPType            ksp_type;
748cb32e2e7SValeria Barra     KSPConvergedReason reason;
749cb32e2e7SValeria Barra     PetscReal          rnorm;
750cb32e2e7SValeria Barra     PetscInt           its;
7512b730f8bSJeremy L Thompson     PetscCall(KSPGetType(ksp, &ksp_type));
7522b730f8bSJeremy L Thompson     PetscCall(KSPGetConvergedReason(ksp, &reason));
7532b730f8bSJeremy L Thompson     PetscCall(KSPGetIterationNumber(ksp, &its));
7542b730f8bSJeremy L Thompson     PetscCall(KSPGetResidualNorm(ksp, &rnorm));
755a22c4fb5SJed Brown     if (!test_mode || reason < 0 || rnorm > 1e-8) {
7562b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm,
757cb32e2e7SValeria Barra                             "  KSP:\n"
758cb32e2e7SValeria Barra                             "    KSP Type                           : %s\n"
759cb32e2e7SValeria Barra                             "    KSP Convergence                    : %s\n"
76008140895SJed Brown                             "    Total KSP Iterations               : %" PetscInt_FMT "\n"
761cb32e2e7SValeria Barra                             "    Final rnorm                        : %e\n",
7622b730f8bSJeremy L Thompson                             ksp_type, KSPConvergedReasons[reason], its, (double)rnorm));
763cb32e2e7SValeria Barra     }
7648e87e98bSJed Brown     if (!test_mode) {
7652b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm, "  Performance:\n"));
7668e87e98bSJed Brown     }
7678e87e98bSJed Brown     {
7689b072555Sjeremylt       PetscReal max_error;
7692b730f8bSJeremy L Thompson       PetscCall(ComputeErrorMax(op_apply_ctx, op_error, X, target, &max_error));
7708e87e98bSJed Brown       PetscReal tol = 5e-2;
7719b072555Sjeremylt       if (!test_mode || max_error > tol) {
7722b730f8bSJeremy L Thompson         PetscCall(MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm));
7732b730f8bSJeremy L Thompson         PetscCall(MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm));
7742b730f8bSJeremy L Thompson         PetscCall(PetscPrintf(comm,
7758e87e98bSJed Brown                               "    Pointwise Error (max)              : %e\n"
7768e87e98bSJed Brown                               "    CG Solve Time                      : %g (%g) sec\n",
7772b730f8bSJeremy L Thompson                               (double)max_error, rt_max, rt_min));
778cb32e2e7SValeria Barra       }
779cb32e2e7SValeria Barra     }
7808d0bb2bbSvaleriabarra     if (!test_mode) {
7812b730f8bSJeremy L Thompson       PetscCall(
7822b730f8bSJeremy L Thompson           PetscPrintf(comm, "    DoFs/Sec in CG                     : %g (%g) million\n", 1e-6 * gsize * its / rt_max, 1e-6 * gsize * its / rt_min));
783cb32e2e7SValeria Barra     }
784cb32e2e7SValeria Barra   }
785cb32e2e7SValeria Barra 
786cb32e2e7SValeria Barra   if (write_solution) {
7879b072555Sjeremylt     PetscViewer vtk_viewer_soln;
788cb32e2e7SValeria Barra 
7892b730f8bSJeremy L Thompson     PetscCall(PetscViewerCreate(comm, &vtk_viewer_soln));
7902b730f8bSJeremy L Thompson     PetscCall(PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK));
7912b730f8bSJeremy L Thompson     PetscCall(PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"));
7922b730f8bSJeremy L Thompson     PetscCall(VecView(X, vtk_viewer_soln));
7932b730f8bSJeremy L Thompson     PetscCall(PetscViewerDestroy(&vtk_viewer_soln));
794cb32e2e7SValeria Barra   }
795cb32e2e7SValeria Barra 
7962b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&rhs));
7972b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&rhs_loc));
7982b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&X));
7992b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&op_apply_ctx->X_loc));
8002b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&op_apply_ctx->Y_loc));
8012b730f8bSJeremy L Thompson   PetscCall(VecScatterDestroy(&l_to_g));
8022b730f8bSJeremy L Thompson   PetscCall(VecScatterDestroy(&l_to_g_0));
8032b730f8bSJeremy L Thompson   PetscCall(VecScatterDestroy(&g_to_g_D));
8042b730f8bSJeremy L Thompson   PetscCall(MatDestroy(&mat));
8052b730f8bSJeremy L Thompson   PetscCall(KSPDestroy(&ksp));
806cb32e2e7SValeria Barra 
807d4d45553Srezgarshakeri   CeedVectorDestroy(&op_apply_ctx->x_ceed);
808d4d45553Srezgarshakeri   CeedVectorDestroy(&op_apply_ctx->y_ceed);
809d4d45553Srezgarshakeri   CeedVectorDestroy(&op_apply_ctx->q_data);
810cb32e2e7SValeria Barra   CeedVectorDestroy(&target);
8119b072555Sjeremylt   CeedOperatorDestroy(&op_setup_geo);
8129b072555Sjeremylt   CeedOperatorDestroy(&op_setup_rhs);
8139b072555Sjeremylt   CeedOperatorDestroy(&op_apply);
8149b072555Sjeremylt   CeedOperatorDestroy(&op_error);
8159b072555Sjeremylt   CeedElemRestrictionDestroy(&elem_restr_u);
8169b072555Sjeremylt   CeedElemRestrictionDestroy(&elem_restr_x);
8179b072555Sjeremylt   CeedElemRestrictionDestroy(&elem_restr_u_i);
8189b072555Sjeremylt   CeedElemRestrictionDestroy(&elem_restr_qd_i);
8199b072555Sjeremylt   CeedQFunctionDestroy(&qf_setup_geo);
8209b072555Sjeremylt   CeedQFunctionDestroy(&qf_setup_rhs);
8219b072555Sjeremylt   CeedQFunctionDestroy(&qf_apply);
8229b072555Sjeremylt   CeedQFunctionDestroy(&qf_error);
8239b072555Sjeremylt   CeedBasisDestroy(&basis_u);
8249b072555Sjeremylt   CeedBasisDestroy(&basis_x);
825cb32e2e7SValeria Barra   CeedDestroy(&ceed);
8262b730f8bSJeremy L Thompson   PetscCall(PetscFree(op_apply_ctx));
827cb32e2e7SValeria Barra   return PetscFinalize();
828cb32e2e7SValeria Barra }
829