xref: /libCEED/examples/petsc/multigrid.c (revision ea61e9ac44808524e4667c1525a05976f536c19c)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
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 //
10*ea61e9acSJeremy L Thompson // This example demonstrates a simple usage of libCEED with PETSc to solve the CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps.
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 //
256c5df90dSjeremylt //TESTARGS -ceed {ceed_resource} -test -problem bp3 -degree 3
26a9b2c5ddSrezgarshakeri //TESTARGS -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 
48636cccdbSjeremylt #if PETSC_VERSION_LT(3, 12, 0)
49636cccdbSjeremylt #ifdef PETSC_HAVE_CUDA
50636cccdbSjeremylt #include <petsccuda.h>
51*ea61e9acSJeremy L Thompson // Note: With PETSc prior to version 3.12.0, providing the source path to include 'cublas_v2.h' will be needed to use 'petsccuda.h'.
52636cccdbSjeremylt #endif
53636cccdbSjeremylt #endif
546c5df90dSjeremylt 
556c5df90dSjeremylt int main(int argc, char **argv) {
566c5df90dSjeremylt   MPI_Comm comm;
572b730f8bSJeremy L Thompson   char     filename[PETSC_MAX_PATH_LEN], ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self";
586c5df90dSjeremylt   double   my_rt_start, my_rt, rt_min, rt_max;
592b730f8bSJeremy 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,
602b730f8bSJeremy L Thompson            *level_degrees;
616c5df90dSjeremylt   PetscScalar          *r;
62cfa59c5bSRey   PetscScalar           eps = 1.0;
63de1229c5Srezgarshakeri   PetscBool             test_mode, benchmark_mode, read_mesh, write_solution, simplex;
649b072555Sjeremylt   PetscLogStage         solve_stage;
6505b9c820SJed Brown   PetscLogEvent         assemble_event;
669b072555Sjeremylt   DM                   *dm, dm_orig;
676c5df90dSjeremylt   KSP                   ksp;
686c5df90dSjeremylt   PC                    pc;
699b072555Sjeremylt   Mat                  *mat_O, *mat_pr, mat_coarse;
709b072555Sjeremylt   Vec                  *X, *X_loc, *mult, rhs, rhs_loc;
719b072555Sjeremylt   PetscMemType          mem_type;
726c88e6a2Srezgarshakeri   OperatorApplyContext *op_apply_ctx, op_error_ctx;
73d4d45553Srezgarshakeri   ProlongRestrContext  *pr_restr_ctx;
746c5df90dSjeremylt   Ceed                  ceed;
759b072555Sjeremylt   CeedData             *ceed_data;
769b072555Sjeremylt   CeedVector            rhs_ceed, target;
771c9a79dbSrezgarshakeri   CeedQFunction         qf_error;
789b072555Sjeremylt   CeedOperator          op_error;
799b072555Sjeremylt   BPType                bp_choice;
809b072555Sjeremylt   CoarsenType           coarsen;
816c5df90dSjeremylt 
822b730f8bSJeremy L Thompson   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
836c5df90dSjeremylt   comm = PETSC_COMM_WORLD;
846c5df90dSjeremylt 
856c5df90dSjeremylt   // Parse command line options
8667490bc6SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
879b072555Sjeremylt   bp_choice = CEED_BP3;
882b730f8bSJeremy L Thompson   PetscCall(PetscOptionsEnum("-problem", "CEED benchmark problem to solve", NULL, bp_types, (PetscEnum)bp_choice, (PetscEnum *)&bp_choice, NULL));
899b072555Sjeremylt   num_comp_u = bp_options[bp_choice].num_comp_u;
906c5df90dSjeremylt   test_mode  = PETSC_FALSE;
912b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, test_mode, &test_mode, NULL));
926c5df90dSjeremylt   benchmark_mode = PETSC_FALSE;
932b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-benchmark", "Benchmarking mode (prints benchmark statistics)", NULL, benchmark_mode, &benchmark_mode, NULL));
946c5df90dSjeremylt   write_solution = PETSC_FALSE;
952b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-write_solution", "Write solution for visualization", NULL, write_solution, &write_solution, NULL));
96de1229c5Srezgarshakeri   simplex = PETSC_FALSE;
972b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-simplex", "Element topology (default:hex)", NULL, simplex, &simplex, NULL));
98eb6a3b92Srezgarshakeri   if ((bp_choice == CEED_BP5 || bp_choice == CEED_BP6) && (simplex)) {
992b730f8bSJeremy L Thompson     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "BP5/6 is not supported with simplex");
100eb6a3b92Srezgarshakeri   }
1012b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-eps", "Epsilon parameter for Kershaw mesh transformation", NULL, eps, &eps, NULL));
1022b730f8bSJeremy L Thompson   if (eps > 1 || eps <= 0) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "-eps %g must be (0,1]", (double)PetscRealPart(eps));
1036c5df90dSjeremylt   degree = test_mode ? 3 : 2;
1042b730f8bSJeremy L Thompson   PetscCall(PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", NULL, degree, &degree, NULL));
1052b730f8bSJeremy L Thompson   if (degree < 1) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "-degree %" PetscInt_FMT " must be at least 1", degree);
1069b072555Sjeremylt   q_extra = bp_options[bp_choice].q_extra;
1072b730f8bSJeremy L Thompson   PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, q_extra, &q_extra, NULL));
1082b730f8bSJeremy L Thompson   PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, ceed_resource, ceed_resource, sizeof(ceed_resource), NULL));
1096c5df90dSjeremylt   coarsen = COARSEN_UNIFORM;
1102b730f8bSJeremy L Thompson   PetscCall(PetscOptionsEnum("-coarsen", "Coarsening strategy to use", NULL, coarsen_types, (PetscEnum)coarsen, (PetscEnum *)&coarsen, NULL));
111cb32e2e7SValeria Barra   read_mesh = PETSC_FALSE;
1122b730f8bSJeremy L Thompson   PetscCall(PetscOptionsString("-mesh", "Read mesh from file", NULL, filename, filename, sizeof(filename), &read_mesh));
1136c5df90dSjeremylt   if (!read_mesh) {
1146c5df90dSjeremylt     PetscInt tmp = dim;
1152b730f8bSJeremy L Thompson     PetscCall(PetscOptionsIntArray("-cells", "Number of cells per dimension", NULL, mesh_elem, &tmp, NULL));
1166c5df90dSjeremylt   }
11767490bc6SJeremy L Thompson   PetscOptionsEnd();
1186c5df90dSjeremylt 
1199396343dSjeremylt   // Set up libCEED
1209b072555Sjeremylt   CeedInit(ceed_resource, &ceed);
1219b072555Sjeremylt   CeedMemType mem_type_backend;
1229b072555Sjeremylt   CeedGetPreferredMemType(ceed, &mem_type_backend);
1239396343dSjeremylt 
1246c5df90dSjeremylt   // Setup DM
1256c5df90dSjeremylt   if (read_mesh) {
1262b730f8bSJeremy L Thompson     PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, NULL, PETSC_TRUE, &dm_orig));
1276c5df90dSjeremylt   } else {
1282b730f8bSJeremy L Thompson     PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, simplex, mesh_elem, NULL, NULL, NULL, PETSC_TRUE, &dm_orig));
1296c5df90dSjeremylt   }
1306c5df90dSjeremylt 
1319b072555Sjeremylt   VecType vec_type;
1329b072555Sjeremylt   switch (mem_type_backend) {
1332b730f8bSJeremy L Thompson     case CEED_MEM_HOST:
1342b730f8bSJeremy L Thompson       vec_type = VECSTANDARD;
1352b730f8bSJeremy L Thompson       break;
136b68a8d79SJed Brown     case CEED_MEM_DEVICE: {
137b68a8d79SJed Brown       const char *resolved;
138b68a8d79SJed Brown       CeedGetResource(ceed, &resolved);
1399b072555Sjeremylt       if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA;
1402b730f8bSJeremy L Thompson       else if (strstr(resolved, "/gpu/hip/occa")) vec_type = VECSTANDARD;  // https://github.com/CEED/libCEED/issues/678
1419b072555Sjeremylt       else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP;
1429b072555Sjeremylt       else vec_type = VECSTANDARD;
143b68a8d79SJed Brown     }
144b68a8d79SJed Brown   }
1452b730f8bSJeremy L Thompson   PetscCall(DMSetVecType(dm_orig, vec_type));
1462b730f8bSJeremy L Thompson   PetscCall(DMSetFromOptions(dm_orig));
1472b730f8bSJeremy L Thompson   PetscCall(DMViewFromOptions(dm_orig, NULL, "-dm_view"));
1483fc8a154SJed Brown 
1493fc8a154SJed Brown   // Apply Kershaw mesh transformation
1502b730f8bSJeremy L Thompson   PetscCall(Kershaw(dm_orig, eps));
151b68a8d79SJed Brown 
1526c5df90dSjeremylt   // Allocate arrays for PETSc objects for each level
1536c5df90dSjeremylt   switch (coarsen) {
1546c5df90dSjeremylt     case COARSEN_UNIFORM:
1559b072555Sjeremylt       num_levels = degree;
1566c5df90dSjeremylt       break;
157dc7d240cSValeria Barra     case COARSEN_LOGARITHMIC:
1589b072555Sjeremylt       num_levels = ceil(log(degree) / log(2)) + 1;
1596c5df90dSjeremylt       break;
1606c5df90dSjeremylt   }
1612b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &level_degrees));
1629b072555Sjeremylt   fine_level = num_levels - 1;
16361608365Sjeremylt 
1646c5df90dSjeremylt   switch (coarsen) {
1656c5df90dSjeremylt     case COARSEN_UNIFORM:
1669b072555Sjeremylt       for (int i = 0; i < num_levels; i++) level_degrees[i] = i + 1;
1676c5df90dSjeremylt       break;
168dc7d240cSValeria Barra     case COARSEN_LOGARITHMIC:
1699b072555Sjeremylt       for (int i = 0; i < num_levels - 1; i++) level_degrees[i] = pow(2, i);
1709b072555Sjeremylt       level_degrees[fine_level] = degree;
1716c5df90dSjeremylt       break;
1726c5df90dSjeremylt   }
1732b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &dm));
1742b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &X));
1752b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &X_loc));
1762b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &mult));
1772b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &op_apply_ctx));
1782b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &pr_restr_ctx));
1792b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &mat_O));
1802b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &mat_pr));
1812b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &l_size));
1822b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &xl_size));
1832b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &g_size));
1846c5df90dSjeremylt 
1851c9a79dbSrezgarshakeri   PetscInt c_start, c_end;
1862b730f8bSJeremy L Thompson   PetscCall(DMPlexGetHeightStratum(dm_orig, 0, &c_start, &c_end));
1871c9a79dbSrezgarshakeri   DMPolytopeType cell_type;
1882b730f8bSJeremy L Thompson   PetscCall(DMPlexGetCellType(dm_orig, c_start, &cell_type));
1891c9a79dbSrezgarshakeri   CeedElemTopology elem_topo = ElemTopologyP2C(cell_type);
1901c9a79dbSrezgarshakeri 
1916c5df90dSjeremylt   // Setup DM and Operator Mat Shells for each level
1929b072555Sjeremylt   for (CeedInt i = 0; i < num_levels; i++) {
1936c5df90dSjeremylt     // Create DM
1942b730f8bSJeremy L Thompson     PetscCall(DMClone(dm_orig, &dm[i]));
1952b730f8bSJeremy L Thompson     PetscCall(DMGetVecType(dm_orig, &vec_type));
1962b730f8bSJeremy L Thompson     PetscCall(DMSetVecType(dm[i], vec_type));
197e83e87a5Sjeremylt     PetscInt dim;
1982b730f8bSJeremy L Thompson     PetscCall(DMGetDimension(dm[i], &dim));
1992b730f8bSJeremy L Thompson     PetscCall(SetupDMByDegree(dm[i], level_degrees[fine_level], q_extra, num_comp_u, dim, bp_options[bp_choice].enforce_bc));
2006c5df90dSjeremylt 
2016c5df90dSjeremylt     // Create vectors
2022b730f8bSJeremy L Thompson     PetscCall(DMCreateGlobalVector(dm[i], &X[i]));
2032b730f8bSJeremy L Thompson     PetscCall(VecGetLocalSize(X[i], &l_size[i]));
2042b730f8bSJeremy L Thompson     PetscCall(VecGetSize(X[i], &g_size[i]));
2052b730f8bSJeremy L Thompson     PetscCall(DMCreateLocalVector(dm[i], &X_loc[i]));
2062b730f8bSJeremy L Thompson     PetscCall(VecGetSize(X_loc[i], &xl_size[i]));
2076c5df90dSjeremylt 
2086c5df90dSjeremylt     // Operator
2092b730f8bSJeremy L Thompson     PetscCall(PetscMalloc1(1, &op_apply_ctx[i]));
2102b730f8bSJeremy L Thompson     PetscCall(PetscMalloc1(1, &op_error_ctx));
2112b730f8bSJeremy L Thompson     PetscCall(MatCreateShell(comm, l_size[i], l_size[i], g_size[i], g_size[i], op_apply_ctx[i], &mat_O[i]));
2122b730f8bSJeremy L Thompson     PetscCall(MatShellSetOperation(mat_O[i], MATOP_MULT, (void (*)(void))MatMult_Ceed));
2132b730f8bSJeremy L Thompson     PetscCall(MatShellSetOperation(mat_O[i], MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag));
2142b730f8bSJeremy L Thompson     PetscCall(MatShellSetVecType(mat_O[i], vec_type));
2156c5df90dSjeremylt 
2166c5df90dSjeremylt     // Level transfers
2176c5df90dSjeremylt     if (i > 0) {
2186c5df90dSjeremylt       // Interp
2192b730f8bSJeremy L Thompson       PetscCall(PetscMalloc1(1, &pr_restr_ctx[i]));
2202b730f8bSJeremy 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]));
2212b730f8bSJeremy L Thompson       PetscCall(MatShellSetOperation(mat_pr[i], MATOP_MULT, (void (*)(void))MatMult_Prolong));
2222b730f8bSJeremy L Thompson       PetscCall(MatShellSetOperation(mat_pr[i], MATOP_MULT_TRANSPOSE, (void (*)(void))MatMult_Restrict));
2232b730f8bSJeremy L Thompson       PetscCall(MatShellSetVecType(mat_pr[i], vec_type));
2246c5df90dSjeremylt     }
2256c5df90dSjeremylt   }
2262b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(X[fine_level], &rhs));
2276c5df90dSjeremylt 
2286c5df90dSjeremylt   // Print global grid information
2296c5df90dSjeremylt   if (!test_mode) {
2309b072555Sjeremylt     PetscInt P = degree + 1, Q = P + q_extra;
2319396343dSjeremylt 
2329b072555Sjeremylt     const char *used_resource;
2339b072555Sjeremylt     CeedGetResource(ceed, &used_resource);
2349396343dSjeremylt 
2352b730f8bSJeremy L Thompson     PetscCall(VecGetType(X[0], &vec_type));
2369396343dSjeremylt 
2372b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm,
238990fdeb6SJeremy L Thompson                           "\n-- CEED Benchmark Problem %" CeedInt_FMT " -- libCEED + PETSc + PCMG --\n"
2399396343dSjeremylt                           "  PETSc:\n"
2409396343dSjeremylt                           "    PETSc Vec Type                          : %s\n"
2416c5df90dSjeremylt                           "  libCEED:\n"
2426c5df90dSjeremylt                           "    libCEED Backend                         : %s\n"
2439396343dSjeremylt                           "    libCEED Backend MemType                 : %s\n"
2446c5df90dSjeremylt                           "  Mesh:\n"
245751eb813Srezgarshakeri                           "    Solution Order (P)                      : %" CeedInt_FMT "\n"
246751eb813Srezgarshakeri                           "    Quadrature  Order (Q)                   : %" CeedInt_FMT "\n"
247de1229c5Srezgarshakeri                           "    Additional quadrature points (q_extra)  : %" CeedInt_FMT "\n"
24808140895SJed Brown                           "    Global Nodes                            : %" PetscInt_FMT "\n"
24908140895SJed Brown                           "    Owned Nodes                             : %" PetscInt_FMT "\n"
25008140895SJed Brown                           "    DoF per node                            : %" PetscInt_FMT "\n"
25151ad7d5bSrezgarshakeri                           "    Element topology                        : %s\n"
2526c5df90dSjeremylt                           "  Multigrid:\n"
253990fdeb6SJeremy L Thompson                           "    Number of Levels                        : %" CeedInt_FMT "\n",
2542b730f8bSJeremy L Thompson                           bp_choice + 1, vec_type, used_resource, CeedMemTypes[mem_type_backend], P, Q, q_extra, g_size[fine_level] / num_comp_u,
2552b730f8bSJeremy L Thompson                           l_size[fine_level] / num_comp_u, num_comp_u, CeedElemTopologies[elem_topo], num_levels));
2566c5df90dSjeremylt   }
2576c5df90dSjeremylt 
2586c5df90dSjeremylt   // Create RHS vector
2592b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(X_loc[fine_level], &rhs_loc));
2602b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(rhs_loc));
2612b730f8bSJeremy L Thompson   PetscCall(VecGetArrayAndMemType(rhs_loc, &r, &mem_type));
2629b072555Sjeremylt   CeedVectorCreate(ceed, xl_size[fine_level], &rhs_ceed);
2639b072555Sjeremylt   CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r);
2646c5df90dSjeremylt 
2656c5df90dSjeremylt   // Set up libCEED operators on each level
2662b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &ceed_data));
267990fdeb6SJeremy L Thompson   for (PetscInt i = 0; i < num_levels; i++) {
2686c5df90dSjeremylt     // Print level information
2699b072555Sjeremylt     if (!test_mode && (i == 0 || i == fine_level)) {
2702b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm,
2712b730f8bSJeremy L Thompson                             "    Level %" PetscInt_FMT " (%s):\n"
272751eb813Srezgarshakeri                             "      Solution Order (P)                    : %" CeedInt_FMT "\n"
27308140895SJed Brown                             "      Global Nodes                          : %" PetscInt_FMT "\n"
27408140895SJed Brown                             "      Owned Nodes                           : %" PetscInt_FMT "\n",
2752b730f8bSJeremy L Thompson                             i, (i ? "fine" : "coarse"), level_degrees[i] + 1, g_size[i] / num_comp_u, l_size[i] / num_comp_u));
2766c5df90dSjeremylt     }
2772b730f8bSJeremy L Thompson     PetscCall(PetscMalloc1(1, &ceed_data[i]));
2782b730f8bSJeremy 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],
2792b730f8bSJeremy L Thompson                                    ceed_data[i], i == (fine_level), rhs_ceed, &target));
2806c5df90dSjeremylt   }
2816c5df90dSjeremylt 
2826c5df90dSjeremylt   // Gather RHS
2839b072555Sjeremylt   CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL);
2842b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(rhs_loc, &r));
2852b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(rhs));
2862b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(dm[fine_level], rhs_loc, ADD_VALUES, rhs));
2879b072555Sjeremylt   CeedVectorDestroy(&rhs_ceed);
2886c5df90dSjeremylt 
28943eb8658SJeremy L Thompson   // Create the error QFunction
2902b730f8bSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].error, bp_options[bp_choice].error_loc, &qf_error);
2919b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP);
2929b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE);
2932b730f8bSJeremy L Thompson   CeedQFunctionAddInput(qf_error, "qdata", ceed_data[fine_level]->q_data_size, CEED_EVAL_NONE);
29438f32c05Srezgarshakeri   CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_INTERP);
2956c5df90dSjeremylt 
2966c5df90dSjeremylt   // Create the error operator
2972b730f8bSJeremy L Thompson   CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_error);
2982b730f8bSJeremy L Thompson   CeedOperatorSetField(op_error, "u", ceed_data[fine_level]->elem_restr_u, ceed_data[fine_level]->basis_u, CEED_VECTOR_ACTIVE);
2992b730f8bSJeremy L Thompson   CeedOperatorSetField(op_error, "true_soln", ceed_data[fine_level]->elem_restr_u_i, CEED_BASIS_COLLOCATED, target);
3002b730f8bSJeremy L Thompson   CeedOperatorSetField(op_error, "qdata", ceed_data[fine_level]->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data[fine_level]->q_data);
3012b730f8bSJeremy L Thompson   CeedOperatorSetField(op_error, "error", ceed_data[fine_level]->elem_restr_u, ceed_data[fine_level]->basis_u, CEED_VECTOR_ACTIVE);
3026c5df90dSjeremylt 
3036c5df90dSjeremylt   // Calculate multiplicity
3049b072555Sjeremylt   for (int i = 0; i < num_levels; i++) {
3056c5df90dSjeremylt     PetscScalar *x;
3066c5df90dSjeremylt 
3076c5df90dSjeremylt     // CEED vector
3082b730f8bSJeremy L Thompson     PetscCall(VecZeroEntries(X_loc[i]));
3092b730f8bSJeremy L Thompson     PetscCall(VecGetArray(X_loc[i], &x));
3109b072555Sjeremylt     CeedVectorSetArray(ceed_data[i]->x_ceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
3116c5df90dSjeremylt 
3126c5df90dSjeremylt     // Multiplicity
3132b730f8bSJeremy L Thompson     CeedElemRestrictionGetMultiplicity(ceed_data[i]->elem_restr_u, ceed_data[i]->x_ceed);
3149b072555Sjeremylt     CeedVectorSyncArray(ceed_data[i]->x_ceed, CEED_MEM_HOST);
3156c5df90dSjeremylt 
3166c5df90dSjeremylt     // Restore vector
3172b730f8bSJeremy L Thompson     PetscCall(VecRestoreArray(X_loc[i], &x));
3186c5df90dSjeremylt 
3196c5df90dSjeremylt     // Creat mult vector
3202b730f8bSJeremy L Thompson     PetscCall(VecDuplicate(X_loc[i], &mult[i]));
3216c5df90dSjeremylt 
3226c5df90dSjeremylt     // Local-to-global
3232b730f8bSJeremy L Thompson     PetscCall(VecZeroEntries(X[i]));
3242b730f8bSJeremy L Thompson     PetscCall(DMLocalToGlobal(dm[i], X_loc[i], ADD_VALUES, X[i]));
3252b730f8bSJeremy L Thompson     PetscCall(VecZeroEntries(X_loc[i]));
3266c5df90dSjeremylt 
3276c5df90dSjeremylt     // Global-to-local
3282b730f8bSJeremy L Thompson     PetscCall(DMGlobalToLocal(dm[i], X[i], INSERT_VALUES, mult[i]));
3292b730f8bSJeremy L Thompson     PetscCall(VecZeroEntries(X[i]));
3306c5df90dSjeremylt 
3316c5df90dSjeremylt     // Multiplicity scaling
3322b730f8bSJeremy L Thompson     PetscCall(VecReciprocal(mult[i]));
3336c5df90dSjeremylt   }
3346c5df90dSjeremylt 
3356c5df90dSjeremylt   // Set up Mat
3369b072555Sjeremylt   for (int i = 0; i < num_levels; i++) {
3376c88e6a2Srezgarshakeri     // Set up apply operator context
3382b730f8bSJeremy L Thompson     PetscCall(SetupApplyOperatorCtx(comm, dm[i], ceed, ceed_data[i], X_loc[i], op_apply_ctx[i]));
3396c5df90dSjeremylt 
3406c5df90dSjeremylt     if (i > 0) {
341a97643b0Sjeremylt       // Prolongation/Restriction Operator
3422b730f8bSJeremy L Thompson       PetscCall(CeedLevelTransferSetup(dm[i - 1], ceed, i, num_comp_u, ceed_data, bp_options[bp_choice], mult[i]));
343d4d45553Srezgarshakeri       pr_restr_ctx[i]->comm        = comm;
344d4d45553Srezgarshakeri       pr_restr_ctx[i]->dmf         = dm[i];
345d4d45553Srezgarshakeri       pr_restr_ctx[i]->dmc         = dm[i - 1];
346d4d45553Srezgarshakeri       pr_restr_ctx[i]->loc_vec_c   = X_loc[i - 1];
347d4d45553Srezgarshakeri       pr_restr_ctx[i]->loc_vec_f   = op_apply_ctx[i]->Y_loc;
348d4d45553Srezgarshakeri       pr_restr_ctx[i]->mult_vec    = mult[i];
349d4d45553Srezgarshakeri       pr_restr_ctx[i]->ceed_vec_c  = op_apply_ctx[i - 1]->x_ceed;
350d4d45553Srezgarshakeri       pr_restr_ctx[i]->ceed_vec_f  = op_apply_ctx[i]->y_ceed;
351d4d45553Srezgarshakeri       pr_restr_ctx[i]->op_prolong  = ceed_data[i]->op_prolong;
352d4d45553Srezgarshakeri       pr_restr_ctx[i]->op_restrict = ceed_data[i]->op_restrict;
353d4d45553Srezgarshakeri       pr_restr_ctx[i]->ceed        = ceed;
3546c5df90dSjeremylt     }
3556c5df90dSjeremylt   }
3566c5df90dSjeremylt 
35753b04fa6SJed Brown   // Assemble coarse grid Jacobian for AMG (or other sparse matrix) solve
3582b730f8bSJeremy L Thompson   PetscCall(DMCreateMatrix(dm[0], &mat_coarse));
35953b04fa6SJed Brown 
3602b730f8bSJeremy L Thompson   PetscCall(PetscLogEventRegister("AssembleMatrix", MAT_CLASSID, &assemble_event));
361cffe6a52SJeremy L Thompson   {
362cffe6a52SJeremy L Thompson     // Assemble matrix analytically
3633047f789SJeremy L Thompson     PetscCount num_entries;
3643047f789SJeremy L Thompson     CeedInt   *rows, *cols;
36553b04fa6SJed Brown     CeedVector coo_values;
3662b730f8bSJeremy L Thompson     CeedOperatorLinearAssembleSymbolic(op_apply_ctx[0]->op, &num_entries, &rows, &cols);
36753b04fa6SJed Brown     ISLocalToGlobalMapping ltog_row, ltog_col;
3682b730f8bSJeremy L Thompson     PetscCall(MatGetLocalToGlobalMapping(mat_coarse, &ltog_row, &ltog_col));
3692b730f8bSJeremy L Thompson     PetscCall(ISLocalToGlobalMappingApply(ltog_row, num_entries, rows, rows));
3702b730f8bSJeremy L Thompson     PetscCall(ISLocalToGlobalMappingApply(ltog_col, num_entries, cols, cols));
3712b730f8bSJeremy L Thompson     PetscCall(MatSetPreallocationCOO(mat_coarse, num_entries, rows, cols));
37253b04fa6SJed Brown     free(rows);
37353b04fa6SJed Brown     free(cols);
37453b04fa6SJed Brown     CeedVectorCreate(ceed, num_entries, &coo_values);
3752b730f8bSJeremy L Thompson     PetscCall(PetscLogEventBegin(assemble_event, mat_coarse, 0, 0, 0));
376d4d45553Srezgarshakeri     CeedOperatorLinearAssemble(op_apply_ctx[0]->op, coo_values);
37753b04fa6SJed Brown     const CeedScalar *values;
37853b04fa6SJed Brown     CeedVectorGetArrayRead(coo_values, CEED_MEM_HOST, &values);
3792b730f8bSJeremy L Thompson     PetscCall(MatSetValuesCOO(mat_coarse, values, ADD_VALUES));
38053b04fa6SJed Brown     CeedVectorRestoreArrayRead(coo_values, &values);
3812b730f8bSJeremy L Thompson     PetscCall(PetscLogEventEnd(assemble_event, mat_coarse, 0, 0, 0));
38253b04fa6SJed Brown     CeedVectorDestroy(&coo_values);
38353b04fa6SJed Brown   }
38415ce0ef0Sjeremylt 
3856c5df90dSjeremylt   // Set up KSP
3862b730f8bSJeremy L Thompson   PetscCall(KSPCreate(comm, &ksp));
3876c5df90dSjeremylt   {
3882b730f8bSJeremy L Thompson     PetscCall(KSPSetType(ksp, KSPCG));
3892b730f8bSJeremy L Thompson     PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
3902b730f8bSJeremy L Thompson     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
3916c5df90dSjeremylt   }
3922b730f8bSJeremy L Thompson   PetscCall(KSPSetFromOptions(ksp));
3932b730f8bSJeremy L Thompson   PetscCall(KSPSetOperators(ksp, mat_O[fine_level], mat_O[fine_level]));
3946c5df90dSjeremylt 
3956c5df90dSjeremylt   // Set up PCMG
3962b730f8bSJeremy L Thompson   PetscCall(KSPGetPC(ksp, &pc));
3979b072555Sjeremylt   PCMGCycleType pcmg_cycle_type = PC_MG_CYCLE_V;
3986c5df90dSjeremylt   {
3992b730f8bSJeremy L Thompson     PetscCall(PCSetType(pc, PCMG));
4006c5df90dSjeremylt 
4016c5df90dSjeremylt     // PCMG levels
4022b730f8bSJeremy L Thompson     PetscCall(PCMGSetLevels(pc, num_levels, NULL));
4039b072555Sjeremylt     for (int i = 0; i < num_levels; i++) {
4046c5df90dSjeremylt       // Smoother
4056c5df90dSjeremylt       KSP smoother;
4066c5df90dSjeremylt       PC  smoother_pc;
4072b730f8bSJeremy L Thompson       PetscCall(PCMGGetSmoother(pc, i, &smoother));
4082b730f8bSJeremy L Thompson       PetscCall(KSPSetType(smoother, KSPCHEBYSHEV));
4092b730f8bSJeremy L Thompson       PetscCall(KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1));
4102b730f8bSJeremy L Thompson       PetscCall(KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE));
4112b730f8bSJeremy L Thompson       PetscCall(KSPSetOperators(smoother, mat_O[i], mat_O[i]));
4122b730f8bSJeremy L Thompson       PetscCall(KSPGetPC(smoother, &smoother_pc));
4132b730f8bSJeremy L Thompson       PetscCall(PCSetType(smoother_pc, PCJACOBI));
4142b730f8bSJeremy L Thompson       PetscCall(PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL));
4156c5df90dSjeremylt 
4166c5df90dSjeremylt       // Work vector
4179b072555Sjeremylt       if (i < num_levels - 1) {
4182b730f8bSJeremy L Thompson         PetscCall(PCMGSetX(pc, i, X[i]));
4196c5df90dSjeremylt       }
4206c5df90dSjeremylt 
4216c5df90dSjeremylt       // Level transfers
4226c5df90dSjeremylt       if (i > 0) {
4236c5df90dSjeremylt         // Interpolation
4242b730f8bSJeremy L Thompson         PetscCall(PCMGSetInterpolation(pc, i, mat_pr[i]));
4256c5df90dSjeremylt       }
4266c5df90dSjeremylt 
4276c5df90dSjeremylt       // Coarse solve
4286c5df90dSjeremylt       KSP coarse;
4296c5df90dSjeremylt       PC  coarse_pc;
4302b730f8bSJeremy L Thompson       PetscCall(PCMGGetCoarseSolve(pc, &coarse));
4312b730f8bSJeremy L Thompson       PetscCall(KSPSetType(coarse, KSPPREONLY));
4322b730f8bSJeremy L Thompson       PetscCall(KSPSetOperators(coarse, mat_coarse, mat_coarse));
43315ce0ef0Sjeremylt 
4342b730f8bSJeremy L Thompson       PetscCall(KSPGetPC(coarse, &coarse_pc));
4352b730f8bSJeremy L Thompson       PetscCall(PCSetType(coarse_pc, PCGAMG));
43615ce0ef0Sjeremylt 
4372b730f8bSJeremy L Thompson       PetscCall(KSPSetOptionsPrefix(coarse, "coarse_"));
4382b730f8bSJeremy L Thompson       PetscCall(PCSetOptionsPrefix(coarse_pc, "coarse_"));
4392b730f8bSJeremy L Thompson       PetscCall(KSPSetFromOptions(coarse));
4402b730f8bSJeremy L Thompson       PetscCall(PCSetFromOptions(coarse_pc));
4416c5df90dSjeremylt     }
4426c5df90dSjeremylt 
4436c5df90dSjeremylt     // PCMG options
4442b730f8bSJeremy L Thompson     PetscCall(PCMGSetType(pc, PC_MG_MULTIPLICATIVE));
4452b730f8bSJeremy L Thompson     PetscCall(PCMGSetNumberSmooth(pc, 3));
4462b730f8bSJeremy L Thompson     PetscCall(PCMGSetCycleType(pc, pcmg_cycle_type));
4476c5df90dSjeremylt   }
4486c5df90dSjeremylt 
4496c5df90dSjeremylt   // First run, if benchmarking
4506c5df90dSjeremylt   if (benchmark_mode) {
4512b730f8bSJeremy L Thompson     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1));
4522b730f8bSJeremy L Thompson     PetscCall(VecZeroEntries(X[fine_level]));
4536c5df90dSjeremylt     my_rt_start = MPI_Wtime();
4542b730f8bSJeremy L Thompson     PetscCall(KSPSolve(ksp, rhs, X[fine_level]));
4556c5df90dSjeremylt     my_rt = MPI_Wtime() - my_rt_start;
4562b730f8bSJeremy L Thompson     PetscCall(MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm));
4576c5df90dSjeremylt     // Set maxits based on first iteration timing
4586c5df90dSjeremylt     if (my_rt > 0.02) {
4592b730f8bSJeremy L Thompson       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5));
4606c5df90dSjeremylt     } else {
4612b730f8bSJeremy L Thompson       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20));
4626c5df90dSjeremylt     }
4636c5df90dSjeremylt   }
4646c5df90dSjeremylt 
4656c5df90dSjeremylt   // Timed solve
4662b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(X[fine_level]));
4672b730f8bSJeremy L Thompson   PetscCall(PetscBarrier((PetscObject)ksp));
46809a940d7Sjeremylt 
46909a940d7Sjeremylt   // -- Performance logging
4702b730f8bSJeremy L Thompson   PetscCall(PetscLogStageRegister("Solve Stage", &solve_stage));
4712b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePush(solve_stage));
47209a940d7Sjeremylt 
47309a940d7Sjeremylt   // -- Solve
4746c5df90dSjeremylt   my_rt_start = MPI_Wtime();
4752b730f8bSJeremy L Thompson   PetscCall(KSPSolve(ksp, rhs, X[fine_level]));
4766c5df90dSjeremylt   my_rt = MPI_Wtime() - my_rt_start;
4776c5df90dSjeremylt 
47809a940d7Sjeremylt   // -- Performance logging
4792b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePop());
48009a940d7Sjeremylt 
4816c5df90dSjeremylt   // Output results
4826c5df90dSjeremylt   {
4839b072555Sjeremylt     KSPType            ksp_type;
4849b072555Sjeremylt     PCMGType           pcmg_type;
4856c5df90dSjeremylt     KSPConvergedReason reason;
4866c5df90dSjeremylt     PetscReal          rnorm;
4876c5df90dSjeremylt     PetscInt           its;
4882b730f8bSJeremy L Thompson     PetscCall(KSPGetType(ksp, &ksp_type));
4892b730f8bSJeremy L Thompson     PetscCall(KSPGetConvergedReason(ksp, &reason));
4902b730f8bSJeremy L Thompson     PetscCall(KSPGetIterationNumber(ksp, &its));
4912b730f8bSJeremy L Thompson     PetscCall(KSPGetResidualNorm(ksp, &rnorm));
4922b730f8bSJeremy L Thompson     PetscCall(PCMGGetType(pc, &pcmg_type));
4936c5df90dSjeremylt     if (!test_mode || reason < 0 || rnorm > 1e-8) {
4942b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm,
4956c5df90dSjeremylt                             "  KSP:\n"
4966c5df90dSjeremylt                             "    KSP Type                                : %s\n"
4976c5df90dSjeremylt                             "    KSP Convergence                         : %s\n"
498a9b2c5ddSrezgarshakeri                             "    Total KSP Iterations                    : %" PetscInt_FMT "\n"
4996c5df90dSjeremylt                             "    Final rnorm                             : %e\n",
5002b730f8bSJeremy L Thompson                             ksp_type, KSPConvergedReasons[reason], its, (double)rnorm));
5012b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm,
5026c5df90dSjeremylt                             "  PCMG:\n"
5036c5df90dSjeremylt                             "    PCMG Type                               : %s\n"
5046c5df90dSjeremylt                             "    PCMG Cycle Type                         : %s\n",
5052b730f8bSJeremy L Thompson                             PCMGTypes[pcmg_type], PCMGCycleTypes[pcmg_cycle_type]));
5066c5df90dSjeremylt     }
5076c5df90dSjeremylt     if (!test_mode) {
5082b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm, "  Performance:\n"));
5096c5df90dSjeremylt     }
5106c5df90dSjeremylt     {
5116c88e6a2Srezgarshakeri       // Set up error operator context
5122b730f8bSJeremy L Thompson       PetscCall(SetupErrorOperatorCtx(comm, dm[fine_level], ceed, ceed_data[fine_level], X_loc[fine_level], op_error, op_error_ctx));
51338f32c05Srezgarshakeri       PetscScalar l2_error;
5142b730f8bSJeremy L Thompson       PetscCall(ComputeL2Error(X[fine_level], &l2_error, op_error_ctx));
5150a8fc04aSrezgarshakeri       PetscReal tol = 5e-2;
51638f32c05Srezgarshakeri       if (!test_mode || l2_error > tol) {
5172b730f8bSJeremy L Thompson         PetscCall(MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm));
5182b730f8bSJeremy L Thompson         PetscCall(MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm));
5192b730f8bSJeremy L Thompson         PetscCall(PetscPrintf(comm,
52038f32c05Srezgarshakeri                               "    L2 Error                                : %e\n"
5216c5df90dSjeremylt                               "    CG Solve Time                           : %g (%g) sec\n",
5222b730f8bSJeremy L Thompson                               (double)l2_error, rt_max, rt_min));
5236c5df90dSjeremylt       }
5246c5df90dSjeremylt     }
5256c5df90dSjeremylt     if (benchmark_mode && (!test_mode)) {
5262b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm, "    DoFs/Sec in CG                            : %g (%g) million\n", 1e-6 * g_size[fine_level] * its / rt_max,
5272b730f8bSJeremy L Thompson                             1e-6 * g_size[fine_level] * its / rt_min));
5286c5df90dSjeremylt     }
5296c5df90dSjeremylt   }
5306c5df90dSjeremylt 
5316c5df90dSjeremylt   if (write_solution) {
5329b072555Sjeremylt     PetscViewer vtk_viewer_soln;
5336c5df90dSjeremylt 
5342b730f8bSJeremy L Thompson     PetscCall(PetscViewerCreate(comm, &vtk_viewer_soln));
5352b730f8bSJeremy L Thompson     PetscCall(PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK));
5362b730f8bSJeremy L Thompson     PetscCall(PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"));
5372b730f8bSJeremy L Thompson     PetscCall(VecView(X[fine_level], vtk_viewer_soln));
5382b730f8bSJeremy L Thompson     PetscCall(PetscViewerDestroy(&vtk_viewer_soln));
5396c5df90dSjeremylt   }
5406c5df90dSjeremylt 
5416c5df90dSjeremylt   // Cleanup
5429b072555Sjeremylt   for (int i = 0; i < num_levels; i++) {
5432b730f8bSJeremy L Thompson     PetscCall(VecDestroy(&X[i]));
5442b730f8bSJeremy L Thompson     PetscCall(VecDestroy(&X_loc[i]));
5452b730f8bSJeremy L Thompson     PetscCall(VecDestroy(&mult[i]));
5462b730f8bSJeremy L Thompson     PetscCall(VecDestroy(&op_apply_ctx[i]->Y_loc));
5472b730f8bSJeremy L Thompson     PetscCall(MatDestroy(&mat_O[i]));
5482b730f8bSJeremy L Thompson     PetscCall(PetscFree(op_apply_ctx[i]));
5496c5df90dSjeremylt     if (i > 0) {
5502b730f8bSJeremy L Thompson       PetscCall(MatDestroy(&mat_pr[i]));
5512b730f8bSJeremy L Thompson       PetscCall(PetscFree(pr_restr_ctx[i]));
5526c5df90dSjeremylt     }
5532b730f8bSJeremy L Thompson     PetscCall(CeedDataDestroy(i, ceed_data[i]));
5542b730f8bSJeremy L Thompson     PetscCall(DMDestroy(&dm[i]));
5556c5df90dSjeremylt   }
5562b730f8bSJeremy L Thompson   PetscCall(PetscFree(level_degrees));
5572b730f8bSJeremy L Thompson   PetscCall(PetscFree(dm));
5582b730f8bSJeremy L Thompson   PetscCall(PetscFree(X));
5592b730f8bSJeremy L Thompson   PetscCall(PetscFree(X_loc));
5602b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&op_error_ctx->Y_loc));
5612b730f8bSJeremy L Thompson   PetscCall(PetscFree(mult));
5622b730f8bSJeremy L Thompson   PetscCall(PetscFree(mat_O));
5632b730f8bSJeremy L Thompson   PetscCall(PetscFree(mat_pr));
5642b730f8bSJeremy L Thompson   PetscCall(PetscFree(ceed_data));
5652b730f8bSJeremy L Thompson   PetscCall(PetscFree(op_apply_ctx));
5662b730f8bSJeremy L Thompson   PetscCall(PetscFree(op_error_ctx));
5672b730f8bSJeremy L Thompson   PetscCall(PetscFree(pr_restr_ctx));
5682b730f8bSJeremy L Thompson   PetscCall(PetscFree(l_size));
5692b730f8bSJeremy L Thompson   PetscCall(PetscFree(xl_size));
5702b730f8bSJeremy L Thompson   PetscCall(PetscFree(g_size));
5712b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&rhs));
5722b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&rhs_loc));
5732b730f8bSJeremy L Thompson   PetscCall(MatDestroy(&mat_coarse));
5742b730f8bSJeremy L Thompson   PetscCall(KSPDestroy(&ksp));
5752b730f8bSJeremy L Thompson   PetscCall(DMDestroy(&dm_orig));
5766c5df90dSjeremylt   CeedVectorDestroy(&target);
5779b072555Sjeremylt   CeedQFunctionDestroy(&qf_error);
5789b072555Sjeremylt   CeedOperatorDestroy(&op_error);
5796c5df90dSjeremylt   CeedDestroy(&ceed);
5806c5df90dSjeremylt   return PetscFinalize();
5816c5df90dSjeremylt }
582