xref: /libCEED/examples/petsc/multigrid.c (revision 5aed82e4fa97acf4ba24a7f10a35f5303a6798e0)
1*5aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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.
36c5df90dSjeremylt //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
56c5df90dSjeremylt //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
76c5df90dSjeremylt 
86c5df90dSjeremylt //                        libCEED + PETSc Example: CEED BPs 3-6 with Multigrid
96c5df90dSjeremylt //
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.
116c5df90dSjeremylt //
126c5df90dSjeremylt // The code uses higher level communication protocols in DMPlex.
136c5df90dSjeremylt //
146c5df90dSjeremylt // Build with:
156c5df90dSjeremylt //
166c5df90dSjeremylt //     make multigrid [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
176c5df90dSjeremylt //
186c5df90dSjeremylt // Sample runs:
196c5df90dSjeremylt //
206c5df90dSjeremylt //     multigrid -problem bp3
2128688798Sjeremylt //     multigrid -problem bp4
2228688798Sjeremylt //     multigrid -problem bp5 -ceed /cpu/self
236c5df90dSjeremylt //     multigrid -problem bp6 -ceed /gpu/cuda
246c5df90dSjeremylt //
25fa5adaf5SJeremy L Thompson //TESTARGS(name="BP3, hex elements") -ceed {ceed_resource} -test -problem bp3 -degree 3
26fa5adaf5SJeremy L Thompson //TESTARGS(name="BP3, tet elements") -ceed {ceed_resource} -test -problem bp3 -degree 3 -simplex
276c5df90dSjeremylt 
286c5df90dSjeremylt /// @file
296c5df90dSjeremylt /// CEED BPs 1-6 multigrid example using PETSc
306c5df90dSjeremylt const char help[] = "Solve CEED BPs using p-multigrid with PETSc and DMPlex\n";
316c5df90dSjeremylt 
32636cccdbSjeremylt #include <ceed.h>
33636cccdbSjeremylt #include <petsc.h>
34636cccdbSjeremylt #include <petscdmplex.h>
35636cccdbSjeremylt #include <petscksp.h>
36636cccdbSjeremylt #include <petscsys.h>
372b730f8bSJeremy L Thompson #include <stdbool.h>
382b730f8bSJeremy L Thompson #include <string.h>
39636cccdbSjeremylt 
40e83e87a5Sjeremylt #include "bps.h"
41636cccdbSjeremylt #include "include/bpsproblemdata.h"
422b730f8bSJeremy L Thompson #include "include/libceedsetup.h"
432b730f8bSJeremy L Thompson #include "include/matops.h"
44636cccdbSjeremylt #include "include/petscutils.h"
45b8962995SJeremy L Thompson #include "include/petscversion.h"
46636cccdbSjeremylt #include "include/structs.h"
47636cccdbSjeremylt 
486c5df90dSjeremylt int main(int argc, char **argv) {
496c5df90dSjeremylt   MPI_Comm comm;
502b730f8bSJeremy L Thompson   char     filename[PETSC_MAX_PATH_LEN], ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self";
516c5df90dSjeremylt   double   my_rt_start, my_rt, rt_min, rt_max;
522b730f8bSJeremy L Thompson   PetscInt degree = 3, q_extra, *l_size, *xl_size, *g_size, dim = 3, fine_level, mesh_elem[3] = {3, 3, 3}, num_comp_u = 1, num_levels = degree,
532b730f8bSJeremy L Thompson            *level_degrees;
54cfa59c5bSRey   PetscScalar           eps = 1.0;
55de1229c5Srezgarshakeri   PetscBool             test_mode, benchmark_mode, read_mesh, write_solution, simplex;
569b072555Sjeremylt   PetscLogStage         solve_stage;
5705b9c820SJed Brown   PetscLogEvent         assemble_event;
589b072555Sjeremylt   DM                   *dm, dm_orig;
596c5df90dSjeremylt   KSP                   ksp;
606c5df90dSjeremylt   PC                    pc;
619b072555Sjeremylt   Mat                  *mat_O, *mat_pr, mat_coarse;
629b072555Sjeremylt   Vec                  *X, *X_loc, *mult, rhs, rhs_loc;
639b072555Sjeremylt   PetscMemType          mem_type;
646c88e6a2Srezgarshakeri   OperatorApplyContext *op_apply_ctx, op_error_ctx;
65d4d45553Srezgarshakeri   ProlongRestrContext  *pr_restr_ctx;
666c5df90dSjeremylt   Ceed                  ceed;
679b072555Sjeremylt   CeedData             *ceed_data;
689b072555Sjeremylt   CeedVector            rhs_ceed, target;
691c9a79dbSrezgarshakeri   CeedQFunction         qf_error;
709b072555Sjeremylt   CeedOperator          op_error;
719b072555Sjeremylt   BPType                bp_choice;
729b072555Sjeremylt   CoarsenType           coarsen;
736c5df90dSjeremylt 
742b730f8bSJeremy L Thompson   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
756c5df90dSjeremylt   comm = PETSC_COMM_WORLD;
766c5df90dSjeremylt 
776c5df90dSjeremylt   // Parse command line options
7867490bc6SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
799b072555Sjeremylt   bp_choice = CEED_BP3;
802b730f8bSJeremy L Thompson   PetscCall(PetscOptionsEnum("-problem", "CEED benchmark problem to solve", NULL, bp_types, (PetscEnum)bp_choice, (PetscEnum *)&bp_choice, NULL));
819b072555Sjeremylt   num_comp_u = bp_options[bp_choice].num_comp_u;
826c5df90dSjeremylt   test_mode  = PETSC_FALSE;
832b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, test_mode, &test_mode, NULL));
846c5df90dSjeremylt   benchmark_mode = PETSC_FALSE;
852b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-benchmark", "Benchmarking mode (prints benchmark statistics)", NULL, benchmark_mode, &benchmark_mode, NULL));
866c5df90dSjeremylt   write_solution = PETSC_FALSE;
872b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-write_solution", "Write solution for visualization", NULL, write_solution, &write_solution, NULL));
88de1229c5Srezgarshakeri   simplex = PETSC_FALSE;
892b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-simplex", "Element topology (default:hex)", NULL, simplex, &simplex, NULL));
90eb6a3b92Srezgarshakeri   if ((bp_choice == CEED_BP5 || bp_choice == CEED_BP6) && (simplex)) {
912b730f8bSJeremy L Thompson     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "BP5/6 is not supported with simplex");
92eb6a3b92Srezgarshakeri   }
932b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-eps", "Epsilon parameter for Kershaw mesh transformation", NULL, eps, &eps, NULL));
942b730f8bSJeremy L Thompson   if (eps > 1 || eps <= 0) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "-eps %g must be (0,1]", (double)PetscRealPart(eps));
956c5df90dSjeremylt   degree = test_mode ? 3 : 2;
962b730f8bSJeremy L Thompson   PetscCall(PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", NULL, degree, &degree, NULL));
972b730f8bSJeremy L Thompson   if (degree < 1) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "-degree %" PetscInt_FMT " must be at least 1", degree);
989b072555Sjeremylt   q_extra = bp_options[bp_choice].q_extra;
992b730f8bSJeremy L Thompson   PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, q_extra, &q_extra, NULL));
1002b730f8bSJeremy L Thompson   PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, ceed_resource, ceed_resource, sizeof(ceed_resource), NULL));
1016c5df90dSjeremylt   coarsen = COARSEN_UNIFORM;
1022b730f8bSJeremy L Thompson   PetscCall(PetscOptionsEnum("-coarsen", "Coarsening strategy to use", NULL, coarsen_types, (PetscEnum)coarsen, (PetscEnum *)&coarsen, NULL));
103cb32e2e7SValeria Barra   read_mesh = PETSC_FALSE;
1042b730f8bSJeremy L Thompson   PetscCall(PetscOptionsString("-mesh", "Read mesh from file", NULL, filename, filename, sizeof(filename), &read_mesh));
1056c5df90dSjeremylt   if (!read_mesh) {
1066c5df90dSjeremylt     PetscInt tmp = dim;
1072b730f8bSJeremy L Thompson     PetscCall(PetscOptionsIntArray("-cells", "Number of cells per dimension", NULL, mesh_elem, &tmp, NULL));
1086c5df90dSjeremylt   }
10967490bc6SJeremy L Thompson   PetscOptionsEnd();
1106c5df90dSjeremylt 
1119396343dSjeremylt   // Set up libCEED
1129b072555Sjeremylt   CeedInit(ceed_resource, &ceed);
1139b072555Sjeremylt   CeedMemType mem_type_backend;
1149b072555Sjeremylt   CeedGetPreferredMemType(ceed, &mem_type_backend);
1159396343dSjeremylt 
1166c5df90dSjeremylt   // Setup DM
1176c5df90dSjeremylt   if (read_mesh) {
1182b730f8bSJeremy L Thompson     PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, NULL, PETSC_TRUE, &dm_orig));
1196c5df90dSjeremylt   } else {
1202b730f8bSJeremy L Thompson     PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, simplex, mesh_elem, NULL, NULL, NULL, PETSC_TRUE, &dm_orig));
1216c5df90dSjeremylt   }
1226c5df90dSjeremylt 
1239b072555Sjeremylt   VecType vec_type;
1249b072555Sjeremylt   switch (mem_type_backend) {
1252b730f8bSJeremy L Thompson     case CEED_MEM_HOST:
1262b730f8bSJeremy L Thompson       vec_type = VECSTANDARD;
1272b730f8bSJeremy L Thompson       break;
128b68a8d79SJed Brown     case CEED_MEM_DEVICE: {
129b68a8d79SJed Brown       const char *resolved;
130b68a8d79SJed Brown       CeedGetResource(ceed, &resolved);
1319b072555Sjeremylt       if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA;
1322b730f8bSJeremy L Thompson       else if (strstr(resolved, "/gpu/hip/occa")) vec_type = VECSTANDARD;  // https://github.com/CEED/libCEED/issues/678
1339b072555Sjeremylt       else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP;
1349b072555Sjeremylt       else vec_type = VECSTANDARD;
135b68a8d79SJed Brown     }
136b68a8d79SJed Brown   }
1372b730f8bSJeremy L Thompson   PetscCall(DMSetVecType(dm_orig, vec_type));
1382b730f8bSJeremy L Thompson   PetscCall(DMSetFromOptions(dm_orig));
1392b730f8bSJeremy L Thompson   PetscCall(DMViewFromOptions(dm_orig, NULL, "-dm_view"));
1403fc8a154SJed Brown 
1413fc8a154SJed Brown   // Apply Kershaw mesh transformation
1422b730f8bSJeremy L Thompson   PetscCall(Kershaw(dm_orig, eps));
143b68a8d79SJed Brown 
1446c5df90dSjeremylt   // Allocate arrays for PETSc objects for each level
1456c5df90dSjeremylt   switch (coarsen) {
1466c5df90dSjeremylt     case COARSEN_UNIFORM:
1479b072555Sjeremylt       num_levels = degree;
1486c5df90dSjeremylt       break;
149dc7d240cSValeria Barra     case COARSEN_LOGARITHMIC:
1509b072555Sjeremylt       num_levels = ceil(log(degree) / log(2)) + 1;
1516c5df90dSjeremylt       break;
1526c5df90dSjeremylt   }
1532b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &level_degrees));
1549b072555Sjeremylt   fine_level = num_levels - 1;
15561608365Sjeremylt 
1566c5df90dSjeremylt   switch (coarsen) {
1576c5df90dSjeremylt     case COARSEN_UNIFORM:
1589b072555Sjeremylt       for (int i = 0; i < num_levels; i++) level_degrees[i] = i + 1;
1596c5df90dSjeremylt       break;
160dc7d240cSValeria Barra     case COARSEN_LOGARITHMIC:
1619b072555Sjeremylt       for (int i = 0; i < num_levels - 1; i++) level_degrees[i] = pow(2, i);
1629b072555Sjeremylt       level_degrees[fine_level] = degree;
1636c5df90dSjeremylt       break;
1646c5df90dSjeremylt   }
1652b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &dm));
1662b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &X));
1672b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &X_loc));
1682b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &mult));
1692b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &op_apply_ctx));
1702b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &pr_restr_ctx));
1712b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &mat_O));
1722b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &mat_pr));
1732b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &l_size));
1742b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &xl_size));
1752b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &g_size));
1766c5df90dSjeremylt 
1771c9a79dbSrezgarshakeri   PetscInt c_start, c_end;
1782b730f8bSJeremy L Thompson   PetscCall(DMPlexGetHeightStratum(dm_orig, 0, &c_start, &c_end));
1791c9a79dbSrezgarshakeri   DMPolytopeType cell_type;
1802b730f8bSJeremy L Thompson   PetscCall(DMPlexGetCellType(dm_orig, c_start, &cell_type));
1811c9a79dbSrezgarshakeri   CeedElemTopology elem_topo = ElemTopologyP2C(cell_type);
1821c9a79dbSrezgarshakeri 
1836c5df90dSjeremylt   // Setup DM and Operator Mat Shells for each level
1849b072555Sjeremylt   for (CeedInt i = 0; i < num_levels; i++) {
1856c5df90dSjeremylt     // Create DM
1862b730f8bSJeremy L Thompson     PetscCall(DMClone(dm_orig, &dm[i]));
1872b730f8bSJeremy L Thompson     PetscCall(DMGetVecType(dm_orig, &vec_type));
1882b730f8bSJeremy L Thompson     PetscCall(DMSetVecType(dm[i], vec_type));
189e83e87a5Sjeremylt     PetscInt dim;
1902b730f8bSJeremy L Thompson     PetscCall(DMGetDimension(dm[i], &dim));
1912b730f8bSJeremy L Thompson     PetscCall(SetupDMByDegree(dm[i], level_degrees[fine_level], q_extra, num_comp_u, dim, bp_options[bp_choice].enforce_bc));
1926c5df90dSjeremylt 
1936c5df90dSjeremylt     // Create vectors
1942b730f8bSJeremy L Thompson     PetscCall(DMCreateGlobalVector(dm[i], &X[i]));
1952b730f8bSJeremy L Thompson     PetscCall(VecGetLocalSize(X[i], &l_size[i]));
1962b730f8bSJeremy L Thompson     PetscCall(VecGetSize(X[i], &g_size[i]));
1972b730f8bSJeremy L Thompson     PetscCall(DMCreateLocalVector(dm[i], &X_loc[i]));
1982b730f8bSJeremy L Thompson     PetscCall(VecGetSize(X_loc[i], &xl_size[i]));
1996c5df90dSjeremylt 
2006c5df90dSjeremylt     // Operator
2012b730f8bSJeremy L Thompson     PetscCall(PetscMalloc1(1, &op_apply_ctx[i]));
2022b730f8bSJeremy L Thompson     PetscCall(PetscMalloc1(1, &op_error_ctx));
2032b730f8bSJeremy L Thompson     PetscCall(MatCreateShell(comm, l_size[i], l_size[i], g_size[i], g_size[i], op_apply_ctx[i], &mat_O[i]));
2042b730f8bSJeremy L Thompson     PetscCall(MatShellSetOperation(mat_O[i], MATOP_MULT, (void (*)(void))MatMult_Ceed));
2052b730f8bSJeremy L Thompson     PetscCall(MatShellSetOperation(mat_O[i], MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag));
2062b730f8bSJeremy L Thompson     PetscCall(MatShellSetVecType(mat_O[i], vec_type));
2076c5df90dSjeremylt 
2086c5df90dSjeremylt     // Level transfers
2096c5df90dSjeremylt     if (i > 0) {
2106c5df90dSjeremylt       // Interp
2112b730f8bSJeremy L Thompson       PetscCall(PetscMalloc1(1, &pr_restr_ctx[i]));
2122b730f8bSJeremy L Thompson       PetscCall(MatCreateShell(comm, l_size[i], l_size[i - 1], g_size[i], g_size[i - 1], pr_restr_ctx[i], &mat_pr[i]));
2132b730f8bSJeremy L Thompson       PetscCall(MatShellSetOperation(mat_pr[i], MATOP_MULT, (void (*)(void))MatMult_Prolong));
2142b730f8bSJeremy L Thompson       PetscCall(MatShellSetOperation(mat_pr[i], MATOP_MULT_TRANSPOSE, (void (*)(void))MatMult_Restrict));
2152b730f8bSJeremy L Thompson       PetscCall(MatShellSetVecType(mat_pr[i], vec_type));
2166c5df90dSjeremylt     }
2176c5df90dSjeremylt   }
2182b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(X[fine_level], &rhs));
2196c5df90dSjeremylt 
2206c5df90dSjeremylt   // Print global grid information
2216c5df90dSjeremylt   if (!test_mode) {
2229b072555Sjeremylt     PetscInt P = degree + 1, Q = P + q_extra;
2239396343dSjeremylt 
2249b072555Sjeremylt     const char *used_resource;
2259b072555Sjeremylt     CeedGetResource(ceed, &used_resource);
2269396343dSjeremylt 
2272b730f8bSJeremy L Thompson     PetscCall(VecGetType(X[0], &vec_type));
2289396343dSjeremylt 
2292b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm,
230990fdeb6SJeremy L Thompson                           "\n-- CEED Benchmark Problem %" CeedInt_FMT " -- libCEED + PETSc + PCMG --\n"
2319396343dSjeremylt                           "  PETSc:\n"
2329396343dSjeremylt                           "    PETSc Vec Type                          : %s\n"
2336c5df90dSjeremylt                           "  libCEED:\n"
2346c5df90dSjeremylt                           "    libCEED Backend                         : %s\n"
2359396343dSjeremylt                           "    libCEED Backend MemType                 : %s\n"
2366c5df90dSjeremylt                           "  Mesh:\n"
237751eb813Srezgarshakeri                           "    Solution Order (P)                      : %" CeedInt_FMT "\n"
238751eb813Srezgarshakeri                           "    Quadrature  Order (Q)                   : %" CeedInt_FMT "\n"
239de1229c5Srezgarshakeri                           "    Additional quadrature points (q_extra)  : %" CeedInt_FMT "\n"
24008140895SJed Brown                           "    Global Nodes                            : %" PetscInt_FMT "\n"
24108140895SJed Brown                           "    Owned Nodes                             : %" PetscInt_FMT "\n"
24208140895SJed Brown                           "    DoF per node                            : %" PetscInt_FMT "\n"
24351ad7d5bSrezgarshakeri                           "    Element topology                        : %s\n"
2446c5df90dSjeremylt                           "  Multigrid:\n"
245990fdeb6SJeremy L Thompson                           "    Number of Levels                        : %" CeedInt_FMT "\n",
2462b730f8bSJeremy L Thompson                           bp_choice + 1, vec_type, used_resource, CeedMemTypes[mem_type_backend], P, Q, q_extra, g_size[fine_level] / num_comp_u,
2472b730f8bSJeremy L Thompson                           l_size[fine_level] / num_comp_u, num_comp_u, CeedElemTopologies[elem_topo], num_levels));
2486c5df90dSjeremylt   }
2496c5df90dSjeremylt 
2506c5df90dSjeremylt   // Create RHS vector
2512b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(X_loc[fine_level], &rhs_loc));
2522b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(rhs_loc));
2539b072555Sjeremylt   CeedVectorCreate(ceed, xl_size[fine_level], &rhs_ceed);
254179e5961SZach Atkins   PetscCall(VecP2C(rhs_loc, &mem_type, rhs_ceed));
2556c5df90dSjeremylt 
2566c5df90dSjeremylt   // Set up libCEED operators on each level
2572b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &ceed_data));
258990fdeb6SJeremy L Thompson   for (PetscInt i = 0; i < num_levels; i++) {
2596c5df90dSjeremylt     // Print level information
2609b072555Sjeremylt     if (!test_mode && (i == 0 || i == fine_level)) {
2612b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm,
2622b730f8bSJeremy L Thompson                             "    Level %" PetscInt_FMT " (%s):\n"
263751eb813Srezgarshakeri                             "      Solution Order (P)                    : %" CeedInt_FMT "\n"
26408140895SJed Brown                             "      Global Nodes                          : %" PetscInt_FMT "\n"
26508140895SJed Brown                             "      Owned Nodes                           : %" PetscInt_FMT "\n",
2662b730f8bSJeremy L Thompson                             i, (i ? "fine" : "coarse"), level_degrees[i] + 1, g_size[i] / num_comp_u, l_size[i] / num_comp_u));
2676c5df90dSjeremylt     }
2682b730f8bSJeremy L Thompson     PetscCall(PetscMalloc1(1, &ceed_data[i]));
2692b730f8bSJeremy L Thompson     PetscCall(SetupLibceedByDegree(dm[i], ceed, level_degrees[i], dim, q_extra, dim, num_comp_u, g_size[i], xl_size[i], bp_options[bp_choice],
2702b730f8bSJeremy L Thompson                                    ceed_data[i], i == (fine_level), rhs_ceed, &target));
2716c5df90dSjeremylt   }
2726c5df90dSjeremylt 
2736c5df90dSjeremylt   // Gather RHS
274179e5961SZach Atkins   PetscCall(VecC2P(rhs_ceed, mem_type, rhs_loc));
2752b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(rhs));
2762b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(dm[fine_level], rhs_loc, ADD_VALUES, rhs));
2779b072555Sjeremylt   CeedVectorDestroy(&rhs_ceed);
2786c5df90dSjeremylt 
27943eb8658SJeremy L Thompson   // Create the error QFunction
2802b730f8bSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].error, bp_options[bp_choice].error_loc, &qf_error);
2819b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP);
2829b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE);
2832b730f8bSJeremy L Thompson   CeedQFunctionAddInput(qf_error, "qdata", ceed_data[fine_level]->q_data_size, CEED_EVAL_NONE);
28438f32c05Srezgarshakeri   CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_INTERP);
2856c5df90dSjeremylt 
2866c5df90dSjeremylt   // Create the error operator
2872b730f8bSJeremy L Thompson   CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_error);
2882b730f8bSJeremy L Thompson   CeedOperatorSetField(op_error, "u", ceed_data[fine_level]->elem_restr_u, ceed_data[fine_level]->basis_u, CEED_VECTOR_ACTIVE);
289356036faSJeremy L Thompson   CeedOperatorSetField(op_error, "true_soln", ceed_data[fine_level]->elem_restr_u_i, CEED_BASIS_NONE, target);
290356036faSJeremy L Thompson   CeedOperatorSetField(op_error, "qdata", ceed_data[fine_level]->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data[fine_level]->q_data);
2912b730f8bSJeremy L Thompson   CeedOperatorSetField(op_error, "error", ceed_data[fine_level]->elem_restr_u, ceed_data[fine_level]->basis_u, CEED_VECTOR_ACTIVE);
2926c5df90dSjeremylt 
2936c5df90dSjeremylt   // Calculate multiplicity
2949b072555Sjeremylt   for (int i = 0; i < num_levels; i++) {
295179e5961SZach Atkins     PetscMemType mem_type;
2966c5df90dSjeremylt 
2976c5df90dSjeremylt     // CEED vector
2982b730f8bSJeremy L Thompson     PetscCall(VecZeroEntries(X_loc[i]));
299179e5961SZach Atkins     PetscCall(VecP2C(X_loc[i], &mem_type, ceed_data[i]->x_ceed));
3006c5df90dSjeremylt 
3016c5df90dSjeremylt     // Multiplicity
3022b730f8bSJeremy L Thompson     CeedElemRestrictionGetMultiplicity(ceed_data[i]->elem_restr_u, ceed_data[i]->x_ceed);
3039b072555Sjeremylt     CeedVectorSyncArray(ceed_data[i]->x_ceed, CEED_MEM_HOST);
3046c5df90dSjeremylt 
3056c5df90dSjeremylt     // Restore vector
306179e5961SZach Atkins     PetscCall(VecC2P(ceed_data[i]->x_ceed, mem_type, X_loc[i]));
3076c5df90dSjeremylt 
3086c5df90dSjeremylt     // Creat mult vector
3092b730f8bSJeremy L Thompson     PetscCall(VecDuplicate(X_loc[i], &mult[i]));
3106c5df90dSjeremylt 
3116c5df90dSjeremylt     // Local-to-global
3122b730f8bSJeremy L Thompson     PetscCall(VecZeroEntries(X[i]));
3132b730f8bSJeremy L Thompson     PetscCall(DMLocalToGlobal(dm[i], X_loc[i], ADD_VALUES, X[i]));
3142b730f8bSJeremy L Thompson     PetscCall(VecZeroEntries(X_loc[i]));
3156c5df90dSjeremylt 
3166c5df90dSjeremylt     // Global-to-local
3172b730f8bSJeremy L Thompson     PetscCall(DMGlobalToLocal(dm[i], X[i], INSERT_VALUES, mult[i]));
3182b730f8bSJeremy L Thompson     PetscCall(VecZeroEntries(X[i]));
3196c5df90dSjeremylt 
3206c5df90dSjeremylt     // Multiplicity scaling
3212b730f8bSJeremy L Thompson     PetscCall(VecReciprocal(mult[i]));
3226c5df90dSjeremylt   }
3236c5df90dSjeremylt 
3246c5df90dSjeremylt   // Set up Mat
3259b072555Sjeremylt   for (int i = 0; i < num_levels; i++) {
3266c88e6a2Srezgarshakeri     // Set up apply operator context
3272b730f8bSJeremy L Thompson     PetscCall(SetupApplyOperatorCtx(comm, dm[i], ceed, ceed_data[i], X_loc[i], op_apply_ctx[i]));
3286c5df90dSjeremylt 
3296c5df90dSjeremylt     if (i > 0) {
330a97643b0Sjeremylt       // Prolongation/Restriction Operator
3312b730f8bSJeremy L Thompson       PetscCall(CeedLevelTransferSetup(dm[i - 1], ceed, i, num_comp_u, ceed_data, bp_options[bp_choice], mult[i]));
332d4d45553Srezgarshakeri       pr_restr_ctx[i]->comm        = comm;
333d4d45553Srezgarshakeri       pr_restr_ctx[i]->dmf         = dm[i];
334d4d45553Srezgarshakeri       pr_restr_ctx[i]->dmc         = dm[i - 1];
335d4d45553Srezgarshakeri       pr_restr_ctx[i]->loc_vec_c   = X_loc[i - 1];
336d4d45553Srezgarshakeri       pr_restr_ctx[i]->loc_vec_f   = op_apply_ctx[i]->Y_loc;
337d4d45553Srezgarshakeri       pr_restr_ctx[i]->mult_vec    = mult[i];
338d4d45553Srezgarshakeri       pr_restr_ctx[i]->ceed_vec_c  = op_apply_ctx[i - 1]->x_ceed;
339d4d45553Srezgarshakeri       pr_restr_ctx[i]->ceed_vec_f  = op_apply_ctx[i]->y_ceed;
340d4d45553Srezgarshakeri       pr_restr_ctx[i]->op_prolong  = ceed_data[i]->op_prolong;
341d4d45553Srezgarshakeri       pr_restr_ctx[i]->op_restrict = ceed_data[i]->op_restrict;
342d4d45553Srezgarshakeri       pr_restr_ctx[i]->ceed        = ceed;
3436c5df90dSjeremylt     }
3446c5df90dSjeremylt   }
3456c5df90dSjeremylt 
34653b04fa6SJed Brown   // Assemble coarse grid Jacobian for AMG (or other sparse matrix) solve
3472b730f8bSJeremy L Thompson   PetscCall(DMCreateMatrix(dm[0], &mat_coarse));
34853b04fa6SJed Brown 
3492b730f8bSJeremy L Thompson   PetscCall(PetscLogEventRegister("AssembleMatrix", MAT_CLASSID, &assemble_event));
350cffe6a52SJeremy L Thompson   {
351cffe6a52SJeremy L Thompson     // Assemble matrix analytically
3523047f789SJeremy L Thompson     PetscCount num_entries;
3533047f789SJeremy L Thompson     CeedInt   *rows, *cols;
35453b04fa6SJed Brown     CeedVector coo_values;
3552b730f8bSJeremy L Thompson     CeedOperatorLinearAssembleSymbolic(op_apply_ctx[0]->op, &num_entries, &rows, &cols);
35653b04fa6SJed Brown     ISLocalToGlobalMapping ltog_row, ltog_col;
3572b730f8bSJeremy L Thompson     PetscCall(MatGetLocalToGlobalMapping(mat_coarse, &ltog_row, &ltog_col));
3582b730f8bSJeremy L Thompson     PetscCall(ISLocalToGlobalMappingApply(ltog_row, num_entries, rows, rows));
3592b730f8bSJeremy L Thompson     PetscCall(ISLocalToGlobalMappingApply(ltog_col, num_entries, cols, cols));
3602b730f8bSJeremy L Thompson     PetscCall(MatSetPreallocationCOO(mat_coarse, num_entries, rows, cols));
36153b04fa6SJed Brown     free(rows);
36253b04fa6SJed Brown     free(cols);
36353b04fa6SJed Brown     CeedVectorCreate(ceed, num_entries, &coo_values);
3642b730f8bSJeremy L Thompson     PetscCall(PetscLogEventBegin(assemble_event, mat_coarse, 0, 0, 0));
365d4d45553Srezgarshakeri     CeedOperatorLinearAssemble(op_apply_ctx[0]->op, coo_values);
36653b04fa6SJed Brown     const CeedScalar *values;
36753b04fa6SJed Brown     CeedVectorGetArrayRead(coo_values, CEED_MEM_HOST, &values);
3682b730f8bSJeremy L Thompson     PetscCall(MatSetValuesCOO(mat_coarse, values, ADD_VALUES));
36953b04fa6SJed Brown     CeedVectorRestoreArrayRead(coo_values, &values);
3702b730f8bSJeremy L Thompson     PetscCall(PetscLogEventEnd(assemble_event, mat_coarse, 0, 0, 0));
37153b04fa6SJed Brown     CeedVectorDestroy(&coo_values);
37253b04fa6SJed Brown   }
37315ce0ef0Sjeremylt 
3746c5df90dSjeremylt   // Set up KSP
3752b730f8bSJeremy L Thompson   PetscCall(KSPCreate(comm, &ksp));
3766c5df90dSjeremylt   {
3772b730f8bSJeremy L Thompson     PetscCall(KSPSetType(ksp, KSPCG));
3782b730f8bSJeremy L Thompson     PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
3792b730f8bSJeremy L Thompson     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
3806c5df90dSjeremylt   }
3812b730f8bSJeremy L Thompson   PetscCall(KSPSetFromOptions(ksp));
3822b730f8bSJeremy L Thompson   PetscCall(KSPSetOperators(ksp, mat_O[fine_level], mat_O[fine_level]));
3836c5df90dSjeremylt 
3846c5df90dSjeremylt   // Set up PCMG
3852b730f8bSJeremy L Thompson   PetscCall(KSPGetPC(ksp, &pc));
3869b072555Sjeremylt   PCMGCycleType pcmg_cycle_type = PC_MG_CYCLE_V;
3876c5df90dSjeremylt   {
3882b730f8bSJeremy L Thompson     PetscCall(PCSetType(pc, PCMG));
3896c5df90dSjeremylt 
3906c5df90dSjeremylt     // PCMG levels
3912b730f8bSJeremy L Thompson     PetscCall(PCMGSetLevels(pc, num_levels, NULL));
3929b072555Sjeremylt     for (int i = 0; i < num_levels; i++) {
3936c5df90dSjeremylt       // Smoother
3946c5df90dSjeremylt       KSP smoother;
3956c5df90dSjeremylt       PC  smoother_pc;
3962b730f8bSJeremy L Thompson       PetscCall(PCMGGetSmoother(pc, i, &smoother));
3972b730f8bSJeremy L Thompson       PetscCall(KSPSetType(smoother, KSPCHEBYSHEV));
3982b730f8bSJeremy L Thompson       PetscCall(KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1));
3992b730f8bSJeremy L Thompson       PetscCall(KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE));
4002b730f8bSJeremy L Thompson       PetscCall(KSPSetOperators(smoother, mat_O[i], mat_O[i]));
4012b730f8bSJeremy L Thompson       PetscCall(KSPGetPC(smoother, &smoother_pc));
4022b730f8bSJeremy L Thompson       PetscCall(PCSetType(smoother_pc, PCJACOBI));
4032b730f8bSJeremy L Thompson       PetscCall(PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL));
4046c5df90dSjeremylt 
4056c5df90dSjeremylt       // Work vector
4069b072555Sjeremylt       if (i < num_levels - 1) {
4072b730f8bSJeremy L Thompson         PetscCall(PCMGSetX(pc, i, X[i]));
4086c5df90dSjeremylt       }
4096c5df90dSjeremylt 
4106c5df90dSjeremylt       // Level transfers
4116c5df90dSjeremylt       if (i > 0) {
4126c5df90dSjeremylt         // Interpolation
4132b730f8bSJeremy L Thompson         PetscCall(PCMGSetInterpolation(pc, i, mat_pr[i]));
4146c5df90dSjeremylt       }
4156c5df90dSjeremylt 
4166c5df90dSjeremylt       // Coarse solve
4176c5df90dSjeremylt       KSP coarse;
4186c5df90dSjeremylt       PC  coarse_pc;
4192b730f8bSJeremy L Thompson       PetscCall(PCMGGetCoarseSolve(pc, &coarse));
4202b730f8bSJeremy L Thompson       PetscCall(KSPSetType(coarse, KSPPREONLY));
4212b730f8bSJeremy L Thompson       PetscCall(KSPSetOperators(coarse, mat_coarse, mat_coarse));
42215ce0ef0Sjeremylt 
4232b730f8bSJeremy L Thompson       PetscCall(KSPGetPC(coarse, &coarse_pc));
4242b730f8bSJeremy L Thompson       PetscCall(PCSetType(coarse_pc, PCGAMG));
42515ce0ef0Sjeremylt 
4262b730f8bSJeremy L Thompson       PetscCall(KSPSetOptionsPrefix(coarse, "coarse_"));
4272b730f8bSJeremy L Thompson       PetscCall(PCSetOptionsPrefix(coarse_pc, "coarse_"));
4282b730f8bSJeremy L Thompson       PetscCall(KSPSetFromOptions(coarse));
4292b730f8bSJeremy L Thompson       PetscCall(PCSetFromOptions(coarse_pc));
4306c5df90dSjeremylt     }
4316c5df90dSjeremylt 
4326c5df90dSjeremylt     // PCMG options
4332b730f8bSJeremy L Thompson     PetscCall(PCMGSetType(pc, PC_MG_MULTIPLICATIVE));
4342b730f8bSJeremy L Thompson     PetscCall(PCMGSetNumberSmooth(pc, 3));
4352b730f8bSJeremy L Thompson     PetscCall(PCMGSetCycleType(pc, pcmg_cycle_type));
4366c5df90dSjeremylt   }
4376c5df90dSjeremylt 
4386c5df90dSjeremylt   // First run, if benchmarking
4396c5df90dSjeremylt   if (benchmark_mode) {
4402b730f8bSJeremy L Thompson     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1));
4412b730f8bSJeremy L Thompson     PetscCall(VecZeroEntries(X[fine_level]));
4426c5df90dSjeremylt     my_rt_start = MPI_Wtime();
4432b730f8bSJeremy L Thompson     PetscCall(KSPSolve(ksp, rhs, X[fine_level]));
4446c5df90dSjeremylt     my_rt = MPI_Wtime() - my_rt_start;
4452b730f8bSJeremy L Thompson     PetscCall(MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm));
4466c5df90dSjeremylt     // Set maxits based on first iteration timing
4476c5df90dSjeremylt     if (my_rt > 0.02) {
4482b730f8bSJeremy L Thompson       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5));
4496c5df90dSjeremylt     } else {
4502b730f8bSJeremy L Thompson       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20));
4516c5df90dSjeremylt     }
4526c5df90dSjeremylt   }
4536c5df90dSjeremylt 
4546c5df90dSjeremylt   // Timed solve
4552b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(X[fine_level]));
4562b730f8bSJeremy L Thompson   PetscCall(PetscBarrier((PetscObject)ksp));
45709a940d7Sjeremylt 
45809a940d7Sjeremylt   // -- Performance logging
4592b730f8bSJeremy L Thompson   PetscCall(PetscLogStageRegister("Solve Stage", &solve_stage));
4602b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePush(solve_stage));
46109a940d7Sjeremylt 
46209a940d7Sjeremylt   // -- Solve
4636c5df90dSjeremylt   my_rt_start = MPI_Wtime();
4642b730f8bSJeremy L Thompson   PetscCall(KSPSolve(ksp, rhs, X[fine_level]));
4656c5df90dSjeremylt   my_rt = MPI_Wtime() - my_rt_start;
4666c5df90dSjeremylt 
46709a940d7Sjeremylt   // -- Performance logging
4682b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePop());
46909a940d7Sjeremylt 
4706c5df90dSjeremylt   // Output results
4716c5df90dSjeremylt   {
4729b072555Sjeremylt     KSPType            ksp_type;
4739b072555Sjeremylt     PCMGType           pcmg_type;
4746c5df90dSjeremylt     KSPConvergedReason reason;
4756c5df90dSjeremylt     PetscReal          rnorm;
4766c5df90dSjeremylt     PetscInt           its;
4772b730f8bSJeremy L Thompson     PetscCall(KSPGetType(ksp, &ksp_type));
4782b730f8bSJeremy L Thompson     PetscCall(KSPGetConvergedReason(ksp, &reason));
4792b730f8bSJeremy L Thompson     PetscCall(KSPGetIterationNumber(ksp, &its));
4802b730f8bSJeremy L Thompson     PetscCall(KSPGetResidualNorm(ksp, &rnorm));
4812b730f8bSJeremy L Thompson     PetscCall(PCMGGetType(pc, &pcmg_type));
4826c5df90dSjeremylt     if (!test_mode || reason < 0 || rnorm > 1e-8) {
4832b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm,
4846c5df90dSjeremylt                             "  KSP:\n"
4856c5df90dSjeremylt                             "    KSP Type                                : %s\n"
4866c5df90dSjeremylt                             "    KSP Convergence                         : %s\n"
487a9b2c5ddSrezgarshakeri                             "    Total KSP Iterations                    : %" PetscInt_FMT "\n"
4886c5df90dSjeremylt                             "    Final rnorm                             : %e\n",
4892b730f8bSJeremy L Thompson                             ksp_type, KSPConvergedReasons[reason], its, (double)rnorm));
4902b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm,
4916c5df90dSjeremylt                             "  PCMG:\n"
4926c5df90dSjeremylt                             "    PCMG Type                               : %s\n"
4936c5df90dSjeremylt                             "    PCMG Cycle Type                         : %s\n",
4942b730f8bSJeremy L Thompson                             PCMGTypes[pcmg_type], PCMGCycleTypes[pcmg_cycle_type]));
4956c5df90dSjeremylt     }
4966c5df90dSjeremylt     if (!test_mode) {
4972b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm, "  Performance:\n"));
4986c5df90dSjeremylt     }
4996c5df90dSjeremylt     {
5006c88e6a2Srezgarshakeri       // Set up error operator context
5012b730f8bSJeremy L Thompson       PetscCall(SetupErrorOperatorCtx(comm, dm[fine_level], ceed, ceed_data[fine_level], X_loc[fine_level], op_error, op_error_ctx));
50238f32c05Srezgarshakeri       PetscScalar l2_error;
5032b730f8bSJeremy L Thompson       PetscCall(ComputeL2Error(X[fine_level], &l2_error, op_error_ctx));
5040a8fc04aSrezgarshakeri       PetscReal tol = 5e-2;
50538f32c05Srezgarshakeri       if (!test_mode || l2_error > tol) {
5062b730f8bSJeremy L Thompson         PetscCall(MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm));
5072b730f8bSJeremy L Thompson         PetscCall(MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm));
5082b730f8bSJeremy L Thompson         PetscCall(PetscPrintf(comm,
50938f32c05Srezgarshakeri                               "    L2 Error                                : %e\n"
5106c5df90dSjeremylt                               "    CG Solve Time                           : %g (%g) sec\n",
5112b730f8bSJeremy L Thompson                               (double)l2_error, rt_max, rt_min));
5126c5df90dSjeremylt       }
5136c5df90dSjeremylt     }
5146c5df90dSjeremylt     if (benchmark_mode && (!test_mode)) {
5152b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm, "    DoFs/Sec in CG                            : %g (%g) million\n", 1e-6 * g_size[fine_level] * its / rt_max,
5162b730f8bSJeremy L Thompson                             1e-6 * g_size[fine_level] * its / rt_min));
5176c5df90dSjeremylt     }
5186c5df90dSjeremylt   }
5196c5df90dSjeremylt 
5206c5df90dSjeremylt   if (write_solution) {
5219b072555Sjeremylt     PetscViewer vtk_viewer_soln;
5226c5df90dSjeremylt 
5232b730f8bSJeremy L Thompson     PetscCall(PetscViewerCreate(comm, &vtk_viewer_soln));
5242b730f8bSJeremy L Thompson     PetscCall(PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK));
5252b730f8bSJeremy L Thompson     PetscCall(PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"));
5262b730f8bSJeremy L Thompson     PetscCall(VecView(X[fine_level], vtk_viewer_soln));
5272b730f8bSJeremy L Thompson     PetscCall(PetscViewerDestroy(&vtk_viewer_soln));
5286c5df90dSjeremylt   }
5296c5df90dSjeremylt 
5306c5df90dSjeremylt   // Cleanup
5319b072555Sjeremylt   for (int i = 0; i < num_levels; i++) {
5322b730f8bSJeremy L Thompson     PetscCall(VecDestroy(&X[i]));
5332b730f8bSJeremy L Thompson     PetscCall(VecDestroy(&X_loc[i]));
5342b730f8bSJeremy L Thompson     PetscCall(VecDestroy(&mult[i]));
5352b730f8bSJeremy L Thompson     PetscCall(VecDestroy(&op_apply_ctx[i]->Y_loc));
5362b730f8bSJeremy L Thompson     PetscCall(MatDestroy(&mat_O[i]));
5372b730f8bSJeremy L Thompson     PetscCall(PetscFree(op_apply_ctx[i]));
5386c5df90dSjeremylt     if (i > 0) {
5392b730f8bSJeremy L Thompson       PetscCall(MatDestroy(&mat_pr[i]));
5402b730f8bSJeremy L Thompson       PetscCall(PetscFree(pr_restr_ctx[i]));
5416c5df90dSjeremylt     }
5422b730f8bSJeremy L Thompson     PetscCall(CeedDataDestroy(i, ceed_data[i]));
5432b730f8bSJeremy L Thompson     PetscCall(DMDestroy(&dm[i]));
5446c5df90dSjeremylt   }
5452b730f8bSJeremy L Thompson   PetscCall(PetscFree(level_degrees));
5462b730f8bSJeremy L Thompson   PetscCall(PetscFree(dm));
5472b730f8bSJeremy L Thompson   PetscCall(PetscFree(X));
5482b730f8bSJeremy L Thompson   PetscCall(PetscFree(X_loc));
5492b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&op_error_ctx->Y_loc));
5502b730f8bSJeremy L Thompson   PetscCall(PetscFree(mult));
5512b730f8bSJeremy L Thompson   PetscCall(PetscFree(mat_O));
5522b730f8bSJeremy L Thompson   PetscCall(PetscFree(mat_pr));
5532b730f8bSJeremy L Thompson   PetscCall(PetscFree(ceed_data));
5542b730f8bSJeremy L Thompson   PetscCall(PetscFree(op_apply_ctx));
5552b730f8bSJeremy L Thompson   PetscCall(PetscFree(op_error_ctx));
5562b730f8bSJeremy L Thompson   PetscCall(PetscFree(pr_restr_ctx));
5572b730f8bSJeremy L Thompson   PetscCall(PetscFree(l_size));
5582b730f8bSJeremy L Thompson   PetscCall(PetscFree(xl_size));
5592b730f8bSJeremy L Thompson   PetscCall(PetscFree(g_size));
5602b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&rhs));
5612b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&rhs_loc));
5622b730f8bSJeremy L Thompson   PetscCall(MatDestroy(&mat_coarse));
5632b730f8bSJeremy L Thompson   PetscCall(KSPDestroy(&ksp));
5642b730f8bSJeremy L Thompson   PetscCall(DMDestroy(&dm_orig));
5656c5df90dSjeremylt   CeedVectorDestroy(&target);
5669b072555Sjeremylt   CeedQFunctionDestroy(&qf_error);
5679b072555Sjeremylt   CeedOperatorDestroy(&op_error);
5686c5df90dSjeremylt   CeedDestroy(&ceed);
5696c5df90dSjeremylt   return PetscFinalize();
5706c5df90dSjeremylt }
571