1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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 //
10ea61e9acSJeremy 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 //
12ea61e9acSJeremy 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>
3549aac155SJeremy L Thompson #include <petscdm.h>
363d576824SJeremy L Thompson #include <petscksp.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
MemTypeP2C(PetscMemType mem_type)462b730f8bSJeremy L Thompson static CeedMemType MemTypeP2C(PetscMemType mem_type) { return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST; }
47b68a8d79SJed Brown
Split3(PetscInt size,PetscInt m[3],bool reverse)48cb32e2e7SValeria Barra static void Split3(PetscInt size, PetscInt m[3], bool reverse) {
499b072555Sjeremylt for (PetscInt d = 0, size_left = size; d < 3; d++) {
509b072555Sjeremylt PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(size_left, 1. / (3 - d)));
519b072555Sjeremylt while (try * (size_left / try) != size_left) try++;
52cb32e2e7SValeria Barra m[reverse ? 2 - d : d] = try;
539b072555Sjeremylt size_left /= try;
54cb32e2e7SValeria Barra }
55cb32e2e7SValeria Barra }
56cb32e2e7SValeria Barra
Max3(const PetscInt a[3])572b730f8bSJeremy L Thompson static PetscInt Max3(const PetscInt a[3]) { return PetscMax(a[0], PetscMax(a[1], a[2])); }
Min3(const PetscInt a[3])582b730f8bSJeremy L Thompson static PetscInt Min3(const PetscInt a[3]) { return PetscMin(a[0], PetscMin(a[1], a[2])); }
GlobalNodes(const PetscInt p[3],const PetscInt i_rank[3],PetscInt degree,const PetscInt mesh_elem[3],PetscInt m_nodes[3])592b730f8bSJeremy 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]) {
602b730f8bSJeremy L Thompson for (int d = 0; d < 3; d++) m_nodes[d] = degree * mesh_elem[d] + (i_rank[d] == p[d] - 1);
61cb32e2e7SValeria Barra }
GlobalStart(const PetscInt p[3],const PetscInt i_rank[3],PetscInt degree,const PetscInt mesh_elem[3])622b730f8bSJeremy L Thompson static PetscInt GlobalStart(const PetscInt p[3], const PetscInt i_rank[3], PetscInt degree, const PetscInt mesh_elem[3]) {
63cb32e2e7SValeria Barra PetscInt start = 0;
64cb32e2e7SValeria Barra // Dumb brute-force is easier to read
65cb32e2e7SValeria Barra for (PetscInt i = 0; i < p[0]; i++) {
66cb32e2e7SValeria Barra for (PetscInt j = 0; j < p[1]; j++) {
67cb32e2e7SValeria Barra for (PetscInt k = 0; k < p[2]; k++) {
689b072555Sjeremylt PetscInt m_nodes[3], ijk_rank[] = {i, j, k};
699b072555Sjeremylt if (i == i_rank[0] && j == i_rank[1] && k == i_rank[2]) return start;
709b072555Sjeremylt GlobalNodes(p, ijk_rank, degree, mesh_elem, m_nodes);
719b072555Sjeremylt start += m_nodes[0] * m_nodes[1] * m_nodes[2];
72cb32e2e7SValeria Barra }
73cb32e2e7SValeria Barra }
74cb32e2e7SValeria Barra }
75cb32e2e7SValeria Barra return -1;
76cb32e2e7SValeria Barra }
CreateRestriction(Ceed ceed,const PetscInt mesh_elem[3],CeedInt P,CeedInt num_comp,CeedElemRestriction * elem_restr)774d00b080SJeremy L Thompson static PetscErrorCode CreateRestriction(Ceed ceed, const PetscInt mesh_elem[3], CeedInt P, CeedInt num_comp, CeedElemRestriction *elem_restr) {
789b072555Sjeremylt const PetscInt num_elem = mesh_elem[0] * mesh_elem[1] * mesh_elem[2];
794d00b080SJeremy L Thompson CeedInt m_nodes[3], *idx, *idx_p;
80cb32e2e7SValeria Barra
815dfaedb8SJed Brown PetscFunctionBeginUser;
82cb32e2e7SValeria Barra // Get indicies
839b072555Sjeremylt for (int d = 0; d < 3; d++) m_nodes[d] = mesh_elem[d] * (P - 1) + 1;
849b072555Sjeremylt idx_p = idx = malloc(num_elem * P * P * P * sizeof idx[0]);
852b730f8bSJeremy L Thompson for (CeedInt i = 0; i < mesh_elem[0]; i++) {
862b730f8bSJeremy L Thompson for (CeedInt j = 0; j < mesh_elem[1]; j++) {
872b730f8bSJeremy L Thompson for (CeedInt k = 0; k < mesh_elem[2]; k++, idx_p += P * P * P) {
882b730f8bSJeremy L Thompson for (CeedInt ii = 0; ii < P; ii++) {
892b730f8bSJeremy L Thompson for (CeedInt jj = 0; jj < P; jj++) {
90cb32e2e7SValeria Barra for (CeedInt kk = 0; kk < P; kk++) {
91cb32e2e7SValeria Barra if (0) { // This is the C-style (i,j,k) ordering that I prefer
922b730f8bSJeremy 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));
93cb32e2e7SValeria Barra } else { // (k,j,i) ordering for consistency with MFEM example
942b730f8bSJeremy 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));
952b730f8bSJeremy L Thompson }
962b730f8bSJeremy L Thompson }
972b730f8bSJeremy L Thompson }
982b730f8bSJeremy L Thompson }
992b730f8bSJeremy L Thompson }
100cb32e2e7SValeria Barra }
101cb32e2e7SValeria Barra }
102cb32e2e7SValeria Barra
103cb32e2e7SValeria Barra // Setup CEED restriction
1042b730f8bSJeremy 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,
1052b730f8bSJeremy L Thompson idx, elem_restr);
106ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
107cb32e2e7SValeria Barra }
108cb32e2e7SValeria Barra
109cb32e2e7SValeria Barra // Data for PETSc
110d4d45553Srezgarshakeri typedef struct OperatorApplyContext_ *OperatorApplyContext;
111d4d45553Srezgarshakeri struct OperatorApplyContext_ {
112cb32e2e7SValeria Barra MPI_Comm comm;
1139b072555Sjeremylt VecScatter l_to_g; // Scatter for all entries
1149b072555Sjeremylt VecScatter l_to_g_0; // Skip Dirichlet values
1159b072555Sjeremylt VecScatter g_to_g_D; // global-to-global; only Dirichlet values
1169b072555Sjeremylt Vec X_loc, Y_loc;
1179b072555Sjeremylt CeedVector x_ceed, y_ceed;
118cb32e2e7SValeria Barra CeedOperator op;
1199b072555Sjeremylt CeedVector q_data;
120cb32e2e7SValeria Barra Ceed ceed;
121cb32e2e7SValeria Barra };
122cb32e2e7SValeria Barra
123cb32e2e7SValeria Barra // BP Options
1242b730f8bSJeremy L Thompson typedef enum { CEED_BP1 = 0, CEED_BP2 = 1, CEED_BP3 = 2, CEED_BP4 = 3, CEED_BP5 = 4, CEED_BP6 = 5 } BPType;
1252b730f8bSJeremy L Thompson static const char *const bp_types[] = {"bp1", "bp2", "bp3", "bp4", "bp5", "bp6", "BPType", "CEED_BP", 0};
126cb32e2e7SValeria Barra
127cb32e2e7SValeria Barra // BP specific data
128cb32e2e7SValeria Barra typedef struct {
1299b072555Sjeremylt CeedInt num_comp_u, q_data_size, q_extra;
1309b072555Sjeremylt CeedQFunctionUser setup_geo, setup_rhs, apply, error;
1319b072555Sjeremylt const char *setup_geo_loc, *setup_rhs_loc, *apply_loc, *error_loc;
1329b072555Sjeremylt CeedEvalMode in_mode, out_mode;
1339b072555Sjeremylt CeedQuadMode q_mode;
1349b072555Sjeremylt } BPData;
135cb32e2e7SValeria Barra
1369b072555Sjeremylt BPData bp_options[6] = {
1372b730f8bSJeremy L Thompson [CEED_BP1] = {.num_comp_u = 1,
1389b072555Sjeremylt .q_data_size = 1,
1399b072555Sjeremylt .q_extra = 1,
1409b072555Sjeremylt .setup_geo = SetupMassGeo,
1419b072555Sjeremylt .setup_rhs = SetupMassRhs,
142cb32e2e7SValeria Barra .apply = Mass,
143cb32e2e7SValeria Barra .error = Error,
1449b072555Sjeremylt .setup_geo_loc = SetupMassGeo_loc,
1459b072555Sjeremylt .setup_rhs_loc = SetupMassRhs_loc,
1469b072555Sjeremylt .apply_loc = Mass_loc,
1479b072555Sjeremylt .error_loc = Error_loc,
1489b072555Sjeremylt .in_mode = CEED_EVAL_INTERP,
1499b072555Sjeremylt .out_mode = CEED_EVAL_INTERP,
1502b730f8bSJeremy L Thompson .q_mode = CEED_GAUSS },
1512b730f8bSJeremy L Thompson [CEED_BP2] = {.num_comp_u = 3,
1529b072555Sjeremylt .q_data_size = 1,
1539b072555Sjeremylt .q_extra = 1,
1549b072555Sjeremylt .setup_geo = SetupMassGeo,
1559b072555Sjeremylt .setup_rhs = SetupMassRhs3,
156cb32e2e7SValeria Barra .apply = Mass3,
157cb32e2e7SValeria Barra .error = Error3,
1589b072555Sjeremylt .setup_geo_loc = SetupMassGeo_loc,
1599b072555Sjeremylt .setup_rhs_loc = SetupMassRhs3_loc,
1609b072555Sjeremylt .apply_loc = Mass3_loc,
1619b072555Sjeremylt .error_loc = Error3_loc,
1629b072555Sjeremylt .in_mode = CEED_EVAL_INTERP,
1639b072555Sjeremylt .out_mode = CEED_EVAL_INTERP,
1642b730f8bSJeremy L Thompson .q_mode = CEED_GAUSS },
1652b730f8bSJeremy L Thompson [CEED_BP3] = {.num_comp_u = 1,
1669b072555Sjeremylt .q_data_size = 7,
1679b072555Sjeremylt .q_extra = 1,
1689b072555Sjeremylt .setup_geo = SetupDiffGeo,
1699b072555Sjeremylt .setup_rhs = SetupDiffRhs,
170cb32e2e7SValeria Barra .apply = Diff,
171cb32e2e7SValeria Barra .error = Error,
1729b072555Sjeremylt .setup_geo_loc = SetupDiffGeo_loc,
1739b072555Sjeremylt .setup_rhs_loc = SetupDiffRhs_loc,
1749b072555Sjeremylt .apply_loc = Diff_loc,
1759b072555Sjeremylt .error_loc = Error_loc,
1769b072555Sjeremylt .in_mode = CEED_EVAL_GRAD,
1779b072555Sjeremylt .out_mode = CEED_EVAL_GRAD,
1782b730f8bSJeremy L Thompson .q_mode = CEED_GAUSS },
1792b730f8bSJeremy L Thompson [CEED_BP4] = {.num_comp_u = 3,
1809b072555Sjeremylt .q_data_size = 7,
1819b072555Sjeremylt .q_extra = 1,
1829b072555Sjeremylt .setup_geo = SetupDiffGeo,
1839b072555Sjeremylt .setup_rhs = SetupDiffRhs3,
184cb32e2e7SValeria Barra .apply = Diff3,
185cb32e2e7SValeria Barra .error = Error3,
1869b072555Sjeremylt .setup_geo_loc = SetupDiffGeo_loc,
1879b072555Sjeremylt .setup_rhs_loc = SetupDiffRhs3_loc,
1889b072555Sjeremylt .apply_loc = Diff3_loc,
1899b072555Sjeremylt .error_loc = Error3_loc,
1909b072555Sjeremylt .in_mode = CEED_EVAL_GRAD,
1919b072555Sjeremylt .out_mode = CEED_EVAL_GRAD,
1922b730f8bSJeremy L Thompson .q_mode = CEED_GAUSS },
1932b730f8bSJeremy L Thompson [CEED_BP5] = {.num_comp_u = 1,
1949b072555Sjeremylt .q_data_size = 7,
1959b072555Sjeremylt .q_extra = 0,
1969b072555Sjeremylt .setup_geo = SetupDiffGeo,
1979b072555Sjeremylt .setup_rhs = SetupDiffRhs,
198cb32e2e7SValeria Barra .apply = Diff,
199cb32e2e7SValeria Barra .error = Error,
2009b072555Sjeremylt .setup_geo_loc = SetupDiffGeo_loc,
2019b072555Sjeremylt .setup_rhs_loc = SetupDiffRhs_loc,
2029b072555Sjeremylt .apply_loc = Diff_loc,
2039b072555Sjeremylt .error_loc = Error_loc,
2049b072555Sjeremylt .in_mode = CEED_EVAL_GRAD,
2059b072555Sjeremylt .out_mode = CEED_EVAL_GRAD,
2062b730f8bSJeremy L Thompson .q_mode = CEED_GAUSS_LOBATTO},
2072b730f8bSJeremy L Thompson [CEED_BP6] = {.num_comp_u = 3,
2089b072555Sjeremylt .q_data_size = 7,
2099b072555Sjeremylt .q_extra = 0,
2109b072555Sjeremylt .setup_geo = SetupDiffGeo,
2119b072555Sjeremylt .setup_rhs = SetupDiffRhs3,
212cb32e2e7SValeria Barra .apply = Diff3,
213cb32e2e7SValeria Barra .error = Error3,
2149b072555Sjeremylt .setup_geo_loc = SetupDiffGeo_loc,
2159b072555Sjeremylt .setup_rhs_loc = SetupDiffRhs3_loc,
2169b072555Sjeremylt .apply_loc = Diff3_loc,
2179b072555Sjeremylt .error_loc = Error3_loc,
2189b072555Sjeremylt .in_mode = CEED_EVAL_GRAD,
2199b072555Sjeremylt .out_mode = CEED_EVAL_GRAD,
2202b730f8bSJeremy L Thompson .q_mode = CEED_GAUSS_LOBATTO}
221cb32e2e7SValeria Barra };
222cb32e2e7SValeria Barra
223cb32e2e7SValeria Barra // This function uses libCEED to compute the action of the mass matrix
MatMult_Mass(Mat A,Vec X,Vec Y)224cb32e2e7SValeria Barra static PetscErrorCode MatMult_Mass(Mat A, Vec X, Vec Y) {
225d4d45553Srezgarshakeri OperatorApplyContext op_apply_ctx;
226cb32e2e7SValeria Barra PetscScalar *x, *y;
2279b072555Sjeremylt PetscMemType x_mem_type, y_mem_type;
228cb32e2e7SValeria Barra
229cb32e2e7SValeria Barra PetscFunctionBeginUser;
2302b730f8bSJeremy L Thompson PetscCall(MatShellGetContext(A, &op_apply_ctx));
2319396343dSjeremylt
2329396343dSjeremylt // Global-to-local
2332b730f8bSJeremy L Thompson PetscCall(VecScatterBegin(op_apply_ctx->l_to_g, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REVERSE));
2342b730f8bSJeremy L Thompson PetscCall(VecScatterEnd(op_apply_ctx->l_to_g, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REVERSE));
235cb32e2e7SValeria Barra
2369396343dSjeremylt // Setup libCEED vectors
2372b730f8bSJeremy L Thompson PetscCall(VecGetArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x, &x_mem_type));
2382b730f8bSJeremy L Thompson PetscCall(VecGetArrayAndMemType(op_apply_ctx->Y_loc, &y, &y_mem_type));
2392b730f8bSJeremy L Thompson CeedVectorSetArray(op_apply_ctx->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x);
2402b730f8bSJeremy L Thompson CeedVectorSetArray(op_apply_ctx->y_ceed, MemTypeP2C(y_mem_type), CEED_USE_POINTER, y);
241cb32e2e7SValeria Barra
2429396343dSjeremylt // Apply libCEED operator
2432b730f8bSJeremy L Thompson CeedOperatorApply(op_apply_ctx->op, op_apply_ctx->x_ceed, op_apply_ctx->y_ceed, CEED_REQUEST_IMMEDIATE);
244cb32e2e7SValeria Barra
2459396343dSjeremylt // Restore PETSc vectors
246d4d45553Srezgarshakeri CeedVectorTakeArray(op_apply_ctx->x_ceed, MemTypeP2C(x_mem_type), NULL);
247d4d45553Srezgarshakeri CeedVectorTakeArray(op_apply_ctx->y_ceed, MemTypeP2C(y_mem_type), NULL);
2482b730f8bSJeremy L Thompson PetscCall(VecRestoreArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x));
2492b730f8bSJeremy L Thompson PetscCall(VecRestoreArrayAndMemType(op_apply_ctx->Y_loc, &y));
250cb32e2e7SValeria Barra
2519396343dSjeremylt // Local-to-global
252cb32e2e7SValeria Barra if (Y) {
2532b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(Y));
2542b730f8bSJeremy L Thompson PetscCall(VecScatterBegin(op_apply_ctx->l_to_g, op_apply_ctx->Y_loc, Y, ADD_VALUES, SCATTER_FORWARD));
2552b730f8bSJeremy L Thompson PetscCall(VecScatterEnd(op_apply_ctx->l_to_g, op_apply_ctx->Y_loc, Y, ADD_VALUES, SCATTER_FORWARD));
256cb32e2e7SValeria Barra }
257ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
258cb32e2e7SValeria Barra }
259cb32e2e7SValeria Barra
260ea61e9acSJeremy L Thompson // This function uses libCEED to compute the action of the Laplacian with Dirichlet boundary conditions
MatMult_Diff(Mat A,Vec X,Vec Y)261cb32e2e7SValeria Barra static PetscErrorCode MatMult_Diff(Mat A, Vec X, Vec Y) {
262d4d45553Srezgarshakeri OperatorApplyContext op_apply_ctx;
263cb32e2e7SValeria Barra PetscScalar *x, *y;
2649b072555Sjeremylt PetscMemType x_mem_type, y_mem_type;
265cb32e2e7SValeria Barra
266cb32e2e7SValeria Barra PetscFunctionBeginUser;
2672b730f8bSJeremy L Thompson PetscCall(MatShellGetContext(A, &op_apply_ctx));
268cb32e2e7SValeria Barra
269cb32e2e7SValeria Barra // Global-to-local
2702b730f8bSJeremy L Thompson PetscCall(VecScatterBegin(op_apply_ctx->l_to_g_0, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REVERSE));
2712b730f8bSJeremy L Thompson PetscCall(VecScatterEnd(op_apply_ctx->l_to_g_0, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REVERSE));
272cb32e2e7SValeria Barra
2739396343dSjeremylt // Setup libCEED vectors
2742b730f8bSJeremy L Thompson PetscCall(VecGetArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x, &x_mem_type));
2752b730f8bSJeremy L Thompson PetscCall(VecGetArrayAndMemType(op_apply_ctx->Y_loc, &y, &y_mem_type));
2762b730f8bSJeremy L Thompson CeedVectorSetArray(op_apply_ctx->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x);
2772b730f8bSJeremy L Thompson CeedVectorSetArray(op_apply_ctx->y_ceed, MemTypeP2C(y_mem_type), CEED_USE_POINTER, y);
278cb32e2e7SValeria Barra
2799396343dSjeremylt // Apply libCEED operator
2802b730f8bSJeremy L Thompson CeedOperatorApply(op_apply_ctx->op, op_apply_ctx->x_ceed, op_apply_ctx->y_ceed, CEED_REQUEST_IMMEDIATE);
281cb32e2e7SValeria Barra
282cb32e2e7SValeria Barra // Restore PETSc vectors
283d4d45553Srezgarshakeri CeedVectorTakeArray(op_apply_ctx->x_ceed, MemTypeP2C(x_mem_type), NULL);
284d4d45553Srezgarshakeri CeedVectorTakeArray(op_apply_ctx->y_ceed, MemTypeP2C(y_mem_type), NULL);
2852b730f8bSJeremy L Thompson PetscCall(VecRestoreArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x));
2862b730f8bSJeremy L Thompson PetscCall(VecRestoreArrayAndMemType(op_apply_ctx->Y_loc, &y));
287cb32e2e7SValeria Barra
288cb32e2e7SValeria Barra // Local-to-global
2892b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(Y));
2902b730f8bSJeremy L Thompson PetscCall(VecScatterBegin(op_apply_ctx->g_to_g_D, X, Y, INSERT_VALUES, SCATTER_FORWARD));
2912b730f8bSJeremy L Thompson PetscCall(VecScatterEnd(op_apply_ctx->g_to_g_D, X, Y, INSERT_VALUES, SCATTER_FORWARD));
2922b730f8bSJeremy L Thompson PetscCall(VecScatterBegin(op_apply_ctx->l_to_g_0, op_apply_ctx->Y_loc, Y, ADD_VALUES, SCATTER_FORWARD));
2932b730f8bSJeremy L Thompson PetscCall(VecScatterEnd(op_apply_ctx->l_to_g_0, op_apply_ctx->Y_loc, Y, ADD_VALUES, SCATTER_FORWARD));
294ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
295cb32e2e7SValeria Barra }
296cb32e2e7SValeria Barra
297cb32e2e7SValeria Barra // This function calculates the error in the final solution
ComputeErrorMax(OperatorApplyContext op_apply_ctx,CeedOperator op_error,Vec X,CeedVector target,PetscReal * max_error)2982b730f8bSJeremy L Thompson static PetscErrorCode ComputeErrorMax(OperatorApplyContext op_apply_ctx, CeedOperator op_error, Vec X, CeedVector target, PetscReal *max_error) {
299cb32e2e7SValeria Barra PetscScalar *x;
3009b072555Sjeremylt PetscMemType mem_type;
301cb32e2e7SValeria Barra CeedVector collocated_error;
3021f9221feSJeremy L Thompson CeedSize length;
303cb32e2e7SValeria Barra
304cb32e2e7SValeria Barra PetscFunctionBeginUser;
305cb32e2e7SValeria Barra CeedVectorGetLength(target, &length);
306d4d45553Srezgarshakeri CeedVectorCreate(op_apply_ctx->ceed, length, &collocated_error);
307cb32e2e7SValeria Barra
308cb32e2e7SValeria Barra // Global-to-local
3092b730f8bSJeremy L Thompson PetscCall(VecScatterBegin(op_apply_ctx->l_to_g, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REVERSE));
3102b730f8bSJeremy L Thompson PetscCall(VecScatterEnd(op_apply_ctx->l_to_g, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REVERSE));
3119396343dSjeremylt
3129396343dSjeremylt // Setup libCEED vector
3132b730f8bSJeremy L Thompson PetscCall(VecGetArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x, &mem_type));
3142b730f8bSJeremy L Thompson CeedVectorSetArray(op_apply_ctx->x_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, x);
315cb32e2e7SValeria Barra
3169396343dSjeremylt // Apply libCEED operator
3172b730f8bSJeremy L Thompson CeedOperatorApply(op_error, op_apply_ctx->x_ceed, collocated_error, CEED_REQUEST_IMMEDIATE);
318cb32e2e7SValeria Barra
319cb32e2e7SValeria Barra // Restore PETSc vector
320d4d45553Srezgarshakeri CeedVectorTakeArray(op_apply_ctx->x_ceed, MemTypeP2C(mem_type), NULL);
3212b730f8bSJeremy L Thompson PetscCall(VecRestoreArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x));
322cb32e2e7SValeria Barra
323cb32e2e7SValeria Barra // Reduce max error
3249b072555Sjeremylt *max_error = 0;
325cb32e2e7SValeria Barra const CeedScalar *e;
326cb32e2e7SValeria Barra CeedVectorGetArrayRead(collocated_error, CEED_MEM_HOST, &e);
327cb32e2e7SValeria Barra for (CeedInt i = 0; i < length; i++) {
3289b072555Sjeremylt *max_error = PetscMax(*max_error, PetscAbsScalar(e[i]));
329cb32e2e7SValeria Barra }
330cb32e2e7SValeria Barra CeedVectorRestoreArrayRead(collocated_error, &e);
3312b730f8bSJeremy L Thompson PetscCall(MPI_Allreduce(MPI_IN_PLACE, max_error, 1, MPIU_REAL, MPIU_MAX, op_apply_ctx->comm));
332cb32e2e7SValeria Barra
333cb32e2e7SValeria Barra // Cleanup
334cb32e2e7SValeria Barra CeedVectorDestroy(&collocated_error);
335ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
336cb32e2e7SValeria Barra }
337cb32e2e7SValeria Barra
main(int argc,char ** argv)338cb32e2e7SValeria Barra int main(int argc, char **argv) {
339cb32e2e7SValeria Barra MPI_Comm comm;
3409b072555Sjeremylt char ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self";
341cb32e2e7SValeria Barra double my_rt_start, my_rt, rt_min, rt_max;
3422b730f8bSJeremy 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,
3432b730f8bSJeremy L Thompson ksp_max_it_clip[2];
344cb32e2e7SValeria Barra PetscScalar *r;
345cb32e2e7SValeria Barra PetscBool test_mode, benchmark_mode, write_solution;
346cb32e2e7SValeria Barra PetscMPIInt size, rank;
3479b072555Sjeremylt PetscLogStage solve_stage;
3489b072555Sjeremylt Vec X, X_loc, rhs, rhs_loc;
349cb32e2e7SValeria Barra Mat mat;
350cb32e2e7SValeria Barra KSP ksp;
3519b072555Sjeremylt VecScatter l_to_g, l_to_g_0, g_to_g_D;
3529b072555Sjeremylt PetscMemType mem_type;
353d4d45553Srezgarshakeri OperatorApplyContext op_apply_ctx;
354cb32e2e7SValeria Barra Ceed ceed;
3559b072555Sjeremylt CeedBasis basis_x, basis_u;
3569b072555Sjeremylt CeedElemRestriction elem_restr_x, elem_restr_u, elem_restr_u_i, elem_restr_qd_i;
3579b072555Sjeremylt CeedQFunction qf_setup_geo, qf_setup_rhs, qf_apply, qf_error;
3589b072555Sjeremylt CeedOperator op_setup_geo, op_setup_rhs, op_apply, op_error;
3599b072555Sjeremylt CeedVector x_coord, q_data, rhs_ceed, target;
360cb32e2e7SValeria Barra CeedInt P, Q;
3619b072555Sjeremylt const CeedInt dim = 3, num_comp_x = 3;
3629b072555Sjeremylt BPType bp_choice;
363cb32e2e7SValeria Barra
3642b730f8bSJeremy L Thompson PetscCall(PetscInitialize(&argc, &argv, NULL, help));
365cb32e2e7SValeria Barra comm = PETSC_COMM_WORLD;
36632d2ee49SValeria Barra
36732d2ee49SValeria Barra // Read command line options
36867490bc6SJeremy L Thompson PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
3699b072555Sjeremylt bp_choice = CEED_BP1;
3702b730f8bSJeremy L Thompson PetscCall(PetscOptionsEnum("-problem", "CEED benchmark problem to solve", NULL, bp_types, (PetscEnum)bp_choice, (PetscEnum *)&bp_choice, NULL));
3719b072555Sjeremylt num_comp_u = bp_options[bp_choice].num_comp_u;
372cb32e2e7SValeria Barra test_mode = PETSC_FALSE;
3732b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, test_mode, &test_mode, NULL));
374cb32e2e7SValeria Barra benchmark_mode = PETSC_FALSE;
3752b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-benchmark", "Benchmarking mode (prints benchmark statistics)", NULL, benchmark_mode, &benchmark_mode, NULL));
376cb32e2e7SValeria Barra write_solution = PETSC_FALSE;
3772b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-write_solution", "Write solution for visualization", NULL, write_solution, &write_solution, NULL));
378cb32e2e7SValeria Barra degree = test_mode ? 3 : 1;
3792b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", NULL, degree, °ree, NULL));
3809b072555Sjeremylt q_extra = bp_options[bp_choice].q_extra;
3812b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, q_extra, &q_extra, NULL));
3822b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, ceed_resource, ceed_resource, sizeof(ceed_resource), NULL));
3839b072555Sjeremylt local_nodes = 1000;
3842b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-local", "Target number of locally owned nodes per process", NULL, local_nodes, &local_nodes, NULL));
3852fbc6e41SJeremy L Thompson PetscInt two = 2;
3862fbc6e41SJeremy L Thompson ksp_max_it_clip[0] = 5;
3872fbc6e41SJeremy L Thompson ksp_max_it_clip[1] = 20;
3881a8516d0SJames Wright PetscCall(PetscOptionsIntArray("-ksp_max_it_clip", "Min and max number of iterations to use during benchmarking", NULL, ksp_max_it_clip, &two,
3891a8516d0SJames Wright NULL));
39067490bc6SJeremy L Thompson PetscOptionsEnd();
391cb32e2e7SValeria Barra P = degree + 1;
3929b072555Sjeremylt Q = P + q_extra;
393cb32e2e7SValeria Barra
3949396343dSjeremylt // Set up libCEED
3959b072555Sjeremylt CeedInit(ceed_resource, &ceed);
3969b072555Sjeremylt CeedMemType mem_type_backend;
3979b072555Sjeremylt CeedGetPreferredMemType(ceed, &mem_type_backend);
3989396343dSjeremylt
3999b072555Sjeremylt VecType default_vec_type = NULL, vec_type;
4009b072555Sjeremylt switch (mem_type_backend) {
4012b730f8bSJeremy L Thompson case CEED_MEM_HOST:
4022b730f8bSJeremy L Thompson default_vec_type = VECSTANDARD;
4032b730f8bSJeremy L Thompson break;
404b68a8d79SJed Brown case CEED_MEM_DEVICE: {
405b68a8d79SJed Brown const char *resolved;
4068c03e814SJeremy L Thompson
407b68a8d79SJed Brown CeedGetResource(ceed, &resolved);
4089b072555Sjeremylt if (strstr(resolved, "/gpu/cuda")) default_vec_type = VECCUDA;
4099b072555Sjeremylt else if (strstr(resolved, "/gpu/hip")) default_vec_type = VECHIP;
4109b072555Sjeremylt else default_vec_type = VECSTANDARD;
411b68a8d79SJed Brown }
412b68a8d79SJed Brown }
4139396343dSjeremylt
414cb32e2e7SValeria Barra // Determine size of process grid
4152b730f8bSJeremy L Thompson PetscCall(MPI_Comm_size(comm, &size));
416cb32e2e7SValeria Barra Split3(size, p, false);
417cb32e2e7SValeria Barra
4189b072555Sjeremylt // Find a nicely composite number of elements no less than local_nodes
4192b730f8bSJeremy L Thompson for (local_elem = PetscMax(1, local_nodes / (degree * degree * degree));; local_elem++) {
4209b072555Sjeremylt Split3(local_elem, mesh_elem, true);
4219b072555Sjeremylt if (Max3(mesh_elem) / Min3(mesh_elem) <= 2) break;
422cb32e2e7SValeria Barra }
423cb32e2e7SValeria Barra
424cb32e2e7SValeria Barra // Find my location in the process grid
4252b730f8bSJeremy L Thompson PetscCall(MPI_Comm_rank(comm, &rank));
4269b072555Sjeremylt for (int d = 0, rank_left = rank; d < dim; d++) {
427cb32e2e7SValeria Barra const int pstride[3] = {p[1] * p[2], p[2], 1};
4289b072555Sjeremylt i_rank[d] = rank_left / pstride[d];
4299b072555Sjeremylt rank_left -= i_rank[d] * pstride[d];
430cb32e2e7SValeria Barra }
431cb32e2e7SValeria Barra
4329b072555Sjeremylt GlobalNodes(p, i_rank, degree, mesh_elem, m_nodes);
433cb32e2e7SValeria Barra
434cb32e2e7SValeria Barra // Setup global vector
4352b730f8bSJeremy L Thompson PetscCall(VecCreate(comm, &X));
4362b730f8bSJeremy L Thompson PetscCall(VecSetType(X, default_vec_type));
4372b730f8bSJeremy L Thompson PetscCall(VecSetSizes(X, m_nodes[0] * m_nodes[1] * m_nodes[2] * num_comp_u, PETSC_DECIDE));
4382b730f8bSJeremy L Thompson PetscCall(VecSetFromOptions(X));
4392b730f8bSJeremy L Thompson PetscCall(VecSetUp(X));
440cb32e2e7SValeria Barra
441cb32e2e7SValeria Barra // Print summary
4424d00b080SJeremy L Thompson PetscInt gsize;
4434d00b080SJeremy L Thompson
4442b730f8bSJeremy L Thompson PetscCall(VecGetSize(X, &gsize));
445dc7d240cSValeria Barra if (!test_mode) {
4469b072555Sjeremylt const char *used_resource;
4479b072555Sjeremylt CeedGetResource(ceed, &used_resource);
4489396343dSjeremylt
4492b730f8bSJeremy L Thompson PetscCall(VecGetType(X, &vec_type));
4509396343dSjeremylt
4512b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm,
452990fdeb6SJeremy L Thompson "\n-- CEED Benchmark Problem %" CeedInt_FMT " -- libCEED + PETSc --\n"
4539396343dSjeremylt " PETSc:\n"
4549396343dSjeremylt " PETSc Vec Type : %s\n"
455cb32e2e7SValeria Barra " libCEED:\n"
456cb32e2e7SValeria Barra " libCEED Backend : %s\n"
4579396343dSjeremylt " libCEED Backend MemType : %s\n"
458cb32e2e7SValeria Barra " Mesh:\n"
459751eb813Srezgarshakeri " Solution Order (P) : %" CeedInt_FMT "\n"
460751eb813Srezgarshakeri " Quadrature Order (Q) : %" CeedInt_FMT "\n"
46108140895SJed Brown " Global nodes : %" PetscInt_FMT "\n"
4622b730f8bSJeremy L Thompson " Process Decomposition : %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n"
4632b730f8bSJeremy L Thompson " Local Elements : %" PetscInt_FMT " = %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n"
4642b730f8bSJeremy L Thompson " Owned nodes : %" PetscInt_FMT " = %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n"
46508140895SJed Brown " DoF per node : %" PetscInt_FMT "\n",
4662b730f8bSJeremy 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],
4672b730f8bSJeremy 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],
4682b730f8bSJeremy L Thompson m_nodes[2], num_comp_u));
469cb32e2e7SValeria Barra }
470cb32e2e7SValeria Barra
471cb32e2e7SValeria Barra {
4729b072555Sjeremylt l_size = 1;
473cb32e2e7SValeria Barra for (int d = 0; d < dim; d++) {
4749b072555Sjeremylt l_nodes[d] = mesh_elem[d] * degree + 1;
4759b072555Sjeremylt l_size *= l_nodes[d];
476cb32e2e7SValeria Barra }
4772b730f8bSJeremy L Thompson PetscCall(VecCreate(PETSC_COMM_SELF, &X_loc));
4782b730f8bSJeremy L Thompson PetscCall(VecSetType(X_loc, default_vec_type));
4792b730f8bSJeremy L Thompson PetscCall(VecSetSizes(X_loc, l_size * num_comp_u, PETSC_DECIDE));
4802b730f8bSJeremy L Thompson PetscCall(VecSetFromOptions(X_loc));
4812b730f8bSJeremy L Thompson PetscCall(VecSetUp(X_loc));
482cb32e2e7SValeria Barra
483cb32e2e7SValeria Barra // Create local-to-global scatter
4849b072555Sjeremylt PetscInt *l_to_g_ind, *l_to_g_ind_0, *loc_ind, l_0_count;
4859b072555Sjeremylt IS l_to_g_is, l_to_g_is_0, loc_is;
4869b072555Sjeremylt PetscInt g_start[2][2][2], g_m_nodes[2][2][2][dim];
487cb32e2e7SValeria Barra
488cb32e2e7SValeria Barra for (int i = 0; i < 2; i++) {
489cb32e2e7SValeria Barra for (int j = 0; j < 2; j++) {
490cb32e2e7SValeria Barra for (int k = 0; k < 2; k++) {
4919b072555Sjeremylt PetscInt ijk_rank[3] = {i_rank[0] + i, i_rank[1] + j, i_rank[2] + k};
4929b072555Sjeremylt g_start[i][j][k] = GlobalStart(p, ijk_rank, degree, mesh_elem);
4939b072555Sjeremylt GlobalNodes(p, ijk_rank, degree, mesh_elem, g_m_nodes[i][j][k]);
494cb32e2e7SValeria Barra }
495cb32e2e7SValeria Barra }
496cb32e2e7SValeria Barra }
497cb32e2e7SValeria Barra
4982b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(l_size, &l_to_g_ind));
4992b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(l_size, &l_to_g_ind_0));
5002b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(l_size, &loc_ind));
5019b072555Sjeremylt l_0_count = 0;
5022b730f8bSJeremy L Thompson for (PetscInt i = 0, ir, ii; ir = i >= m_nodes[0], ii = i - ir * m_nodes[0], i < l_nodes[0]; i++) {
5032b730f8bSJeremy L Thompson for (PetscInt j = 0, jr, jj; jr = j >= m_nodes[1], jj = j - jr * m_nodes[1], j < l_nodes[1]; j++) {
5042b730f8bSJeremy L Thompson for (PetscInt k = 0, kr, kk; kr = k >= m_nodes[2], kk = k - kr * m_nodes[2], k < l_nodes[2]; k++) {
5059b072555Sjeremylt PetscInt here = (i * l_nodes[1] + j) * l_nodes[2] + k;
5062b730f8bSJeremy 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;
5072b730f8bSJeremy L Thompson if ((i_rank[0] == 0 && i == 0) || (i_rank[1] == 0 && j == 0) || (i_rank[2] == 0 && k == 0) ||
5082b730f8bSJeremy L Thompson (i_rank[0] + 1 == p[0] && i + 1 == l_nodes[0]) || (i_rank[1] + 1 == p[1] && j + 1 == l_nodes[1]) ||
5096c10af5dSJeremy L Thompson (i_rank[2] + 1 == p[2] && k + 1 == l_nodes[2])) {
510cb32e2e7SValeria Barra continue;
5116c10af5dSJeremy L Thompson }
5129b072555Sjeremylt l_to_g_ind_0[l_0_count] = l_to_g_ind[here];
5139b072555Sjeremylt loc_ind[l_0_count++] = here;
514cb32e2e7SValeria Barra }
5152b730f8bSJeremy L Thompson }
5162b730f8bSJeremy L Thompson }
5172b730f8bSJeremy L Thompson PetscCall(ISCreateBlock(comm, num_comp_u, l_size, l_to_g_ind, PETSC_OWN_POINTER, &l_to_g_is));
5182b730f8bSJeremy L Thompson PetscCall(VecScatterCreate(X_loc, NULL, X, l_to_g_is, &l_to_g));
5192b730f8bSJeremy L Thompson PetscCall(ISCreateBlock(comm, num_comp_u, l_0_count, l_to_g_ind_0, PETSC_OWN_POINTER, &l_to_g_is_0));
5202b730f8bSJeremy L Thompson PetscCall(ISCreateBlock(comm, num_comp_u, l_0_count, loc_ind, PETSC_OWN_POINTER, &loc_is));
5212b730f8bSJeremy L Thompson PetscCall(VecScatterCreate(X_loc, loc_is, X, l_to_g_is_0, &l_to_g_0));
522cb32e2e7SValeria Barra {
523ea61e9acSJeremy 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)
5249b072555Sjeremylt PetscInt x_start, x_end, *ind_D, count_D = 0;
5259b072555Sjeremylt IS is_D;
526cb32e2e7SValeria Barra const PetscScalar *x;
5272b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(X_loc));
5282b730f8bSJeremy L Thompson PetscCall(VecSet(X, 1.0));
5292b730f8bSJeremy L Thompson PetscCall(VecScatterBegin(l_to_g_0, X_loc, X, INSERT_VALUES, SCATTER_FORWARD));
5302b730f8bSJeremy L Thompson PetscCall(VecScatterEnd(l_to_g_0, X_loc, X, INSERT_VALUES, SCATTER_FORWARD));
5312b730f8bSJeremy L Thompson PetscCall(VecGetOwnershipRange(X, &x_start, &x_end));
5322b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(x_end - x_start, &ind_D));
5332b730f8bSJeremy L Thompson PetscCall(VecGetArrayRead(X, &x));
5349b072555Sjeremylt for (PetscInt i = 0; i < x_end - x_start; i++) {
5359b072555Sjeremylt if (x[i] == 1.) ind_D[count_D++] = x_start + i;
536cb32e2e7SValeria Barra }
5372b730f8bSJeremy L Thompson PetscCall(VecRestoreArrayRead(X, &x));
5382b730f8bSJeremy L Thompson PetscCall(ISCreateGeneral(comm, count_D, ind_D, PETSC_COPY_VALUES, &is_D));
5392b730f8bSJeremy L Thompson PetscCall(PetscFree(ind_D));
5402b730f8bSJeremy L Thompson PetscCall(VecScatterCreate(X, is_D, X, is_D, &g_to_g_D));
5412b730f8bSJeremy L Thompson PetscCall(ISDestroy(&is_D));
542cb32e2e7SValeria Barra }
5432b730f8bSJeremy L Thompson PetscCall(ISDestroy(&l_to_g_is));
5442b730f8bSJeremy L Thompson PetscCall(ISDestroy(&l_to_g_is_0));
5452b730f8bSJeremy L Thompson PetscCall(ISDestroy(&loc_is));
546cb32e2e7SValeria Barra }
547cb32e2e7SValeria Barra
548cb32e2e7SValeria Barra // CEED bases
5492b730f8bSJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_u, P, Q, bp_options[bp_choice].q_mode, &basis_u);
5502b730f8bSJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, bp_options[bp_choice].q_mode, &basis_x);
551cb32e2e7SValeria Barra
552cb32e2e7SValeria Barra // CEED restrictions
553179e5961SZach Atkins PetscCall(CreateRestriction(ceed, mesh_elem, P, num_comp_u, &elem_restr_u));
554179e5961SZach Atkins PetscCall(CreateRestriction(ceed, mesh_elem, 2, dim, &elem_restr_x));
5559b072555Sjeremylt CeedInt num_elem = mesh_elem[0] * mesh_elem[1] * mesh_elem[2];
5564d00b080SJeremy L Thompson
5572b730f8bSJeremy 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);
5582b730f8bSJeremy L Thompson CeedElemRestrictionCreateStrided(ceed, num_elem, Q * Q * Q, bp_options[bp_choice].q_data_size,
5592b730f8bSJeremy L Thompson bp_options[bp_choice].q_data_size * num_elem * Q * Q * Q, CEED_STRIDES_BACKEND, &elem_restr_qd_i);
560cb32e2e7SValeria Barra {
5619b072555Sjeremylt CeedScalar *x_loc;
5622b730f8bSJeremy L Thompson CeedInt shape[3] = {mesh_elem[0] + 1, mesh_elem[1] + 1, mesh_elem[2] + 1}, len = shape[0] * shape[1] * shape[2];
5634d00b080SJeremy L Thompson
5649b072555Sjeremylt x_loc = malloc(len * num_comp_x * sizeof x_loc[0]);
565cb32e2e7SValeria Barra for (CeedInt i = 0; i < shape[0]; i++) {
566cb32e2e7SValeria Barra for (CeedInt j = 0; j < shape[1]; j++) {
567cb32e2e7SValeria Barra for (CeedInt k = 0; k < shape[2]; k++) {
5682b730f8bSJeremy 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]);
5692b730f8bSJeremy 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]);
5702b730f8bSJeremy 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]);
571cb32e2e7SValeria Barra }
572cb32e2e7SValeria Barra }
573cb32e2e7SValeria Barra }
5749b072555Sjeremylt CeedVectorCreate(ceed, len * num_comp_x, &x_coord);
5759b072555Sjeremylt CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_OWN_POINTER, x_loc);
576cb32e2e7SValeria Barra }
577cb32e2e7SValeria Barra
5789b072555Sjeremylt // Create the QFunction that builds the operator quadrature data
5792b730f8bSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].setup_geo, bp_options[bp_choice].setup_geo_loc, &qf_setup_geo);
5809b072555Sjeremylt CeedQFunctionAddInput(qf_setup_geo, "x", num_comp_x, CEED_EVAL_INTERP);
5819b072555Sjeremylt CeedQFunctionAddInput(qf_setup_geo, "dx", num_comp_x * dim, CEED_EVAL_GRAD);
5829b072555Sjeremylt CeedQFunctionAddInput(qf_setup_geo, "weight", 1, CEED_EVAL_WEIGHT);
5832b730f8bSJeremy L Thompson CeedQFunctionAddOutput(qf_setup_geo, "q_data", bp_options[bp_choice].q_data_size, CEED_EVAL_NONE);
584cb32e2e7SValeria Barra
5859b072555Sjeremylt // Create the QFunction that sets up the RHS and true solution
5862b730f8bSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].setup_rhs, bp_options[bp_choice].setup_rhs_loc, &qf_setup_rhs);
5879b072555Sjeremylt CeedQFunctionAddInput(qf_setup_rhs, "x", num_comp_x, CEED_EVAL_INTERP);
5882b730f8bSJeremy L Thompson CeedQFunctionAddInput(qf_setup_rhs, "q_data", bp_options[bp_choice].q_data_size, CEED_EVAL_NONE);
5899b072555Sjeremylt CeedQFunctionAddOutput(qf_setup_rhs, "true_soln", num_comp_u, CEED_EVAL_NONE);
5909b072555Sjeremylt CeedQFunctionAddOutput(qf_setup_rhs, "rhs", num_comp_u, CEED_EVAL_INTERP);
591cb32e2e7SValeria Barra
592cb32e2e7SValeria Barra // Set up PDE operator
5932b730f8bSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].apply, bp_options[bp_choice].apply_loc, &qf_apply);
594cb32e2e7SValeria Barra // Add inputs and outputs
5959b072555Sjeremylt CeedInt in_scale = bp_options[bp_choice].in_mode == CEED_EVAL_GRAD ? 3 : 1;
5969b072555Sjeremylt CeedInt out_scale = bp_options[bp_choice].out_mode == CEED_EVAL_GRAD ? 3 : 1;
5972b730f8bSJeremy L Thompson CeedQFunctionAddInput(qf_apply, "u", num_comp_u * in_scale, bp_options[bp_choice].in_mode);
5982b730f8bSJeremy L Thompson CeedQFunctionAddInput(qf_apply, "q_data", bp_options[bp_choice].q_data_size, CEED_EVAL_NONE);
5992b730f8bSJeremy L Thompson CeedQFunctionAddOutput(qf_apply, "v", num_comp_u * out_scale, bp_options[bp_choice].out_mode);
600cb32e2e7SValeria Barra
601cb32e2e7SValeria Barra // Create the error qfunction
6022b730f8bSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].error, bp_options[bp_choice].error_loc, &qf_error);
6039b072555Sjeremylt CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP);
6049b072555Sjeremylt CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE);
6052b730f8bSJeremy L Thompson CeedQFunctionAddInput(qf_error, "qdata", bp_options[bp_choice].q_data_size, CEED_EVAL_NONE);
6069b072555Sjeremylt CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_NONE);
607cb32e2e7SValeria Barra
608cb32e2e7SValeria Barra // Create the persistent vectors that will be needed in setup
6099b072555Sjeremylt CeedInt num_qpts;
6109b072555Sjeremylt CeedBasisGetNumQuadraturePoints(basis_u, &num_qpts);
6112b730f8bSJeremy L Thompson CeedVectorCreate(ceed, bp_options[bp_choice].q_data_size * num_elem * num_qpts, &q_data);
6129b072555Sjeremylt CeedVectorCreate(ceed, num_elem * num_qpts * num_comp_u, &target);
6139b072555Sjeremylt CeedVectorCreate(ceed, l_size * num_comp_u, &rhs_ceed);
614cb32e2e7SValeria Barra
615cb32e2e7SValeria Barra // Create the operator that builds the quadrature data for the ceed operator
6162b730f8bSJeremy L Thompson CeedOperatorCreate(ceed, qf_setup_geo, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setup_geo);
6172b730f8bSJeremy L Thompson CeedOperatorSetField(op_setup_geo, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
6182b730f8bSJeremy L Thompson CeedOperatorSetField(op_setup_geo, "dx", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
6192b730f8bSJeremy L Thompson CeedOperatorSetField(op_setup_geo, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE);
620356036faSJeremy L Thompson CeedOperatorSetField(op_setup_geo, "q_data", elem_restr_qd_i, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
621cb32e2e7SValeria Barra
622cb32e2e7SValeria Barra // Create the operator that builds the RHS and true solution
6232b730f8bSJeremy L Thompson CeedOperatorCreate(ceed, qf_setup_rhs, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setup_rhs);
6242b730f8bSJeremy L Thompson CeedOperatorSetField(op_setup_rhs, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
625356036faSJeremy L Thompson CeedOperatorSetField(op_setup_rhs, "q_data", elem_restr_qd_i, CEED_BASIS_NONE, q_data);
626356036faSJeremy L Thompson CeedOperatorSetField(op_setup_rhs, "true_soln", elem_restr_u_i, CEED_BASIS_NONE, target);
6272b730f8bSJeremy L Thompson CeedOperatorSetField(op_setup_rhs, "rhs", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
628cb32e2e7SValeria Barra
629cb32e2e7SValeria Barra // Create the mass or diff operator
6302b730f8bSJeremy L Thompson CeedOperatorCreate(ceed, qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_apply);
6319b072555Sjeremylt CeedOperatorSetField(op_apply, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
632356036faSJeremy L Thompson CeedOperatorSetField(op_apply, "q_data", elem_restr_qd_i, CEED_BASIS_NONE, q_data);
6339b072555Sjeremylt CeedOperatorSetField(op_apply, "v", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
634cb32e2e7SValeria Barra
635cb32e2e7SValeria Barra // Create the error operator
6362b730f8bSJeremy L Thompson CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_error);
6379b072555Sjeremylt CeedOperatorSetField(op_error, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
638356036faSJeremy L Thompson CeedOperatorSetField(op_error, "true_soln", elem_restr_u_i, CEED_BASIS_NONE, target);
639356036faSJeremy L Thompson CeedOperatorSetField(op_error, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data);
640356036faSJeremy L Thompson CeedOperatorSetField(op_error, "error", elem_restr_u_i, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
641cb32e2e7SValeria Barra
642cb32e2e7SValeria Barra // Set up Mat
6432b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &op_apply_ctx));
644d4d45553Srezgarshakeri op_apply_ctx->comm = comm;
645d4d45553Srezgarshakeri op_apply_ctx->l_to_g = l_to_g;
6469b072555Sjeremylt if (bp_choice != CEED_BP1 && bp_choice != CEED_BP2) {
647d4d45553Srezgarshakeri op_apply_ctx->l_to_g_0 = l_to_g_0;
648d4d45553Srezgarshakeri op_apply_ctx->g_to_g_D = g_to_g_D;
649cb32e2e7SValeria Barra }
650d4d45553Srezgarshakeri op_apply_ctx->X_loc = X_loc;
6512b730f8bSJeremy L Thompson PetscCall(VecDuplicate(X_loc, &op_apply_ctx->Y_loc));
652d4d45553Srezgarshakeri CeedVectorCreate(ceed, l_size * num_comp_u, &op_apply_ctx->x_ceed);
653d4d45553Srezgarshakeri CeedVectorCreate(ceed, l_size * num_comp_u, &op_apply_ctx->y_ceed);
654d4d45553Srezgarshakeri op_apply_ctx->op = op_apply;
655d4d45553Srezgarshakeri op_apply_ctx->q_data = q_data;
656d4d45553Srezgarshakeri op_apply_ctx->ceed = ceed;
657cb32e2e7SValeria Barra
6582b730f8bSJeremy 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,
6592b730f8bSJeremy L Thompson PETSC_DECIDE, op_apply_ctx, &mat));
6609b072555Sjeremylt if (bp_choice == CEED_BP1 || bp_choice == CEED_BP2) {
6612b730f8bSJeremy L Thompson PetscCall(MatShellSetOperation(mat, MATOP_MULT, (void (*)(void))MatMult_Mass));
662cb32e2e7SValeria Barra } else {
6632b730f8bSJeremy L Thompson PetscCall(MatShellSetOperation(mat, MATOP_MULT, (void (*)(void))MatMult_Diff));
664cb32e2e7SValeria Barra }
6652b730f8bSJeremy L Thompson PetscCall(VecGetType(X, &vec_type));
6662b730f8bSJeremy L Thompson PetscCall(MatShellSetVecType(mat, vec_type));
667cb32e2e7SValeria Barra
668cb32e2e7SValeria Barra // Get RHS vector
6692b730f8bSJeremy L Thompson PetscCall(VecDuplicate(X, &rhs));
6702b730f8bSJeremy L Thompson PetscCall(VecDuplicate(X_loc, &rhs_loc));
6712b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(rhs_loc));
6722b730f8bSJeremy L Thompson PetscCall(VecGetArrayAndMemType(rhs_loc, &r, &mem_type));
6739b072555Sjeremylt CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r);
674cb32e2e7SValeria Barra
6759b072555Sjeremylt // Setup q_data, rhs, and target
6769b072555Sjeremylt CeedOperatorApply(op_setup_geo, x_coord, q_data, CEED_REQUEST_IMMEDIATE);
6779b072555Sjeremylt CeedOperatorApply(op_setup_rhs, x_coord, rhs_ceed, CEED_REQUEST_IMMEDIATE);
6789b072555Sjeremylt CeedVectorDestroy(&x_coord);
679cb32e2e7SValeria Barra
680cb32e2e7SValeria Barra // Gather RHS
6812b730f8bSJeremy L Thompson PetscCall(CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL));
6822b730f8bSJeremy L Thompson PetscCall(VecRestoreArrayAndMemType(rhs_loc, &r));
6832b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(rhs));
6842b730f8bSJeremy L Thompson PetscCall(VecScatterBegin(l_to_g, rhs_loc, rhs, ADD_VALUES, SCATTER_FORWARD));
6852b730f8bSJeremy L Thompson PetscCall(VecScatterEnd(l_to_g, rhs_loc, rhs, ADD_VALUES, SCATTER_FORWARD));
6869b072555Sjeremylt CeedVectorDestroy(&rhs_ceed);
687cb32e2e7SValeria Barra
6882b730f8bSJeremy L Thompson PetscCall(KSPCreate(comm, &ksp));
689cb32e2e7SValeria Barra {
690cb32e2e7SValeria Barra PC pc;
6912b730f8bSJeremy L Thompson PetscCall(KSPGetPC(ksp, &pc));
6929b072555Sjeremylt if (bp_choice == CEED_BP1 || bp_choice == CEED_BP2) {
6932b730f8bSJeremy L Thompson PetscCall(PCSetType(pc, PCJACOBI));
6942b730f8bSJeremy L Thompson PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
695cb32e2e7SValeria Barra } else {
6962b730f8bSJeremy L Thompson PetscCall(PCSetType(pc, PCNONE));
697cb32e2e7SValeria Barra }
6982b730f8bSJeremy L Thompson PetscCall(KSPSetType(ksp, KSPCG));
6992b730f8bSJeremy L Thompson PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
7002b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
701cb32e2e7SValeria Barra }
7022b730f8bSJeremy L Thompson PetscCall(KSPSetOperators(ksp, mat, mat));
703da9108adSvaleriabarra // First run's performance log is not considered for benchmarking purposes
7042b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1));
705cb32e2e7SValeria Barra my_rt_start = MPI_Wtime();
7062b730f8bSJeremy L Thompson PetscCall(KSPSolve(ksp, rhs, X));
707cb32e2e7SValeria Barra my_rt = MPI_Wtime() - my_rt_start;
7082b730f8bSJeremy L Thompson PetscCall(MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm));
709cb32e2e7SValeria Barra // Set maxits based on first iteration timing
710cb32e2e7SValeria Barra if (my_rt > 0.02) {
7112b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, ksp_max_it_clip[0]));
712cb32e2e7SValeria Barra } else {
7132b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, ksp_max_it_clip[1]));
714cb32e2e7SValeria Barra }
7152b730f8bSJeremy L Thompson PetscCall(KSPSetFromOptions(ksp));
71609a940d7Sjeremylt
717cb32e2e7SValeria Barra // Timed solve
7182b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(X));
7192b730f8bSJeremy L Thompson PetscCall(PetscBarrier((PetscObject)ksp));
72009a940d7Sjeremylt
72109a940d7Sjeremylt // -- Performance logging
7222b730f8bSJeremy L Thompson PetscCall(PetscLogStageRegister("Solve Stage", &solve_stage));
7232b730f8bSJeremy L Thompson PetscCall(PetscLogStagePush(solve_stage));
72409a940d7Sjeremylt
72509a940d7Sjeremylt // -- Solve
726cb32e2e7SValeria Barra my_rt_start = MPI_Wtime();
7272b730f8bSJeremy L Thompson PetscCall(KSPSolve(ksp, rhs, X));
728cb32e2e7SValeria Barra my_rt = MPI_Wtime() - my_rt_start;
72909a940d7Sjeremylt
73009a940d7Sjeremylt // -- Performance logging
7312b730f8bSJeremy L Thompson PetscCall(PetscLogStagePop());
73209a940d7Sjeremylt
7338e87e98bSJed Brown // Output results
734cb32e2e7SValeria Barra {
7359b072555Sjeremylt KSPType ksp_type;
736cb32e2e7SValeria Barra KSPConvergedReason reason;
737cb32e2e7SValeria Barra PetscReal rnorm;
738cb32e2e7SValeria Barra PetscInt its;
7392b730f8bSJeremy L Thompson PetscCall(KSPGetType(ksp, &ksp_type));
7402b730f8bSJeremy L Thompson PetscCall(KSPGetConvergedReason(ksp, &reason));
7412b730f8bSJeremy L Thompson PetscCall(KSPGetIterationNumber(ksp, &its));
7422b730f8bSJeremy L Thompson PetscCall(KSPGetResidualNorm(ksp, &rnorm));
743a22c4fb5SJed Brown if (!test_mode || reason < 0 || rnorm > 1e-8) {
7442b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm,
745cb32e2e7SValeria Barra " KSP:\n"
746cb32e2e7SValeria Barra " KSP Type : %s\n"
747cb32e2e7SValeria Barra " KSP Convergence : %s\n"
74808140895SJed Brown " Total KSP Iterations : %" PetscInt_FMT "\n"
749cb32e2e7SValeria Barra " Final rnorm : %e\n",
7502b730f8bSJeremy L Thompson ksp_type, KSPConvergedReasons[reason], its, (double)rnorm));
751cb32e2e7SValeria Barra }
7528e87e98bSJed Brown if (!test_mode) {
7532b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, " Performance:\n"));
7548e87e98bSJed Brown }
7558e87e98bSJed Brown {
7569b072555Sjeremylt PetscReal max_error;
7572b730f8bSJeremy L Thompson PetscCall(ComputeErrorMax(op_apply_ctx, op_error, X, target, &max_error));
7588e87e98bSJed Brown PetscReal tol = 5e-2;
7599b072555Sjeremylt if (!test_mode || max_error > tol) {
7602b730f8bSJeremy L Thompson PetscCall(MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm));
7612b730f8bSJeremy L Thompson PetscCall(MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm));
7622b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm,
7638e87e98bSJed Brown " Pointwise Error (max) : %e\n"
7648e87e98bSJed Brown " CG Solve Time : %g (%g) sec\n",
7652b730f8bSJeremy L Thompson (double)max_error, rt_max, rt_min));
766cb32e2e7SValeria Barra }
767cb32e2e7SValeria Barra }
7688d0bb2bbSvaleriabarra if (!test_mode) {
7691a8516d0SJames Wright PetscCall(PetscPrintf(comm, " DoFs/Sec in CG : %g (%g) million\n", 1e-6 * gsize * its / rt_max,
7701a8516d0SJames Wright 1e-6 * gsize * its / rt_min));
771cb32e2e7SValeria Barra }
772cb32e2e7SValeria Barra }
773cb32e2e7SValeria Barra
774cb32e2e7SValeria Barra if (write_solution) {
7759b072555Sjeremylt PetscViewer vtk_viewer_soln;
776cb32e2e7SValeria Barra
7772b730f8bSJeremy L Thompson PetscCall(PetscViewerCreate(comm, &vtk_viewer_soln));
7782b730f8bSJeremy L Thompson PetscCall(PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK));
7792b730f8bSJeremy L Thompson PetscCall(PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"));
7802b730f8bSJeremy L Thompson PetscCall(VecView(X, vtk_viewer_soln));
7812b730f8bSJeremy L Thompson PetscCall(PetscViewerDestroy(&vtk_viewer_soln));
782cb32e2e7SValeria Barra }
783cb32e2e7SValeria Barra
7842b730f8bSJeremy L Thompson PetscCall(VecDestroy(&rhs));
7852b730f8bSJeremy L Thompson PetscCall(VecDestroy(&rhs_loc));
7862b730f8bSJeremy L Thompson PetscCall(VecDestroy(&X));
7872b730f8bSJeremy L Thompson PetscCall(VecDestroy(&op_apply_ctx->X_loc));
7882b730f8bSJeremy L Thompson PetscCall(VecDestroy(&op_apply_ctx->Y_loc));
7892b730f8bSJeremy L Thompson PetscCall(VecScatterDestroy(&l_to_g));
7902b730f8bSJeremy L Thompson PetscCall(VecScatterDestroy(&l_to_g_0));
7912b730f8bSJeremy L Thompson PetscCall(VecScatterDestroy(&g_to_g_D));
7922b730f8bSJeremy L Thompson PetscCall(MatDestroy(&mat));
7932b730f8bSJeremy L Thompson PetscCall(KSPDestroy(&ksp));
794cb32e2e7SValeria Barra
795d4d45553Srezgarshakeri CeedVectorDestroy(&op_apply_ctx->x_ceed);
796d4d45553Srezgarshakeri CeedVectorDestroy(&op_apply_ctx->y_ceed);
797d4d45553Srezgarshakeri CeedVectorDestroy(&op_apply_ctx->q_data);
798cb32e2e7SValeria Barra CeedVectorDestroy(&target);
7999b072555Sjeremylt CeedElemRestrictionDestroy(&elem_restr_u);
8009b072555Sjeremylt CeedElemRestrictionDestroy(&elem_restr_x);
8019b072555Sjeremylt CeedElemRestrictionDestroy(&elem_restr_u_i);
8029b072555Sjeremylt CeedElemRestrictionDestroy(&elem_restr_qd_i);
80342a37a0bSJeremy L Thompson CeedBasisDestroy(&basis_u);
80442a37a0bSJeremy L Thompson CeedBasisDestroy(&basis_x);
8059b072555Sjeremylt CeedQFunctionDestroy(&qf_setup_geo);
8069b072555Sjeremylt CeedQFunctionDestroy(&qf_setup_rhs);
8079b072555Sjeremylt CeedQFunctionDestroy(&qf_apply);
8089b072555Sjeremylt CeedQFunctionDestroy(&qf_error);
80942a37a0bSJeremy L Thompson CeedOperatorDestroy(&op_setup_geo);
81042a37a0bSJeremy L Thompson CeedOperatorDestroy(&op_setup_rhs);
81142a37a0bSJeremy L Thompson CeedOperatorDestroy(&op_apply);
81242a37a0bSJeremy L Thompson CeedOperatorDestroy(&op_error);
813cb32e2e7SValeria Barra CeedDestroy(&ceed);
81442a37a0bSJeremy L Thompson
8152b730f8bSJeremy L Thompson PetscCall(PetscFree(op_apply_ctx));
816cb32e2e7SValeria Barra return PetscFinalize();
817cb32e2e7SValeria Barra }
818