xref: /libCEED/examples/petsc/multigrid.c (revision a9b2c5dd65841e8c467007aaab2d4e9349426d5a)
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 //
106c5df90dSjeremylt // This example demonstrates a simple usage of libCEED with PETSc to solve the
116c5df90dSjeremylt // CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps.
126c5df90dSjeremylt //
136c5df90dSjeremylt // The code uses higher level communication protocols in DMPlex.
146c5df90dSjeremylt //
156c5df90dSjeremylt // Build with:
166c5df90dSjeremylt //
176c5df90dSjeremylt //     make multigrid [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
186c5df90dSjeremylt //
196c5df90dSjeremylt // Sample runs:
206c5df90dSjeremylt //
216c5df90dSjeremylt //     multigrid -problem bp3
2228688798Sjeremylt //     multigrid -problem bp4
2328688798Sjeremylt //     multigrid -problem bp5 -ceed /cpu/self
246c5df90dSjeremylt //     multigrid -problem bp6 -ceed /gpu/cuda
256c5df90dSjeremylt //
266c5df90dSjeremylt //TESTARGS -ceed {ceed_resource} -test -problem bp3 -degree 3
27*a9b2c5ddSrezgarshakeri //TESTARGS -ceed {ceed_resource} -test -problem bp3 -degree 3 -simplex
286c5df90dSjeremylt 
296c5df90dSjeremylt /// @file
306c5df90dSjeremylt /// CEED BPs 1-6 multigrid example using PETSc
316c5df90dSjeremylt const char help[] = "Solve CEED BPs using p-multigrid with PETSc and DMPlex\n";
326c5df90dSjeremylt 
33636cccdbSjeremylt #include <stdbool.h>
34636cccdbSjeremylt #include <string.h>
35636cccdbSjeremylt #include <ceed.h>
36636cccdbSjeremylt #include <petsc.h>
37636cccdbSjeremylt #include <petscdmplex.h>
38636cccdbSjeremylt #include <petscksp.h>
39636cccdbSjeremylt #include <petscsys.h>
40636cccdbSjeremylt 
41e83e87a5Sjeremylt #include "bps.h"
42636cccdbSjeremylt #include "include/bpsproblemdata.h"
43636cccdbSjeremylt #include "include/petscutils.h"
44b8962995SJeremy L Thompson #include "include/petscversion.h"
45636cccdbSjeremylt #include "include/matops.h"
46636cccdbSjeremylt #include "include/structs.h"
47636cccdbSjeremylt #include "include/libceedsetup.h"
48636cccdbSjeremylt 
49636cccdbSjeremylt #if PETSC_VERSION_LT(3,12,0)
50636cccdbSjeremylt #ifdef PETSC_HAVE_CUDA
51636cccdbSjeremylt #include <petsccuda.h>
52636cccdbSjeremylt // Note: With PETSc prior to version 3.12.0, providing the source path to
53636cccdbSjeremylt //       include 'cublas_v2.h' will be needed to use 'petsccuda.h'.
54636cccdbSjeremylt #endif
55636cccdbSjeremylt #endif
566c5df90dSjeremylt 
576c5df90dSjeremylt int main(int argc, char **argv) {
586c5df90dSjeremylt   PetscInt ierr;
596c5df90dSjeremylt   MPI_Comm comm;
60cb0b5415Sjeremylt   char filename[PETSC_MAX_PATH_LEN],
619b072555Sjeremylt        ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self";
626c5df90dSjeremylt   double my_rt_start, my_rt, rt_min, rt_max;
639b072555Sjeremylt   PetscInt degree = 3, q_extra, *l_size, *xl_size, *g_size, dim = 3, fine_level,
649b072555Sjeremylt            mesh_elem[3] = {3, 3, 3}, num_comp_u = 1, num_levels = degree, *level_degrees;
656c5df90dSjeremylt   PetscScalar *r;
66cfa59c5bSRey   PetscScalar eps = 1.0;
67de1229c5Srezgarshakeri   PetscBool test_mode, benchmark_mode, read_mesh, write_solution, simplex;
689b072555Sjeremylt   PetscLogStage solve_stage;
6905b9c820SJed Brown   PetscLogEvent assemble_event;
709b072555Sjeremylt   DM  *dm, dm_orig;
716c5df90dSjeremylt   KSP ksp;
726c5df90dSjeremylt   PC pc;
739b072555Sjeremylt   Mat *mat_O, *mat_pr, mat_coarse;
749b072555Sjeremylt   Vec *X, *X_loc, *mult, rhs, rhs_loc;
759b072555Sjeremylt   PetscMemType mem_type;
769b072555Sjeremylt   UserO *user_O;
779b072555Sjeremylt   UserProlongRestr *user_pr;
786c5df90dSjeremylt   Ceed ceed;
799b072555Sjeremylt   CeedData *ceed_data;
809b072555Sjeremylt   CeedVector rhs_ceed, target;
811c9a79dbSrezgarshakeri   CeedQFunction qf_error;
829b072555Sjeremylt   CeedOperator op_error;
839b072555Sjeremylt   BPType bp_choice;
849b072555Sjeremylt   CoarsenType coarsen;
856c5df90dSjeremylt 
866c5df90dSjeremylt   ierr = PetscInitialize(&argc, &argv, NULL, help);
876c5df90dSjeremylt   if (ierr) return ierr;
886c5df90dSjeremylt   comm = PETSC_COMM_WORLD;
896c5df90dSjeremylt 
906c5df90dSjeremylt   // Parse command line options
9167490bc6SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
929b072555Sjeremylt   bp_choice = CEED_BP3;
936c5df90dSjeremylt   ierr = PetscOptionsEnum("-problem",
946c5df90dSjeremylt                           "CEED benchmark problem to solve", NULL,
959b072555Sjeremylt                           bp_types, (PetscEnum)bp_choice, (PetscEnum *)&bp_choice,
966c5df90dSjeremylt                           NULL); CHKERRQ(ierr);
979b072555Sjeremylt   num_comp_u = bp_options[bp_choice].num_comp_u;
986c5df90dSjeremylt   test_mode = PETSC_FALSE;
996c5df90dSjeremylt   ierr = PetscOptionsBool("-test",
1006c5df90dSjeremylt                           "Testing mode (do not print unless error is large)",
1016c5df90dSjeremylt                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
1026c5df90dSjeremylt   benchmark_mode = PETSC_FALSE;
1036c5df90dSjeremylt   ierr = PetscOptionsBool("-benchmark",
1046c5df90dSjeremylt                           "Benchmarking mode (prints benchmark statistics)",
1056c5df90dSjeremylt                           NULL, benchmark_mode, &benchmark_mode, NULL);
1066c5df90dSjeremylt   CHKERRQ(ierr);
1076c5df90dSjeremylt   write_solution = PETSC_FALSE;
1086c5df90dSjeremylt   ierr = PetscOptionsBool("-write_solution",
1096c5df90dSjeremylt                           "Write solution for visualization",
1106c5df90dSjeremylt                           NULL, write_solution, &write_solution, NULL);
1116c5df90dSjeremylt   CHKERRQ(ierr);
112de1229c5Srezgarshakeri   simplex = PETSC_FALSE;
113de1229c5Srezgarshakeri   ierr = PetscOptionsBool("-simplex",
114de1229c5Srezgarshakeri                           "Element topology (default:hex)",
115de1229c5Srezgarshakeri                           NULL, simplex, &simplex, NULL);
116de1229c5Srezgarshakeri   CHKERRQ(ierr);
117cfa59c5bSRey   ierr = PetscOptionsScalar("-eps",
118cfa59c5bSRey                             "Epsilon parameter for Kershaw mesh transformation",
119cfa59c5bSRey                             NULL, eps, &eps, NULL);
1207578c821SJeremy L Thompson   if (eps > 1 || eps <= 0) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
12108140895SJed Brown                                      "-eps %g must be (0,1]", (double)PetscRealPart(eps));
1226c5df90dSjeremylt   degree = test_mode ? 3 : 2;
1236c5df90dSjeremylt   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
1246c5df90dSjeremylt                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
1257578c821SJeremy L Thompson   if (degree < 1) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
12608140895SJed Brown                             "-degree %" PetscInt_FMT " must be at least 1", degree);
1279b072555Sjeremylt   q_extra = bp_options[bp_choice].q_extra;
1289b072555Sjeremylt   ierr = PetscOptionsInt("-q_extra", "Number of extra quadrature points",
1299b072555Sjeremylt                          NULL, q_extra, &q_extra, NULL); CHKERRQ(ierr);
1306c5df90dSjeremylt   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
1319b072555Sjeremylt                             NULL, ceed_resource, ceed_resource,
1329b072555Sjeremylt                             sizeof(ceed_resource), NULL); CHKERRQ(ierr);
1336c5df90dSjeremylt   coarsen = COARSEN_UNIFORM;
1346c5df90dSjeremylt   ierr = PetscOptionsEnum("-coarsen",
1356c5df90dSjeremylt                           "Coarsening strategy to use", NULL,
1369b072555Sjeremylt                           coarsen_types, (PetscEnum)coarsen,
1376c5df90dSjeremylt                           (PetscEnum *)&coarsen, NULL); CHKERRQ(ierr);
138cb32e2e7SValeria Barra   read_mesh = PETSC_FALSE;
1396c5df90dSjeremylt   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL,
1406c5df90dSjeremylt                             filename, filename, sizeof(filename), &read_mesh);
1416c5df90dSjeremylt   CHKERRQ(ierr);
1426c5df90dSjeremylt   if (!read_mesh) {
1436c5df90dSjeremylt     PetscInt tmp = dim;
1446c5df90dSjeremylt     ierr = PetscOptionsIntArray("-cells","Number of cells per dimension", NULL,
1459b072555Sjeremylt                                 mesh_elem, &tmp, NULL); CHKERRQ(ierr);
1466c5df90dSjeremylt   }
14767490bc6SJeremy L Thompson   PetscOptionsEnd();
1486c5df90dSjeremylt 
1499396343dSjeremylt   // Set up libCEED
1509b072555Sjeremylt   CeedInit(ceed_resource, &ceed);
1519b072555Sjeremylt   CeedMemType mem_type_backend;
1529b072555Sjeremylt   CeedGetPreferredMemType(ceed, &mem_type_backend);
1539396343dSjeremylt 
1546c5df90dSjeremylt   // Setup DM
1556c5df90dSjeremylt   if (read_mesh) {
1567ed3e4cdSJeremy L Thompson     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, NULL, PETSC_TRUE,
1577ed3e4cdSJeremy L Thompson                                 &dm_orig);
1586c5df90dSjeremylt     CHKERRQ(ierr);
1596c5df90dSjeremylt   } else {
160de1229c5Srezgarshakeri     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, simplex, mesh_elem, NULL,
1619b072555Sjeremylt                                NULL, NULL, PETSC_TRUE, &dm_orig); CHKERRQ(ierr);
1626c5df90dSjeremylt   }
1636c5df90dSjeremylt 
1649b072555Sjeremylt   VecType vec_type;
1659b072555Sjeremylt   switch (mem_type_backend) {
1669b072555Sjeremylt   case CEED_MEM_HOST: vec_type = VECSTANDARD; break;
167b68a8d79SJed Brown   case CEED_MEM_DEVICE: {
168b68a8d79SJed Brown     const char *resolved;
169b68a8d79SJed Brown     CeedGetResource(ceed, &resolved);
1709b072555Sjeremylt     if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA;
171ca2d516cSJed Brown     else if (strstr(resolved, "/gpu/hip/occa"))
1729b072555Sjeremylt       vec_type = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678
1739b072555Sjeremylt     else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP;
1749b072555Sjeremylt     else vec_type = VECSTANDARD;
175b68a8d79SJed Brown   }
176b68a8d79SJed Brown   }
1779b072555Sjeremylt   ierr = DMSetVecType(dm_orig, vec_type); CHKERRQ(ierr);
1789b072555Sjeremylt   ierr = DMSetFromOptions(dm_orig); CHKERRQ(ierr);
1793fc8a154SJed Brown   ierr = DMViewFromOptions(dm_orig, NULL, "-dm_view"); CHKERRQ(ierr);
1803fc8a154SJed Brown 
1813fc8a154SJed Brown   // Apply Kershaw mesh transformation
1823fc8a154SJed Brown   ierr = Kershaw(dm_orig, eps); CHKERRQ(ierr);
183b68a8d79SJed Brown 
1846c5df90dSjeremylt   // Allocate arrays for PETSc objects for each level
1856c5df90dSjeremylt   switch (coarsen) {
1866c5df90dSjeremylt   case COARSEN_UNIFORM:
1879b072555Sjeremylt     num_levels = degree;
1886c5df90dSjeremylt     break;
189dc7d240cSValeria Barra   case COARSEN_LOGARITHMIC:
1909b072555Sjeremylt     num_levels = ceil(log(degree)/log(2)) + 1;
1916c5df90dSjeremylt     break;
1926c5df90dSjeremylt   }
1939b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &level_degrees); CHKERRQ(ierr);
1949b072555Sjeremylt   fine_level = num_levels - 1;
19561608365Sjeremylt 
1966c5df90dSjeremylt   switch (coarsen) {
1976c5df90dSjeremylt   case COARSEN_UNIFORM:
1989b072555Sjeremylt     for (int i=0; i<num_levels; i++) level_degrees[i] = i + 1;
1996c5df90dSjeremylt     break;
200dc7d240cSValeria Barra   case COARSEN_LOGARITHMIC:
2019b072555Sjeremylt     for (int i=0; i<num_levels - 1; i++) level_degrees[i] = pow(2,i);
2029b072555Sjeremylt     level_degrees[fine_level] = degree;
2036c5df90dSjeremylt     break;
2046c5df90dSjeremylt   }
2059b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &dm); CHKERRQ(ierr);
2069b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &X); CHKERRQ(ierr);
2079b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &X_loc); CHKERRQ(ierr);
2089b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &mult); CHKERRQ(ierr);
2099b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &user_O); CHKERRQ(ierr);
2109b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &user_pr); CHKERRQ(ierr);
2119b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &mat_O); CHKERRQ(ierr);
2129b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &mat_pr); CHKERRQ(ierr);
2139b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &l_size); CHKERRQ(ierr);
2149b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &xl_size); CHKERRQ(ierr);
2159b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &g_size); CHKERRQ(ierr);
2166c5df90dSjeremylt 
2171c9a79dbSrezgarshakeri   PetscInt c_start, c_end;
2181c9a79dbSrezgarshakeri   ierr = DMPlexGetHeightStratum(dm_orig, 0, &c_start, &c_end); CHKERRQ(ierr);
2191c9a79dbSrezgarshakeri   DMPolytopeType  cell_type;
2201c9a79dbSrezgarshakeri   ierr = DMPlexGetCellType(dm_orig, c_start, &cell_type); CHKERRQ(ierr);
2211c9a79dbSrezgarshakeri   CeedElemTopology elem_topo = ElemTopologyP2C(cell_type);
2221c9a79dbSrezgarshakeri 
2236c5df90dSjeremylt   // Setup DM and Operator Mat Shells for each level
2249b072555Sjeremylt   for (CeedInt i=0; i<num_levels; i++) {
2256c5df90dSjeremylt     // Create DM
2269b072555Sjeremylt     ierr = DMClone(dm_orig, &dm[i]); CHKERRQ(ierr);
2279b072555Sjeremylt     ierr = DMGetVecType(dm_orig, &vec_type); CHKERRQ(ierr);
2289b072555Sjeremylt     ierr = DMSetVecType(dm[i], vec_type); CHKERRQ(ierr);
229e83e87a5Sjeremylt     PetscInt dim;
230e83e87a5Sjeremylt     ierr = DMGetDimension(dm[i], &dim); CHKERRQ(ierr);
231de1229c5Srezgarshakeri     ierr = SetupDMByDegree(dm[i], level_degrees[i], q_extra, num_comp_u, dim,
2329b072555Sjeremylt                            bp_options[bp_choice].enforce_bc, bp_options[bp_choice].bc_func);
2336c5df90dSjeremylt     CHKERRQ(ierr);
2346c5df90dSjeremylt 
2356c5df90dSjeremylt     // Create vectors
2366c5df90dSjeremylt     ierr = DMCreateGlobalVector(dm[i], &X[i]); CHKERRQ(ierr);
2379b072555Sjeremylt     ierr = VecGetLocalSize(X[i], &l_size[i]); CHKERRQ(ierr);
2389b072555Sjeremylt     ierr = VecGetSize(X[i], &g_size[i]); CHKERRQ(ierr);
2399b072555Sjeremylt     ierr = DMCreateLocalVector(dm[i], &X_loc[i]); CHKERRQ(ierr);
2409b072555Sjeremylt     ierr = VecGetSize(X_loc[i], &xl_size[i]); CHKERRQ(ierr);
2416c5df90dSjeremylt 
2426c5df90dSjeremylt     // Operator
2439b072555Sjeremylt     ierr = PetscMalloc1(1, &user_O[i]); CHKERRQ(ierr);
2449b072555Sjeremylt     ierr = MatCreateShell(comm, l_size[i], l_size[i], g_size[i], g_size[i],
2459b072555Sjeremylt                           user_O[i], &mat_O[i]); CHKERRQ(ierr);
2469b072555Sjeremylt     ierr = MatShellSetOperation(mat_O[i], MATOP_MULT,
247ce74dcefSjeremylt                                 (void(*)(void))MatMult_Ceed); CHKERRQ(ierr);
2489b072555Sjeremylt     ierr = MatShellSetOperation(mat_O[i], MATOP_GET_DIAGONAL,
249ce74dcefSjeremylt                                 (void(*)(void))MatGetDiag); CHKERRQ(ierr);
2509b072555Sjeremylt     ierr = MatShellSetVecType(mat_O[i], vec_type); CHKERRQ(ierr);
2516c5df90dSjeremylt 
2526c5df90dSjeremylt     // Level transfers
2536c5df90dSjeremylt     if (i > 0) {
2546c5df90dSjeremylt       // Interp
2559b072555Sjeremylt       ierr = PetscMalloc1(1, &user_pr[i]); CHKERRQ(ierr);
2569b072555Sjeremylt       ierr = MatCreateShell(comm, l_size[i], l_size[i-1], g_size[i], g_size[i-1],
2579b072555Sjeremylt                             user_pr[i], &mat_pr[i]); CHKERRQ(ierr);
2589b072555Sjeremylt       ierr = MatShellSetOperation(mat_pr[i], MATOP_MULT,
259a97643b0Sjeremylt                                   (void(*)(void))MatMult_Prolong);
2606c5df90dSjeremylt       CHKERRQ(ierr);
2619b072555Sjeremylt       ierr = MatShellSetOperation(mat_pr[i], MATOP_MULT_TRANSPOSE,
2626c5df90dSjeremylt                                   (void(*)(void))MatMult_Restrict);
2636c5df90dSjeremylt       CHKERRQ(ierr);
2649b072555Sjeremylt       ierr = MatShellSetVecType(mat_pr[i], vec_type); CHKERRQ(ierr);
2656c5df90dSjeremylt     }
2666c5df90dSjeremylt   }
2679b072555Sjeremylt   ierr = VecDuplicate(X[fine_level], &rhs); CHKERRQ(ierr);
2686c5df90dSjeremylt 
2696c5df90dSjeremylt   // Print global grid information
2706c5df90dSjeremylt   if (!test_mode) {
2719b072555Sjeremylt     PetscInt P = degree + 1, Q = P + q_extra;
2729396343dSjeremylt 
2739b072555Sjeremylt     const char *used_resource;
2749b072555Sjeremylt     CeedGetResource(ceed, &used_resource);
2759396343dSjeremylt 
2769b072555Sjeremylt     ierr = VecGetType(X[0], &vec_type); CHKERRQ(ierr);
2779396343dSjeremylt 
2786c5df90dSjeremylt     ierr = PetscPrintf(comm,
279990fdeb6SJeremy L Thompson                        "\n-- CEED Benchmark Problem %" CeedInt_FMT " -- libCEED + PETSc + PCMG --\n"
2809396343dSjeremylt                        "  PETSc:\n"
2819396343dSjeremylt                        "    PETSc Vec Type                          : %s\n"
2826c5df90dSjeremylt                        "  libCEED:\n"
2836c5df90dSjeremylt                        "    libCEED Backend                         : %s\n"
2849396343dSjeremylt                        "    libCEED Backend MemType                 : %s\n"
2856c5df90dSjeremylt                        "  Mesh:\n"
286990fdeb6SJeremy L Thompson                        "    Number of 1D Basis Nodes (p)            : %" CeedInt_FMT "\n"
287990fdeb6SJeremy L Thompson                        "    Number of 1D Quadrature Points (q)      : %" CeedInt_FMT "\n"
288de1229c5Srezgarshakeri                        "    Additional quadrature points (q_extra)  : %" CeedInt_FMT "\n"
28908140895SJed Brown                        "    Global Nodes                            : %" PetscInt_FMT "\n"
29008140895SJed Brown                        "    Owned Nodes                             : %" PetscInt_FMT "\n"
29108140895SJed Brown                        "    DoF per node                            : %" PetscInt_FMT "\n"
2926c5df90dSjeremylt                        "  Multigrid:\n"
293990fdeb6SJeremy L Thompson                        "    Number of Levels                        : %" CeedInt_FMT "\n",
2949b072555Sjeremylt                        bp_choice+1, vec_type, used_resource,
295de1229c5Srezgarshakeri                        CeedMemTypes[mem_type_backend], P, Q, q_extra,
296de1229c5Srezgarshakeri                        g_size[fine_level]/num_comp_u, l_size[fine_level]/num_comp_u,
2971c9a79dbSrezgarshakeri                        num_comp_u, CeedElemTopologies[elem_topo],
2981c9a79dbSrezgarshakeri                        num_levels); CHKERRQ(ierr);
2996c5df90dSjeremylt   }
3006c5df90dSjeremylt 
3016c5df90dSjeremylt   // Create RHS vector
3029b072555Sjeremylt   ierr = VecDuplicate(X_loc[fine_level], &rhs_loc); CHKERRQ(ierr);
3039b072555Sjeremylt   ierr = VecZeroEntries(rhs_loc); CHKERRQ(ierr);
3049b072555Sjeremylt   ierr = VecGetArrayAndMemType(rhs_loc, &r, &mem_type); CHKERRQ(ierr);
3059b072555Sjeremylt   CeedVectorCreate(ceed, xl_size[fine_level], &rhs_ceed);
3069b072555Sjeremylt   CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r);
3076c5df90dSjeremylt 
3086c5df90dSjeremylt   // Set up libCEED operators on each level
3099b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &ceed_data); CHKERRQ(ierr);
310990fdeb6SJeremy L Thompson   for (PetscInt i=0; i<num_levels; i++) {
3116c5df90dSjeremylt     // Print level information
3129b072555Sjeremylt     if (!test_mode && (i == 0 || i == fine_level)) {
31308140895SJed Brown       ierr = PetscPrintf(comm,"    Level %" PetscInt_FMT " (%s):\n"
314990fdeb6SJeremy L Thompson                          "      Number of 1D Basis Nodes (p)     : %" CeedInt_FMT "\n"
31508140895SJed Brown                          "      Global Nodes                     : %" PetscInt_FMT "\n"
31608140895SJed Brown                          "      Owned Nodes                      : %" PetscInt_FMT "\n",
3179b072555Sjeremylt                          i, (i? "fine" : "coarse"), level_degrees[i] + 1,
3189b072555Sjeremylt                          g_size[i]/num_comp_u, l_size[i]/num_comp_u); CHKERRQ(ierr);
3196c5df90dSjeremylt     }
3209b072555Sjeremylt     ierr = PetscMalloc1(1, &ceed_data[i]); CHKERRQ(ierr);
3219b072555Sjeremylt     ierr = SetupLibceedByDegree(dm[i], ceed, level_degrees[i], dim, q_extra,
3229b072555Sjeremylt                                 dim, num_comp_u, g_size[i], xl_size[i], bp_options[bp_choice],
3239b072555Sjeremylt                                 ceed_data[i], i==(fine_level), rhs_ceed, &target);
324f9342789SJeremy L Thompson     CHKERRQ(ierr);
3256c5df90dSjeremylt   }
3266c5df90dSjeremylt 
3276c5df90dSjeremylt   // Gather RHS
3289b072555Sjeremylt   CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL);
3299b072555Sjeremylt   ierr = VecRestoreArrayAndMemType(rhs_loc, &r); CHKERRQ(ierr);
3306c5df90dSjeremylt   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
3319b072555Sjeremylt   ierr = DMLocalToGlobal(dm[fine_level], rhs_loc, ADD_VALUES, rhs); CHKERRQ(ierr);
3329b072555Sjeremylt   CeedVectorDestroy(&rhs_ceed);
3336c5df90dSjeremylt 
33443eb8658SJeremy L Thompson   // Create the error QFunction
3359b072555Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].error,
3369b072555Sjeremylt                               bp_options[bp_choice].error_loc, &qf_error);
3379b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP);
3389b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE);
3399b072555Sjeremylt   CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_NONE);
3406c5df90dSjeremylt 
3416c5df90dSjeremylt   // Create the error operator
3429b072555Sjeremylt   CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
3439b072555Sjeremylt                      &op_error);
3449b072555Sjeremylt   CeedOperatorSetField(op_error, "u", ceed_data[fine_level]->elem_restr_u,
3459b072555Sjeremylt                        ceed_data[fine_level]->basis_u, CEED_VECTOR_ACTIVE);
3469b072555Sjeremylt   CeedOperatorSetField(op_error, "true_soln",
3479b072555Sjeremylt                        ceed_data[fine_level]->elem_restr_u_i,
348a8d32208Sjeremylt                        CEED_BASIS_COLLOCATED, target);
3499b072555Sjeremylt   CeedOperatorSetField(op_error, "error", ceed_data[fine_level]->elem_restr_u_i,
350a8d32208Sjeremylt                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
3516c5df90dSjeremylt 
3526c5df90dSjeremylt   // Calculate multiplicity
3539b072555Sjeremylt   for (int i=0; i<num_levels; i++) {
3546c5df90dSjeremylt     PetscScalar *x;
3556c5df90dSjeremylt 
3566c5df90dSjeremylt     // CEED vector
3579b072555Sjeremylt     ierr = VecZeroEntries(X_loc[i]); CHKERRQ(ierr);
3589b072555Sjeremylt     ierr = VecGetArray(X_loc[i], &x); CHKERRQ(ierr);
3599b072555Sjeremylt     CeedVectorSetArray(ceed_data[i]->x_ceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
3606c5df90dSjeremylt 
3616c5df90dSjeremylt     // Multiplicity
3629b072555Sjeremylt     CeedElemRestrictionGetMultiplicity(ceed_data[i]->elem_restr_u,
3639b072555Sjeremylt                                        ceed_data[i]->x_ceed);
3649b072555Sjeremylt     CeedVectorSyncArray(ceed_data[i]->x_ceed, CEED_MEM_HOST);
3656c5df90dSjeremylt 
3666c5df90dSjeremylt     // Restore vector
3679b072555Sjeremylt     ierr = VecRestoreArray(X_loc[i], &x); CHKERRQ(ierr);
3686c5df90dSjeremylt 
3696c5df90dSjeremylt     // Creat mult vector
3709b072555Sjeremylt     ierr = VecDuplicate(X_loc[i], &mult[i]); CHKERRQ(ierr);
3716c5df90dSjeremylt 
3726c5df90dSjeremylt     // Local-to-global
3736c5df90dSjeremylt     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
3749b072555Sjeremylt     ierr = DMLocalToGlobal(dm[i], X_loc[i], ADD_VALUES, X[i]);
3756c5df90dSjeremylt     CHKERRQ(ierr);
3769b072555Sjeremylt     ierr = VecZeroEntries(X_loc[i]); CHKERRQ(ierr);
3776c5df90dSjeremylt 
3786c5df90dSjeremylt     // Global-to-local
379483f8b0dSjeremylt     ierr = DMGlobalToLocal(dm[i], X[i], INSERT_VALUES, mult[i]);
3806c5df90dSjeremylt     CHKERRQ(ierr);
3816c5df90dSjeremylt     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
3826c5df90dSjeremylt 
3836c5df90dSjeremylt     // Multiplicity scaling
3846c5df90dSjeremylt     ierr = VecReciprocal(mult[i]);
3856c5df90dSjeremylt   }
3866c5df90dSjeremylt 
3876c5df90dSjeremylt   // Set up Mat
3889b072555Sjeremylt   for (int i=0; i<num_levels; i++) {
389226c3a8fSjeremylt     // User Operator
3909b072555Sjeremylt     user_O[i]->comm = comm;
3919b072555Sjeremylt     user_O[i]->dm = dm[i];
3929b072555Sjeremylt     user_O[i]->X_loc = X_loc[i];
3939b072555Sjeremylt     ierr = VecDuplicate(X_loc[i], &user_O[i]->Y_loc); CHKERRQ(ierr);
3949b072555Sjeremylt     user_O[i]->x_ceed = ceed_data[i]->x_ceed;
3959b072555Sjeremylt     user_O[i]->y_ceed = ceed_data[i]->y_ceed;
3969b072555Sjeremylt     user_O[i]->op = ceed_data[i]->op_apply;
3979b072555Sjeremylt     user_O[i]->ceed = ceed;
3986c5df90dSjeremylt 
3996c5df90dSjeremylt     if (i > 0) {
400a97643b0Sjeremylt       // Prolongation/Restriction Operator
4011c9a79dbSrezgarshakeri       ierr = CeedLevelTransferSetup(dm[i], ceed, i, num_comp_u, ceed_data,
4021c9a79dbSrezgarshakeri                                     mult[i]); CHKERRQ(ierr);
4039b072555Sjeremylt       user_pr[i]->comm = comm;
4049b072555Sjeremylt       user_pr[i]->dmf = dm[i];
4059b072555Sjeremylt       user_pr[i]->dmc = dm[i-1];
4069b072555Sjeremylt       user_pr[i]->loc_vec_c = X_loc[i-1];
4079b072555Sjeremylt       user_pr[i]->loc_vec_f = user_O[i]->Y_loc;
4089b072555Sjeremylt       user_pr[i]->mult_vec = mult[i];
4099b072555Sjeremylt       user_pr[i]->ceed_vec_c = user_O[i-1]->x_ceed;
4109b072555Sjeremylt       user_pr[i]->ceed_vec_f = user_O[i]->y_ceed;
4119b072555Sjeremylt       user_pr[i]->op_prolong = ceed_data[i]->op_prolong;
4129b072555Sjeremylt       user_pr[i]->op_restrict = ceed_data[i]->op_restrict;
4139b072555Sjeremylt       user_pr[i]->ceed = ceed;
4146c5df90dSjeremylt     }
4156c5df90dSjeremylt   }
4166c5df90dSjeremylt 
41753b04fa6SJed Brown   // Assemble coarse grid Jacobian for AMG (or other sparse matrix) solve
41853b04fa6SJed Brown   ierr = DMCreateMatrix(dm[0], &mat_coarse); CHKERRQ(ierr);
41953b04fa6SJed Brown 
42005b9c820SJed Brown   ierr = PetscLogEventRegister("AssembleMatrix", MAT_CLASSID, &assemble_event);
42105b9c820SJed Brown   CHKERRQ(ierr);
422cffe6a52SJeremy L Thompson   {
423cffe6a52SJeremy L Thompson     // Assemble matrix analytically
4243047f789SJeremy L Thompson     PetscCount num_entries;
4253047f789SJeremy L Thompson     CeedInt *rows, *cols;
42653b04fa6SJed Brown     CeedVector coo_values;
42753b04fa6SJed Brown     CeedOperatorLinearAssembleSymbolic(user_O[0]->op, &num_entries, &rows, &cols);
42853b04fa6SJed Brown     ISLocalToGlobalMapping ltog_row, ltog_col;
42953b04fa6SJed Brown     ierr = MatGetLocalToGlobalMapping(mat_coarse, &ltog_row, &ltog_col);
43053b04fa6SJed Brown     CHKERRQ(ierr);
43153b04fa6SJed Brown     ierr = ISLocalToGlobalMappingApply(ltog_row, num_entries, rows, rows);
43253b04fa6SJed Brown     CHKERRQ(ierr);
43353b04fa6SJed Brown     ierr = ISLocalToGlobalMappingApply(ltog_col, num_entries, cols, cols);
43453b04fa6SJed Brown     CHKERRQ(ierr);
43553b04fa6SJed Brown     ierr = MatSetPreallocationCOO(mat_coarse, num_entries, rows, cols);
43653b04fa6SJed Brown     CHKERRQ(ierr);
43753b04fa6SJed Brown     free(rows);
43853b04fa6SJed Brown     free(cols);
43953b04fa6SJed Brown     CeedVectorCreate(ceed, num_entries, &coo_values);
44005b9c820SJed Brown     ierr = PetscLogEventBegin(assemble_event, mat_coarse, 0, 0, 0); CHKERRQ(ierr);
44153b04fa6SJed Brown     CeedOperatorLinearAssemble(user_O[0]->op, coo_values);
44253b04fa6SJed Brown     const CeedScalar *values;
44353b04fa6SJed Brown     CeedVectorGetArrayRead(coo_values, CEED_MEM_HOST, &values);
44453b04fa6SJed Brown     ierr = MatSetValuesCOO(mat_coarse, values, ADD_VALUES); CHKERRQ(ierr);
44553b04fa6SJed Brown     CeedVectorRestoreArrayRead(coo_values, &values);
44605b9c820SJed Brown     ierr = PetscLogEventEnd(assemble_event, mat_coarse, 0, 0, 0); CHKERRQ(ierr);
44753b04fa6SJed Brown     CeedVectorDestroy(&coo_values);
44853b04fa6SJed Brown   }
44915ce0ef0Sjeremylt 
4506c5df90dSjeremylt   // Set up KSP
4516c5df90dSjeremylt   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
4526c5df90dSjeremylt   {
4536c5df90dSjeremylt     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
4546c5df90dSjeremylt     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
4556c5df90dSjeremylt     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
4566c5df90dSjeremylt                             PETSC_DEFAULT); CHKERRQ(ierr);
4576c5df90dSjeremylt   }
4586c5df90dSjeremylt   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
4599b072555Sjeremylt   ierr = KSPSetOperators(ksp, mat_O[fine_level], mat_O[fine_level]);
4606c5df90dSjeremylt   CHKERRQ(ierr);
4616c5df90dSjeremylt 
4626c5df90dSjeremylt   // Set up PCMG
4636c5df90dSjeremylt   ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
4649b072555Sjeremylt   PCMGCycleType pcmg_cycle_type = PC_MG_CYCLE_V;
4656c5df90dSjeremylt   {
4666c5df90dSjeremylt     ierr = PCSetType(pc, PCMG); CHKERRQ(ierr);
4676c5df90dSjeremylt 
4686c5df90dSjeremylt     // PCMG levels
4699b072555Sjeremylt     ierr = PCMGSetLevels(pc, num_levels, NULL); CHKERRQ(ierr);
4709b072555Sjeremylt     for (int i=0; i<num_levels; i++) {
4716c5df90dSjeremylt       // Smoother
4726c5df90dSjeremylt       KSP smoother;
4736c5df90dSjeremylt       PC smoother_pc;
4746c5df90dSjeremylt       ierr = PCMGGetSmoother(pc, i, &smoother); CHKERRQ(ierr);
4756c5df90dSjeremylt       ierr = KSPSetType(smoother, KSPCHEBYSHEV); CHKERRQ(ierr);
4766c5df90dSjeremylt       ierr = KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1); CHKERRQ(ierr);
4776c5df90dSjeremylt       ierr = KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE); CHKERRQ(ierr);
4789b072555Sjeremylt       ierr = KSPSetOperators(smoother, mat_O[i], mat_O[i]); CHKERRQ(ierr);
4796c5df90dSjeremylt       ierr = KSPGetPC(smoother, &smoother_pc); CHKERRQ(ierr);
4806c5df90dSjeremylt       ierr = PCSetType(smoother_pc, PCJACOBI); CHKERRQ(ierr);
4816c5df90dSjeremylt       ierr = PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
4826c5df90dSjeremylt 
4836c5df90dSjeremylt       // Work vector
4849b072555Sjeremylt       if (i < num_levels - 1) {
4856c5df90dSjeremylt         ierr = PCMGSetX(pc, i, X[i]); CHKERRQ(ierr);
4866c5df90dSjeremylt       }
4876c5df90dSjeremylt 
4886c5df90dSjeremylt       // Level transfers
4896c5df90dSjeremylt       if (i > 0) {
4906c5df90dSjeremylt         // Interpolation
4919b072555Sjeremylt         ierr = PCMGSetInterpolation(pc, i, mat_pr[i]); CHKERRQ(ierr);
4926c5df90dSjeremylt       }
4936c5df90dSjeremylt 
4946c5df90dSjeremylt       // Coarse solve
4956c5df90dSjeremylt       KSP coarse;
4966c5df90dSjeremylt       PC coarse_pc;
4976c5df90dSjeremylt       ierr = PCMGGetCoarseSolve(pc, &coarse); CHKERRQ(ierr);
49815ce0ef0Sjeremylt       ierr = KSPSetType(coarse, KSPPREONLY); CHKERRQ(ierr);
4999b072555Sjeremylt       ierr = KSPSetOperators(coarse, mat_coarse, mat_coarse); CHKERRQ(ierr);
50015ce0ef0Sjeremylt 
5016c5df90dSjeremylt       ierr = KSPGetPC(coarse, &coarse_pc); CHKERRQ(ierr);
50215ce0ef0Sjeremylt       ierr = PCSetType(coarse_pc, PCGAMG); CHKERRQ(ierr);
50315ce0ef0Sjeremylt 
50415ce0ef0Sjeremylt       ierr = KSPSetOptionsPrefix(coarse, "coarse_"); CHKERRQ(ierr);
50515ce0ef0Sjeremylt       ierr = PCSetOptionsPrefix(coarse_pc, "coarse_"); CHKERRQ(ierr);
50615ce0ef0Sjeremylt       ierr = KSPSetFromOptions(coarse); CHKERRQ(ierr);
50715ce0ef0Sjeremylt       ierr = PCSetFromOptions(coarse_pc); CHKERRQ(ierr);
5086c5df90dSjeremylt     }
5096c5df90dSjeremylt 
5106c5df90dSjeremylt     // PCMG options
5116c5df90dSjeremylt     ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr);
5126c5df90dSjeremylt     ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr);
5139b072555Sjeremylt     ierr = PCMGSetCycleType(pc, pcmg_cycle_type); CHKERRQ(ierr);
5146c5df90dSjeremylt   }
5156c5df90dSjeremylt 
5166c5df90dSjeremylt   // First run, if benchmarking
5176c5df90dSjeremylt   if (benchmark_mode) {
5186c5df90dSjeremylt     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
5196c5df90dSjeremylt     CHKERRQ(ierr);
5209b072555Sjeremylt     ierr = VecZeroEntries(X[fine_level]); CHKERRQ(ierr);
5216c5df90dSjeremylt     my_rt_start = MPI_Wtime();
5229b072555Sjeremylt     ierr = KSPSolve(ksp, rhs, X[fine_level]); CHKERRQ(ierr);
5236c5df90dSjeremylt     my_rt = MPI_Wtime() - my_rt_start;
5246c5df90dSjeremylt     ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
5256c5df90dSjeremylt     CHKERRQ(ierr);
5266c5df90dSjeremylt     // Set maxits based on first iteration timing
5276c5df90dSjeremylt     if (my_rt > 0.02) {
5286c5df90dSjeremylt       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
5296c5df90dSjeremylt       CHKERRQ(ierr);
5306c5df90dSjeremylt     } else {
5316c5df90dSjeremylt       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
5326c5df90dSjeremylt       CHKERRQ(ierr);
5336c5df90dSjeremylt     }
5346c5df90dSjeremylt   }
5356c5df90dSjeremylt 
5366c5df90dSjeremylt   // Timed solve
5379b072555Sjeremylt   ierr = VecZeroEntries(X[fine_level]); CHKERRQ(ierr);
5386c5df90dSjeremylt   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
53909a940d7Sjeremylt 
54009a940d7Sjeremylt   // -- Performance logging
5419b072555Sjeremylt   ierr = PetscLogStageRegister("Solve Stage", &solve_stage); CHKERRQ(ierr);
5429b072555Sjeremylt   ierr = PetscLogStagePush(solve_stage); CHKERRQ(ierr);
54309a940d7Sjeremylt 
54409a940d7Sjeremylt   // -- Solve
5456c5df90dSjeremylt   my_rt_start = MPI_Wtime();
5469b072555Sjeremylt   ierr = KSPSolve(ksp, rhs, X[fine_level]); CHKERRQ(ierr);
5476c5df90dSjeremylt   my_rt = MPI_Wtime() - my_rt_start;
5486c5df90dSjeremylt 
54909a940d7Sjeremylt 
55009a940d7Sjeremylt   // -- Performance logging
55109a940d7Sjeremylt   ierr = PetscLogStagePop();
55209a940d7Sjeremylt 
5536c5df90dSjeremylt   // Output results
5546c5df90dSjeremylt   {
5559b072555Sjeremylt     KSPType ksp_type;
5569b072555Sjeremylt     PCMGType pcmg_type;
5576c5df90dSjeremylt     KSPConvergedReason reason;
5586c5df90dSjeremylt     PetscReal rnorm;
5596c5df90dSjeremylt     PetscInt its;
5609b072555Sjeremylt     ierr = KSPGetType(ksp, &ksp_type); CHKERRQ(ierr);
5616c5df90dSjeremylt     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
5626c5df90dSjeremylt     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
5636c5df90dSjeremylt     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
5649b072555Sjeremylt     ierr = PCMGGetType(pc, &pcmg_type); CHKERRQ(ierr);
5656c5df90dSjeremylt     if (!test_mode || reason < 0 || rnorm > 1e-8) {
5666c5df90dSjeremylt       ierr = PetscPrintf(comm,
5676c5df90dSjeremylt                          "  KSP:\n"
5686c5df90dSjeremylt                          "    KSP Type                                : %s\n"
5696c5df90dSjeremylt                          "    KSP Convergence                         : %s\n"
570*a9b2c5ddSrezgarshakeri                          "    Total KSP Iterations                    : %" PetscInt_FMT "\n"
5716c5df90dSjeremylt                          "    Final rnorm                             : %e\n",
5729b072555Sjeremylt                          ksp_type, KSPConvergedReasons[reason], its,
5736c5df90dSjeremylt                          (double)rnorm); CHKERRQ(ierr);
5746c5df90dSjeremylt       ierr = PetscPrintf(comm,
5756c5df90dSjeremylt                          "  PCMG:\n"
5766c5df90dSjeremylt                          "    PCMG Type                               : %s\n"
5776c5df90dSjeremylt                          "    PCMG Cycle Type                         : %s\n",
5789b072555Sjeremylt                          PCMGTypes[pcmg_type],
5799b072555Sjeremylt                          PCMGCycleTypes[pcmg_cycle_type]); CHKERRQ(ierr);
5806c5df90dSjeremylt     }
5816c5df90dSjeremylt     if (!test_mode) {
5826c5df90dSjeremylt       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
5836c5df90dSjeremylt     }
5846c5df90dSjeremylt     {
5859b072555Sjeremylt       PetscReal max_error;
5869b072555Sjeremylt       ierr = ComputeErrorMax(user_O[fine_level], op_error, X[fine_level], target,
5879b072555Sjeremylt                              &max_error); CHKERRQ(ierr);
5886c5df90dSjeremylt       PetscReal tol = 5e-2;
5899b072555Sjeremylt       if (!test_mode || max_error > tol) {
5906c5df90dSjeremylt         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
5916c5df90dSjeremylt         CHKERRQ(ierr);
5926c5df90dSjeremylt         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
5936c5df90dSjeremylt         CHKERRQ(ierr);
5946c5df90dSjeremylt         ierr = PetscPrintf(comm,
5956c5df90dSjeremylt                            "    Pointwise Error (max)                   : %e\n"
5966c5df90dSjeremylt                            "    CG Solve Time                           : %g (%g) sec\n",
5979b072555Sjeremylt                            (double)max_error, rt_max, rt_min); CHKERRQ(ierr);
5986c5df90dSjeremylt       }
5996c5df90dSjeremylt     }
6006c5df90dSjeremylt     if (benchmark_mode && (!test_mode)) {
6016c5df90dSjeremylt       ierr = PetscPrintf(comm,
6026c5df90dSjeremylt                          "    DoFs/Sec in CG                     : %g (%g) million\n",
6039b072555Sjeremylt                          1e-6*g_size[fine_level]*its/rt_max,
6049b072555Sjeremylt                          1e-6*g_size[fine_level]*its/rt_min);
6056c5df90dSjeremylt       CHKERRQ(ierr);
6066c5df90dSjeremylt     }
6076c5df90dSjeremylt   }
6086c5df90dSjeremylt 
6096c5df90dSjeremylt   if (write_solution) {
6109b072555Sjeremylt     PetscViewer vtk_viewer_soln;
6116c5df90dSjeremylt 
6129b072555Sjeremylt     ierr = PetscViewerCreate(comm, &vtk_viewer_soln); CHKERRQ(ierr);
6139b072555Sjeremylt     ierr = PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK); CHKERRQ(ierr);
6149b072555Sjeremylt     ierr = PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"); CHKERRQ(ierr);
6159b072555Sjeremylt     ierr = VecView(X[fine_level], vtk_viewer_soln); CHKERRQ(ierr);
6169b072555Sjeremylt     ierr = PetscViewerDestroy(&vtk_viewer_soln); CHKERRQ(ierr);
6176c5df90dSjeremylt   }
6186c5df90dSjeremylt 
6196c5df90dSjeremylt   // Cleanup
6209b072555Sjeremylt   for (int i=0; i<num_levels; i++) {
6216c5df90dSjeremylt     ierr = VecDestroy(&X[i]); CHKERRQ(ierr);
6229b072555Sjeremylt     ierr = VecDestroy(&X_loc[i]); CHKERRQ(ierr);
6236c5df90dSjeremylt     ierr = VecDestroy(&mult[i]); CHKERRQ(ierr);
6249b072555Sjeremylt     ierr = VecDestroy(&user_O[i]->Y_loc); CHKERRQ(ierr);
6259b072555Sjeremylt     ierr = MatDestroy(&mat_O[i]); CHKERRQ(ierr);
6269b072555Sjeremylt     ierr = PetscFree(user_O[i]); CHKERRQ(ierr);
6276c5df90dSjeremylt     if (i > 0) {
6289b072555Sjeremylt       ierr = MatDestroy(&mat_pr[i]); CHKERRQ(ierr);
6299b072555Sjeremylt       ierr = PetscFree(user_pr[i]); CHKERRQ(ierr);
6306c5df90dSjeremylt     }
6319b072555Sjeremylt     ierr = CeedDataDestroy(i, ceed_data[i]); CHKERRQ(ierr);
6326c5df90dSjeremylt     ierr = DMDestroy(&dm[i]); CHKERRQ(ierr);
6336c5df90dSjeremylt   }
6349b072555Sjeremylt   ierr = PetscFree(level_degrees); CHKERRQ(ierr);
6356c5df90dSjeremylt   ierr = PetscFree(dm); CHKERRQ(ierr);
6366c5df90dSjeremylt   ierr = PetscFree(X); CHKERRQ(ierr);
6379b072555Sjeremylt   ierr = PetscFree(X_loc); CHKERRQ(ierr);
6386c5df90dSjeremylt   ierr = PetscFree(mult); CHKERRQ(ierr);
6399b072555Sjeremylt   ierr = PetscFree(mat_O); CHKERRQ(ierr);
6409b072555Sjeremylt   ierr = PetscFree(mat_pr); CHKERRQ(ierr);
6419b072555Sjeremylt   ierr = PetscFree(ceed_data); CHKERRQ(ierr);
6429b072555Sjeremylt   ierr = PetscFree(user_O); CHKERRQ(ierr);
6439b072555Sjeremylt   ierr = PetscFree(user_pr); CHKERRQ(ierr);
6449b072555Sjeremylt   ierr = PetscFree(l_size); CHKERRQ(ierr);
6459b072555Sjeremylt   ierr = PetscFree(xl_size); CHKERRQ(ierr);
6469b072555Sjeremylt   ierr = PetscFree(g_size); CHKERRQ(ierr);
6476c5df90dSjeremylt   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
6489b072555Sjeremylt   ierr = VecDestroy(&rhs_loc); CHKERRQ(ierr);
6499b072555Sjeremylt   ierr = MatDestroy(&mat_coarse); CHKERRQ(ierr);
6506c5df90dSjeremylt   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
6519b072555Sjeremylt   ierr = DMDestroy(&dm_orig); CHKERRQ(ierr);
6526c5df90dSjeremylt   CeedVectorDestroy(&target);
6539b072555Sjeremylt   CeedQFunctionDestroy(&qf_error);
6549b072555Sjeremylt   CeedOperatorDestroy(&op_error);
6556c5df90dSjeremylt   CeedDestroy(&ceed);
6566c5df90dSjeremylt   return PetscFinalize();
6576c5df90dSjeremylt }
658