xref: /libCEED/examples/petsc/multigrid.c (revision 9b072555b57804a6f4e0fc2b1ad83be89838f0e5)
16c5df90dSjeremylt // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
26c5df90dSjeremylt // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
36c5df90dSjeremylt // reserved. See files LICENSE and NOTICE for details.
46c5df90dSjeremylt //
56c5df90dSjeremylt // This file is part of CEED, a collection of benchmarks, miniapps, software
66c5df90dSjeremylt // libraries and APIs for efficient high-order finite element and spectral
76c5df90dSjeremylt // element discretizations for exascale applications. For more information and
86c5df90dSjeremylt // source code availability see http://github.com/ceed.
96c5df90dSjeremylt //
106c5df90dSjeremylt // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
116c5df90dSjeremylt // a collaborative effort of two U.S. Department of Energy organizations (Office
126c5df90dSjeremylt // of Science and the National Nuclear Security Administration) responsible for
136c5df90dSjeremylt // the planning and preparation of a capable exascale ecosystem, including
146c5df90dSjeremylt // software, applications, hardware, advanced system engineering and early
156c5df90dSjeremylt // testbed platforms, in support of the nation's exascale computing imperative.
166c5df90dSjeremylt 
176c5df90dSjeremylt //                        libCEED + PETSc Example: CEED BPs 3-6 with Multigrid
186c5df90dSjeremylt //
196c5df90dSjeremylt // This example demonstrates a simple usage of libCEED with PETSc to solve the
206c5df90dSjeremylt // CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps.
216c5df90dSjeremylt //
226c5df90dSjeremylt // The code uses higher level communication protocols in DMPlex.
236c5df90dSjeremylt //
246c5df90dSjeremylt // Build with:
256c5df90dSjeremylt //
266c5df90dSjeremylt //     make multigrid [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
276c5df90dSjeremylt //
286c5df90dSjeremylt // Sample runs:
296c5df90dSjeremylt //
306c5df90dSjeremylt //     multigrid -problem bp3
3128688798Sjeremylt //     multigrid -problem bp4
3228688798Sjeremylt //     multigrid -problem bp5 -ceed /cpu/self
336c5df90dSjeremylt //     multigrid -problem bp6 -ceed /gpu/cuda
346c5df90dSjeremylt //
356c5df90dSjeremylt //TESTARGS -ceed {ceed_resource} -test -problem bp3 -degree 3
366c5df90dSjeremylt 
376c5df90dSjeremylt /// @file
386c5df90dSjeremylt /// CEED BPs 1-6 multigrid example using PETSc
396c5df90dSjeremylt const char help[] = "Solve CEED BPs using p-multigrid with PETSc and DMPlex\n";
406c5df90dSjeremylt 
41e83e87a5Sjeremylt #include "bps.h"
426c5df90dSjeremylt 
436c5df90dSjeremylt int main(int argc, char **argv) {
446c5df90dSjeremylt   PetscInt ierr;
456c5df90dSjeremylt   MPI_Comm comm;
46cb0b5415Sjeremylt   char filename[PETSC_MAX_PATH_LEN],
47*9b072555Sjeremylt        ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self";
486c5df90dSjeremylt   double my_rt_start, my_rt, rt_min, rt_max;
49*9b072555Sjeremylt   PetscInt degree = 3, q_extra, *l_size, *xl_size, *g_size, dim = 3, fine_level,
50*9b072555Sjeremylt            mesh_elem[3] = {3, 3, 3}, num_comp_u = 1, num_levels = degree, *level_degrees;
516c5df90dSjeremylt   PetscScalar *r;
52cfa59c5bSRey   PetscScalar eps = 1.0;
536c5df90dSjeremylt   PetscBool test_mode, benchmark_mode, read_mesh, write_solution;
54*9b072555Sjeremylt   PetscLogStage solve_stage;
55*9b072555Sjeremylt   DM  *dm, dm_orig;
56*9b072555Sjeremylt   SNES snes_dummy;
576c5df90dSjeremylt   KSP ksp;
586c5df90dSjeremylt   PC pc;
59*9b072555Sjeremylt   Mat *mat_O, *mat_pr, mat_coarse;
60*9b072555Sjeremylt   Vec *X, *X_loc, *mult, rhs, rhs_loc;
61*9b072555Sjeremylt   PetscMemType mem_type;
62*9b072555Sjeremylt   UserO *user_O;
63*9b072555Sjeremylt   UserProlongRestr *user_pr;
646c5df90dSjeremylt   Ceed ceed;
65*9b072555Sjeremylt   CeedData *ceed_data;
66*9b072555Sjeremylt   CeedVector rhs_ceed, target;
67*9b072555Sjeremylt   CeedQFunction qf_error, qf_restrict, qf_prolong;
68*9b072555Sjeremylt   CeedOperator op_error;
69*9b072555Sjeremylt   BPType bp_choice;
70*9b072555Sjeremylt   CoarsenType coarsen;
716c5df90dSjeremylt 
726c5df90dSjeremylt   ierr = PetscInitialize(&argc, &argv, NULL, help);
736c5df90dSjeremylt   if (ierr) return ierr;
746c5df90dSjeremylt   comm = PETSC_COMM_WORLD;
756c5df90dSjeremylt 
766c5df90dSjeremylt   // Parse command line options
776c5df90dSjeremylt   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr);
78*9b072555Sjeremylt   bp_choice = CEED_BP3;
796c5df90dSjeremylt   ierr = PetscOptionsEnum("-problem",
806c5df90dSjeremylt                           "CEED benchmark problem to solve", NULL,
81*9b072555Sjeremylt                           bp_types, (PetscEnum)bp_choice, (PetscEnum *)&bp_choice,
826c5df90dSjeremylt                           NULL); CHKERRQ(ierr);
83*9b072555Sjeremylt   num_comp_u = bp_options[bp_choice].num_comp_u;
846c5df90dSjeremylt   test_mode = PETSC_FALSE;
856c5df90dSjeremylt   ierr = PetscOptionsBool("-test",
866c5df90dSjeremylt                           "Testing mode (do not print unless error is large)",
876c5df90dSjeremylt                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
886c5df90dSjeremylt   benchmark_mode = PETSC_FALSE;
896c5df90dSjeremylt   ierr = PetscOptionsBool("-benchmark",
906c5df90dSjeremylt                           "Benchmarking mode (prints benchmark statistics)",
916c5df90dSjeremylt                           NULL, benchmark_mode, &benchmark_mode, NULL);
926c5df90dSjeremylt   CHKERRQ(ierr);
936c5df90dSjeremylt   write_solution = PETSC_FALSE;
946c5df90dSjeremylt   ierr = PetscOptionsBool("-write_solution",
956c5df90dSjeremylt                           "Write solution for visualization",
966c5df90dSjeremylt                           NULL, write_solution, &write_solution, NULL);
976c5df90dSjeremylt   CHKERRQ(ierr);
98cfa59c5bSRey   ierr = PetscOptionsScalar("-eps",
99cfa59c5bSRey                             "Epsilon parameter for Kershaw mesh transformation",
100cfa59c5bSRey                             NULL, eps, &eps, NULL);
101cfa59c5bSRey   if (eps > 1 || eps <= 0) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
102cfa59c5bSRey                                       "-eps %D must be (0,1]", eps);
1036c5df90dSjeremylt   degree = test_mode ? 3 : 2;
1046c5df90dSjeremylt   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
1056c5df90dSjeremylt                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
1066c5df90dSjeremylt   if (degree < 1) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
1076c5df90dSjeremylt                              "-degree %D must be at least 1", degree);
108*9b072555Sjeremylt   q_extra = bp_options[bp_choice].q_extra;
109*9b072555Sjeremylt   ierr = PetscOptionsInt("-q_extra", "Number of extra quadrature points",
110*9b072555Sjeremylt                          NULL, q_extra, &q_extra, NULL); CHKERRQ(ierr);
1116c5df90dSjeremylt   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
112*9b072555Sjeremylt                             NULL, ceed_resource, ceed_resource,
113*9b072555Sjeremylt                             sizeof(ceed_resource), NULL); CHKERRQ(ierr);
1146c5df90dSjeremylt   coarsen = COARSEN_UNIFORM;
1156c5df90dSjeremylt   ierr = PetscOptionsEnum("-coarsen",
1166c5df90dSjeremylt                           "Coarsening strategy to use", NULL,
117*9b072555Sjeremylt                           coarsen_types, (PetscEnum)coarsen,
1186c5df90dSjeremylt                           (PetscEnum *)&coarsen, NULL); CHKERRQ(ierr);
119cb32e2e7SValeria Barra   read_mesh = PETSC_FALSE;
1206c5df90dSjeremylt   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL,
1216c5df90dSjeremylt                             filename, filename, sizeof(filename), &read_mesh);
1226c5df90dSjeremylt   CHKERRQ(ierr);
1236c5df90dSjeremylt   if (!read_mesh) {
1246c5df90dSjeremylt     PetscInt tmp = dim;
1256c5df90dSjeremylt     ierr = PetscOptionsIntArray("-cells","Number of cells per dimension", NULL,
126*9b072555Sjeremylt                                 mesh_elem, &tmp, NULL); CHKERRQ(ierr);
1276c5df90dSjeremylt   }
1286c5df90dSjeremylt   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
1296c5df90dSjeremylt 
1309396343dSjeremylt   // Set up libCEED
131*9b072555Sjeremylt   CeedInit(ceed_resource, &ceed);
132*9b072555Sjeremylt   CeedMemType mem_type_backend;
133*9b072555Sjeremylt   CeedGetPreferredMemType(ceed, &mem_type_backend);
1349396343dSjeremylt 
1356c5df90dSjeremylt   // Setup DM
1366c5df90dSjeremylt   if (read_mesh) {
137*9b072555Sjeremylt     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dm_orig);
1386c5df90dSjeremylt     CHKERRQ(ierr);
1396c5df90dSjeremylt   } else {
140*9b072555Sjeremylt     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_FALSE, mesh_elem, NULL,
141*9b072555Sjeremylt                                NULL, NULL, PETSC_TRUE, &dm_orig); CHKERRQ(ierr);
1426c5df90dSjeremylt   }
1436c5df90dSjeremylt 
1446c5df90dSjeremylt   {
145*9b072555Sjeremylt     DM dm_dist = NULL;
1466c5df90dSjeremylt     PetscPartitioner part;
1476c5df90dSjeremylt 
148*9b072555Sjeremylt     ierr = DMPlexGetPartitioner(dm_orig, &part); CHKERRQ(ierr);
1496c5df90dSjeremylt     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
150*9b072555Sjeremylt     ierr = DMPlexDistribute(dm_orig, 0, NULL, &dm_dist); CHKERRQ(ierr);
151*9b072555Sjeremylt     if (dm_dist) {
152*9b072555Sjeremylt       ierr = DMDestroy(&dm_orig); CHKERRQ(ierr);
153*9b072555Sjeremylt       dm_orig = dm_dist;
1546c5df90dSjeremylt     }
1556c5df90dSjeremylt   }
1566c5df90dSjeremylt 
157*9b072555Sjeremylt   // Apply Kershaw mesh transformation
158*9b072555Sjeremylt   ierr = Kershaw(dm_orig, eps); CHKERRQ(ierr);
159cfa59c5bSRey 
160*9b072555Sjeremylt   VecType vec_type;
161*9b072555Sjeremylt   switch (mem_type_backend) {
162*9b072555Sjeremylt   case CEED_MEM_HOST: vec_type = VECSTANDARD; break;
163b68a8d79SJed Brown   case CEED_MEM_DEVICE: {
164b68a8d79SJed Brown     const char *resolved;
165b68a8d79SJed Brown     CeedGetResource(ceed, &resolved);
166*9b072555Sjeremylt     if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA;
167ca2d516cSJed Brown     else if (strstr(resolved, "/gpu/hip/occa"))
168*9b072555Sjeremylt       vec_type = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678
169*9b072555Sjeremylt     else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP;
170*9b072555Sjeremylt     else vec_type = VECSTANDARD;
171b68a8d79SJed Brown   }
172b68a8d79SJed Brown   }
173*9b072555Sjeremylt   ierr = DMSetVecType(dm_orig, vec_type); CHKERRQ(ierr);
174*9b072555Sjeremylt   ierr = DMSetFromOptions(dm_orig); CHKERRQ(ierr);
175b68a8d79SJed Brown 
1766c5df90dSjeremylt   // Allocate arrays for PETSc objects for each level
1776c5df90dSjeremylt   switch (coarsen) {
1786c5df90dSjeremylt   case COARSEN_UNIFORM:
179*9b072555Sjeremylt     num_levels = degree;
1806c5df90dSjeremylt     break;
181dc7d240cSValeria Barra   case COARSEN_LOGARITHMIC:
182*9b072555Sjeremylt     num_levels = ceil(log(degree)/log(2)) + 1;
1836c5df90dSjeremylt     break;
1846c5df90dSjeremylt   }
185*9b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &level_degrees); CHKERRQ(ierr);
186*9b072555Sjeremylt   fine_level = num_levels - 1;
18761608365Sjeremylt 
1886c5df90dSjeremylt   switch (coarsen) {
1896c5df90dSjeremylt   case COARSEN_UNIFORM:
190*9b072555Sjeremylt     for (int i=0; i<num_levels; i++) level_degrees[i] = i + 1;
1916c5df90dSjeremylt     break;
192dc7d240cSValeria Barra   case COARSEN_LOGARITHMIC:
193*9b072555Sjeremylt     for (int i=0; i<num_levels - 1; i++) level_degrees[i] = pow(2,i);
194*9b072555Sjeremylt     level_degrees[fine_level] = degree;
1956c5df90dSjeremylt     break;
1966c5df90dSjeremylt   }
197*9b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &dm); CHKERRQ(ierr);
198*9b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &X); CHKERRQ(ierr);
199*9b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &X_loc); CHKERRQ(ierr);
200*9b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &mult); CHKERRQ(ierr);
201*9b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &user_O); CHKERRQ(ierr);
202*9b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &user_pr); CHKERRQ(ierr);
203*9b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &mat_O); CHKERRQ(ierr);
204*9b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &mat_pr); CHKERRQ(ierr);
205*9b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &l_size); CHKERRQ(ierr);
206*9b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &xl_size); CHKERRQ(ierr);
207*9b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &g_size); CHKERRQ(ierr);
2086c5df90dSjeremylt 
2096c5df90dSjeremylt   // Setup DM and Operator Mat Shells for each level
210*9b072555Sjeremylt   for (CeedInt i=0; i<num_levels; i++) {
2116c5df90dSjeremylt     // Create DM
212*9b072555Sjeremylt     ierr = DMClone(dm_orig, &dm[i]); CHKERRQ(ierr);
213*9b072555Sjeremylt     ierr = DMGetVecType(dm_orig, &vec_type); CHKERRQ(ierr);
214*9b072555Sjeremylt     ierr = DMSetVecType(dm[i], vec_type); CHKERRQ(ierr);
215e83e87a5Sjeremylt     PetscInt dim;
216e83e87a5Sjeremylt     ierr = DMGetDimension(dm[i], &dim); CHKERRQ(ierr);
217*9b072555Sjeremylt     ierr = SetupDMByDegree(dm[i], level_degrees[i], num_comp_u, dim,
218*9b072555Sjeremylt                            bp_options[bp_choice].enforce_bc, bp_options[bp_choice].bc_func);
2196c5df90dSjeremylt     CHKERRQ(ierr);
2206c5df90dSjeremylt 
2216c5df90dSjeremylt     // Create vectors
2226c5df90dSjeremylt     ierr = DMCreateGlobalVector(dm[i], &X[i]); CHKERRQ(ierr);
223*9b072555Sjeremylt     ierr = VecGetLocalSize(X[i], &l_size[i]); CHKERRQ(ierr);
224*9b072555Sjeremylt     ierr = VecGetSize(X[i], &g_size[i]); CHKERRQ(ierr);
225*9b072555Sjeremylt     ierr = DMCreateLocalVector(dm[i], &X_loc[i]); CHKERRQ(ierr);
226*9b072555Sjeremylt     ierr = VecGetSize(X_loc[i], &xl_size[i]); CHKERRQ(ierr);
2276c5df90dSjeremylt 
2286c5df90dSjeremylt     // Operator
229*9b072555Sjeremylt     ierr = PetscMalloc1(1, &user_O[i]); CHKERRQ(ierr);
230*9b072555Sjeremylt     ierr = MatCreateShell(comm, l_size[i], l_size[i], g_size[i], g_size[i],
231*9b072555Sjeremylt                           user_O[i], &mat_O[i]); CHKERRQ(ierr);
232*9b072555Sjeremylt     ierr = MatShellSetOperation(mat_O[i], MATOP_MULT,
233ce74dcefSjeremylt                                 (void(*)(void))MatMult_Ceed); CHKERRQ(ierr);
234*9b072555Sjeremylt     ierr = MatShellSetOperation(mat_O[i], MATOP_GET_DIAGONAL,
235ce74dcefSjeremylt                                 (void(*)(void))MatGetDiag); CHKERRQ(ierr);
236*9b072555Sjeremylt     ierr = MatShellSetVecType(mat_O[i], vec_type); CHKERRQ(ierr);
2376c5df90dSjeremylt 
2386c5df90dSjeremylt     // Level transfers
2396c5df90dSjeremylt     if (i > 0) {
2406c5df90dSjeremylt       // Interp
241*9b072555Sjeremylt       ierr = PetscMalloc1(1, &user_pr[i]); CHKERRQ(ierr);
242*9b072555Sjeremylt       ierr = MatCreateShell(comm, l_size[i], l_size[i-1], g_size[i], g_size[i-1],
243*9b072555Sjeremylt                             user_pr[i], &mat_pr[i]); CHKERRQ(ierr);
244*9b072555Sjeremylt       ierr = MatShellSetOperation(mat_pr[i], MATOP_MULT,
245a97643b0Sjeremylt                                   (void(*)(void))MatMult_Prolong);
2466c5df90dSjeremylt       CHKERRQ(ierr);
247*9b072555Sjeremylt       ierr = MatShellSetOperation(mat_pr[i], MATOP_MULT_TRANSPOSE,
2486c5df90dSjeremylt                                   (void(*)(void))MatMult_Restrict);
2496c5df90dSjeremylt       CHKERRQ(ierr);
250*9b072555Sjeremylt       ierr = MatShellSetVecType(mat_pr[i], vec_type); CHKERRQ(ierr);
2516c5df90dSjeremylt     }
2526c5df90dSjeremylt   }
253*9b072555Sjeremylt   ierr = VecDuplicate(X[fine_level], &rhs); CHKERRQ(ierr);
2546c5df90dSjeremylt 
2556c5df90dSjeremylt   // Print global grid information
2566c5df90dSjeremylt   if (!test_mode) {
257*9b072555Sjeremylt     PetscInt P = degree + 1, Q = P + q_extra;
2589396343dSjeremylt 
259*9b072555Sjeremylt     const char *used_resource;
260*9b072555Sjeremylt     CeedGetResource(ceed, &used_resource);
2619396343dSjeremylt 
262*9b072555Sjeremylt     ierr = VecGetType(X[0], &vec_type); CHKERRQ(ierr);
2639396343dSjeremylt 
2646c5df90dSjeremylt     ierr = PetscPrintf(comm,
2656c5df90dSjeremylt                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc + PCMG --\n"
2669396343dSjeremylt                        "  PETSc:\n"
2679396343dSjeremylt                        "    PETSc Vec Type                     : %s\n"
2686c5df90dSjeremylt                        "  libCEED:\n"
2696c5df90dSjeremylt                        "    libCEED Backend                    : %s\n"
2709396343dSjeremylt                        "    libCEED Backend MemType            : %s\n"
2716c5df90dSjeremylt                        "  Mesh:\n"
2726c5df90dSjeremylt                        "    Number of 1D Basis Nodes (p)       : %d\n"
2736c5df90dSjeremylt                        "    Number of 1D Quadrature Points (q) : %d\n"
2746c5df90dSjeremylt                        "    Global Nodes                       : %D\n"
2756c5df90dSjeremylt                        "    Owned Nodes                        : %D\n"
276db419314Sjeremylt                        "    DoF per node                       : %D\n"
2776c5df90dSjeremylt                        "  Multigrid:\n"
2786c5df90dSjeremylt                        "    Number of Levels                   : %d\n",
279*9b072555Sjeremylt                        bp_choice+1, vec_type, used_resource,
280*9b072555Sjeremylt                        CeedMemTypes[mem_type_backend],
281*9b072555Sjeremylt                        P, Q, g_size[fine_level]/num_comp_u, l_size[fine_level]/num_comp_u,
282*9b072555Sjeremylt                        num_comp_u, num_levels); CHKERRQ(ierr);
2836c5df90dSjeremylt   }
2846c5df90dSjeremylt 
2856c5df90dSjeremylt   // Create RHS vector
286*9b072555Sjeremylt   ierr = VecDuplicate(X_loc[fine_level], &rhs_loc); CHKERRQ(ierr);
287*9b072555Sjeremylt   ierr = VecZeroEntries(rhs_loc); CHKERRQ(ierr);
288*9b072555Sjeremylt   ierr = VecGetArrayAndMemType(rhs_loc, &r, &mem_type); CHKERRQ(ierr);
289*9b072555Sjeremylt   CeedVectorCreate(ceed, xl_size[fine_level], &rhs_ceed);
290*9b072555Sjeremylt   CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r);
2916c5df90dSjeremylt 
2926c5df90dSjeremylt   // Set up libCEED operators on each level
293*9b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &ceed_data); CHKERRQ(ierr);
294*9b072555Sjeremylt   for (int i=0; i<num_levels; i++) {
2956c5df90dSjeremylt     // Print level information
296*9b072555Sjeremylt     if (!test_mode && (i == 0 || i == fine_level)) {
2976c5df90dSjeremylt       ierr = PetscPrintf(comm,"    Level %D (%s):\n"
2986c5df90dSjeremylt                          "      Number of 1D Basis Nodes (p)     : %d\n"
2996c5df90dSjeremylt                          "      Global Nodes                     : %D\n"
3006c5df90dSjeremylt                          "      Owned Nodes                      : %D\n",
301*9b072555Sjeremylt                          i, (i? "fine" : "coarse"), level_degrees[i] + 1,
302*9b072555Sjeremylt                          g_size[i]/num_comp_u, l_size[i]/num_comp_u); CHKERRQ(ierr);
3036c5df90dSjeremylt     }
304*9b072555Sjeremylt     ierr = PetscMalloc1(1, &ceed_data[i]); CHKERRQ(ierr);
305*9b072555Sjeremylt     ierr = SetupLibceedByDegree(dm[i], ceed, level_degrees[i], dim, q_extra,
306*9b072555Sjeremylt                                 dim, num_comp_u, g_size[i], xl_size[i], bp_options[bp_choice],
307*9b072555Sjeremylt                                 ceed_data[i], i==(fine_level), rhs_ceed, &target);
308f9342789SJeremy L Thompson     CHKERRQ(ierr);
3096c5df90dSjeremylt   }
3106c5df90dSjeremylt 
3116c5df90dSjeremylt   // Gather RHS
312*9b072555Sjeremylt   CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL);
313*9b072555Sjeremylt   ierr = VecRestoreArrayAndMemType(rhs_loc, &r); CHKERRQ(ierr);
3146c5df90dSjeremylt   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
315*9b072555Sjeremylt   ierr = DMLocalToGlobal(dm[fine_level], rhs_loc, ADD_VALUES, rhs); CHKERRQ(ierr);
316*9b072555Sjeremylt   CeedVectorDestroy(&rhs_ceed);
3176c5df90dSjeremylt 
31843eb8658SJeremy L Thompson   // Create the restriction/interpolation QFunction
319*9b072555Sjeremylt   CeedQFunctionCreateIdentity(ceed, num_comp_u, CEED_EVAL_NONE, CEED_EVAL_INTERP,
320*9b072555Sjeremylt                               &qf_restrict);
321*9b072555Sjeremylt   CeedQFunctionCreateIdentity(ceed, num_comp_u, CEED_EVAL_INTERP, CEED_EVAL_NONE,
322*9b072555Sjeremylt                               &qf_prolong);
3236c5df90dSjeremylt 
3246c5df90dSjeremylt   // Set up libCEED level transfer operators
325*9b072555Sjeremylt   ierr = CeedLevelTransferSetup(ceed, num_levels, num_comp_u, ceed_data,
326*9b072555Sjeremylt                                 level_degrees,
327*9b072555Sjeremylt                                 qf_restrict, qf_prolong); CHKERRQ(ierr);
3286c5df90dSjeremylt 
32943eb8658SJeremy L Thompson   // Create the error QFunction
330*9b072555Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].error,
331*9b072555Sjeremylt                               bp_options[bp_choice].error_loc, &qf_error);
332*9b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP);
333*9b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE);
334*9b072555Sjeremylt   CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_NONE);
3356c5df90dSjeremylt 
3366c5df90dSjeremylt   // Create the error operator
337*9b072555Sjeremylt   CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
338*9b072555Sjeremylt                      &op_error);
339*9b072555Sjeremylt   CeedOperatorSetField(op_error, "u", ceed_data[fine_level]->elem_restr_u,
340*9b072555Sjeremylt                        ceed_data[fine_level]->basis_u, CEED_VECTOR_ACTIVE);
341*9b072555Sjeremylt   CeedOperatorSetField(op_error, "true_soln",
342*9b072555Sjeremylt                        ceed_data[fine_level]->elem_restr_u_i,
343a8d32208Sjeremylt                        CEED_BASIS_COLLOCATED, target);
344*9b072555Sjeremylt   CeedOperatorSetField(op_error, "error", ceed_data[fine_level]->elem_restr_u_i,
345a8d32208Sjeremylt                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
3466c5df90dSjeremylt 
3476c5df90dSjeremylt   // Calculate multiplicity
348*9b072555Sjeremylt   for (int i=0; i<num_levels; i++) {
3496c5df90dSjeremylt     PetscScalar *x;
3506c5df90dSjeremylt 
3516c5df90dSjeremylt     // CEED vector
352*9b072555Sjeremylt     ierr = VecZeroEntries(X_loc[i]); CHKERRQ(ierr);
353*9b072555Sjeremylt     ierr = VecGetArray(X_loc[i], &x); CHKERRQ(ierr);
354*9b072555Sjeremylt     CeedVectorSetArray(ceed_data[i]->x_ceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
3556c5df90dSjeremylt 
3566c5df90dSjeremylt     // Multiplicity
357*9b072555Sjeremylt     CeedElemRestrictionGetMultiplicity(ceed_data[i]->elem_restr_u,
358*9b072555Sjeremylt                                        ceed_data[i]->x_ceed);
359*9b072555Sjeremylt     CeedVectorSyncArray(ceed_data[i]->x_ceed, CEED_MEM_HOST);
3606c5df90dSjeremylt 
3616c5df90dSjeremylt     // Restore vector
362*9b072555Sjeremylt     ierr = VecRestoreArray(X_loc[i], &x); CHKERRQ(ierr);
3636c5df90dSjeremylt 
3646c5df90dSjeremylt     // Creat mult vector
365*9b072555Sjeremylt     ierr = VecDuplicate(X_loc[i], &mult[i]); CHKERRQ(ierr);
3666c5df90dSjeremylt 
3676c5df90dSjeremylt     // Local-to-global
3686c5df90dSjeremylt     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
369*9b072555Sjeremylt     ierr = DMLocalToGlobal(dm[i], X_loc[i], ADD_VALUES, X[i]);
3706c5df90dSjeremylt     CHKERRQ(ierr);
371*9b072555Sjeremylt     ierr = VecZeroEntries(X_loc[i]); CHKERRQ(ierr);
3726c5df90dSjeremylt 
3736c5df90dSjeremylt     // Global-to-local
374483f8b0dSjeremylt     ierr = DMGlobalToLocal(dm[i], X[i], INSERT_VALUES, mult[i]);
3756c5df90dSjeremylt     CHKERRQ(ierr);
3766c5df90dSjeremylt     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
3776c5df90dSjeremylt 
3786c5df90dSjeremylt     // Multiplicity scaling
3796c5df90dSjeremylt     ierr = VecReciprocal(mult[i]);
3806c5df90dSjeremylt   }
3816c5df90dSjeremylt 
3826c5df90dSjeremylt   // Set up Mat
383*9b072555Sjeremylt   for (int i=0; i<num_levels; i++) {
384226c3a8fSjeremylt     // User Operator
385*9b072555Sjeremylt     user_O[i]->comm = comm;
386*9b072555Sjeremylt     user_O[i]->dm = dm[i];
387*9b072555Sjeremylt     user_O[i]->X_loc = X_loc[i];
388*9b072555Sjeremylt     ierr = VecDuplicate(X_loc[i], &user_O[i]->Y_loc); CHKERRQ(ierr);
389*9b072555Sjeremylt     user_O[i]->x_ceed = ceed_data[i]->x_ceed;
390*9b072555Sjeremylt     user_O[i]->y_ceed = ceed_data[i]->y_ceed;
391*9b072555Sjeremylt     user_O[i]->op = ceed_data[i]->op_apply;
392*9b072555Sjeremylt     user_O[i]->ceed = ceed;
3936c5df90dSjeremylt 
3946c5df90dSjeremylt     if (i > 0) {
395a97643b0Sjeremylt       // Prolongation/Restriction Operator
396*9b072555Sjeremylt       user_pr[i]->comm = comm;
397*9b072555Sjeremylt       user_pr[i]->dmf = dm[i];
398*9b072555Sjeremylt       user_pr[i]->dmc = dm[i-1];
399*9b072555Sjeremylt       user_pr[i]->loc_vec_c = X_loc[i-1];
400*9b072555Sjeremylt       user_pr[i]->loc_vec_f = user_O[i]->Y_loc;
401*9b072555Sjeremylt       user_pr[i]->mult_vec = mult[i];
402*9b072555Sjeremylt       user_pr[i]->ceed_vec_c = user_O[i-1]->x_ceed;
403*9b072555Sjeremylt       user_pr[i]->ceed_vec_f = user_O[i]->y_ceed;
404*9b072555Sjeremylt       user_pr[i]->op_prolong = ceed_data[i]->op_prolong;
405*9b072555Sjeremylt       user_pr[i]->op_restrict = ceed_data[i]->op_restrict;
406*9b072555Sjeremylt       user_pr[i]->ceed = ceed;
4076c5df90dSjeremylt     }
4086c5df90dSjeremylt   }
4096c5df90dSjeremylt 
41015ce0ef0Sjeremylt   // Setup dummy SNES for AMG coarse solve
411*9b072555Sjeremylt   ierr = SNESCreate(comm, &snes_dummy); CHKERRQ(ierr);
412*9b072555Sjeremylt   ierr = SNESSetDM(snes_dummy, dm[0]); CHKERRQ(ierr);
413*9b072555Sjeremylt   ierr = SNESSetSolution(snes_dummy, X[0]); CHKERRQ(ierr);
41415ce0ef0Sjeremylt 
41515ce0ef0Sjeremylt   // -- Jacobian matrix
41615ce0ef0Sjeremylt   ierr = DMSetMatType(dm[0], MATAIJ); CHKERRQ(ierr);
417*9b072555Sjeremylt   ierr = DMCreateMatrix(dm[0], &mat_coarse); CHKERRQ(ierr);
418*9b072555Sjeremylt   ierr = SNESSetJacobian(snes_dummy, mat_coarse, mat_coarse, NULL,
41915ce0ef0Sjeremylt                          NULL); CHKERRQ(ierr);
42015ce0ef0Sjeremylt 
42115ce0ef0Sjeremylt   // -- Residual evaluation function
422*9b072555Sjeremylt   ierr = SNESSetFunction(snes_dummy, X[0], FormResidual_Ceed,
423*9b072555Sjeremylt                          user_O[0]); CHKERRQ(ierr);
42415ce0ef0Sjeremylt 
42515ce0ef0Sjeremylt   // -- Form Jacobian
426*9b072555Sjeremylt   ierr = SNESComputeJacobianDefaultColor(snes_dummy, X[0], mat_O[0],
427*9b072555Sjeremylt                                          mat_coarse, NULL); CHKERRQ(ierr);
42815ce0ef0Sjeremylt 
4296c5df90dSjeremylt   // Set up KSP
4306c5df90dSjeremylt   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
4316c5df90dSjeremylt   {
4326c5df90dSjeremylt     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
4336c5df90dSjeremylt     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
4346c5df90dSjeremylt     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
4356c5df90dSjeremylt                             PETSC_DEFAULT); CHKERRQ(ierr);
4366c5df90dSjeremylt   }
4376c5df90dSjeremylt   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
438*9b072555Sjeremylt   ierr = KSPSetOperators(ksp, mat_O[fine_level], mat_O[fine_level]);
4396c5df90dSjeremylt   CHKERRQ(ierr);
4406c5df90dSjeremylt 
4416c5df90dSjeremylt   // Set up PCMG
4426c5df90dSjeremylt   ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
443*9b072555Sjeremylt   PCMGCycleType pcmg_cycle_type = PC_MG_CYCLE_V;
4446c5df90dSjeremylt   {
4456c5df90dSjeremylt     ierr = PCSetType(pc, PCMG); CHKERRQ(ierr);
4466c5df90dSjeremylt 
4476c5df90dSjeremylt     // PCMG levels
448*9b072555Sjeremylt     ierr = PCMGSetLevels(pc, num_levels, NULL); CHKERRQ(ierr);
449*9b072555Sjeremylt     for (int i=0; i<num_levels; i++) {
4506c5df90dSjeremylt       // Smoother
4516c5df90dSjeremylt       KSP smoother;
4526c5df90dSjeremylt       PC smoother_pc;
4536c5df90dSjeremylt       ierr = PCMGGetSmoother(pc, i, &smoother); CHKERRQ(ierr);
4546c5df90dSjeremylt       ierr = KSPSetType(smoother, KSPCHEBYSHEV); CHKERRQ(ierr);
4556c5df90dSjeremylt       ierr = KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1); CHKERRQ(ierr);
4566c5df90dSjeremylt       ierr = KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE); CHKERRQ(ierr);
457*9b072555Sjeremylt       ierr = KSPSetOperators(smoother, mat_O[i], mat_O[i]); CHKERRQ(ierr);
4586c5df90dSjeremylt       ierr = KSPGetPC(smoother, &smoother_pc); CHKERRQ(ierr);
4596c5df90dSjeremylt       ierr = PCSetType(smoother_pc, PCJACOBI); CHKERRQ(ierr);
4606c5df90dSjeremylt       ierr = PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
4616c5df90dSjeremylt 
4626c5df90dSjeremylt       // Work vector
463*9b072555Sjeremylt       if (i < num_levels - 1) {
4646c5df90dSjeremylt         ierr = PCMGSetX(pc, i, X[i]); CHKERRQ(ierr);
4656c5df90dSjeremylt       }
4666c5df90dSjeremylt 
4676c5df90dSjeremylt       // Level transfers
4686c5df90dSjeremylt       if (i > 0) {
4696c5df90dSjeremylt         // Interpolation
470*9b072555Sjeremylt         ierr = PCMGSetInterpolation(pc, i, mat_pr[i]); CHKERRQ(ierr);
4716c5df90dSjeremylt       }
4726c5df90dSjeremylt 
4736c5df90dSjeremylt       // Coarse solve
4746c5df90dSjeremylt       KSP coarse;
4756c5df90dSjeremylt       PC coarse_pc;
4766c5df90dSjeremylt       ierr = PCMGGetCoarseSolve(pc, &coarse); CHKERRQ(ierr);
47715ce0ef0Sjeremylt       ierr = KSPSetType(coarse, KSPPREONLY); CHKERRQ(ierr);
478*9b072555Sjeremylt       ierr = KSPSetOperators(coarse, mat_coarse, mat_coarse); CHKERRQ(ierr);
47915ce0ef0Sjeremylt 
4806c5df90dSjeremylt       ierr = KSPGetPC(coarse, &coarse_pc); CHKERRQ(ierr);
48115ce0ef0Sjeremylt       ierr = PCSetType(coarse_pc, PCGAMG); CHKERRQ(ierr);
48215ce0ef0Sjeremylt 
48315ce0ef0Sjeremylt       ierr = KSPSetOptionsPrefix(coarse, "coarse_"); CHKERRQ(ierr);
48415ce0ef0Sjeremylt       ierr = PCSetOptionsPrefix(coarse_pc, "coarse_"); CHKERRQ(ierr);
48515ce0ef0Sjeremylt       ierr = KSPSetFromOptions(coarse); CHKERRQ(ierr);
48615ce0ef0Sjeremylt       ierr = PCSetFromOptions(coarse_pc); CHKERRQ(ierr);
4876c5df90dSjeremylt     }
4886c5df90dSjeremylt 
4896c5df90dSjeremylt     // PCMG options
4906c5df90dSjeremylt     ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr);
4916c5df90dSjeremylt     ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr);
492*9b072555Sjeremylt     ierr = PCMGSetCycleType(pc, pcmg_cycle_type); CHKERRQ(ierr);
4936c5df90dSjeremylt   }
4946c5df90dSjeremylt 
4956c5df90dSjeremylt   // First run, if benchmarking
4966c5df90dSjeremylt   if (benchmark_mode) {
4976c5df90dSjeremylt     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
4986c5df90dSjeremylt     CHKERRQ(ierr);
499*9b072555Sjeremylt     ierr = VecZeroEntries(X[fine_level]); CHKERRQ(ierr);
5006c5df90dSjeremylt     my_rt_start = MPI_Wtime();
501*9b072555Sjeremylt     ierr = KSPSolve(ksp, rhs, X[fine_level]); CHKERRQ(ierr);
5026c5df90dSjeremylt     my_rt = MPI_Wtime() - my_rt_start;
5036c5df90dSjeremylt     ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
5046c5df90dSjeremylt     CHKERRQ(ierr);
5056c5df90dSjeremylt     // Set maxits based on first iteration timing
5066c5df90dSjeremylt     if (my_rt > 0.02) {
5076c5df90dSjeremylt       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
5086c5df90dSjeremylt       CHKERRQ(ierr);
5096c5df90dSjeremylt     } else {
5106c5df90dSjeremylt       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
5116c5df90dSjeremylt       CHKERRQ(ierr);
5126c5df90dSjeremylt     }
5136c5df90dSjeremylt   }
5146c5df90dSjeremylt 
5156c5df90dSjeremylt   // Timed solve
516*9b072555Sjeremylt   ierr = VecZeroEntries(X[fine_level]); CHKERRQ(ierr);
5176c5df90dSjeremylt   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
51809a940d7Sjeremylt 
51909a940d7Sjeremylt   // -- Performance logging
520*9b072555Sjeremylt   ierr = PetscLogStageRegister("Solve Stage", &solve_stage); CHKERRQ(ierr);
521*9b072555Sjeremylt   ierr = PetscLogStagePush(solve_stage); CHKERRQ(ierr);
52209a940d7Sjeremylt 
52309a940d7Sjeremylt   // -- Solve
5246c5df90dSjeremylt   my_rt_start = MPI_Wtime();
525*9b072555Sjeremylt   ierr = KSPSolve(ksp, rhs, X[fine_level]); CHKERRQ(ierr);
5266c5df90dSjeremylt   my_rt = MPI_Wtime() - my_rt_start;
5276c5df90dSjeremylt 
52809a940d7Sjeremylt 
52909a940d7Sjeremylt   // -- Performance logging
53009a940d7Sjeremylt   ierr = PetscLogStagePop();
53109a940d7Sjeremylt 
5326c5df90dSjeremylt   // Output results
5336c5df90dSjeremylt   {
534*9b072555Sjeremylt     KSPType ksp_type;
535*9b072555Sjeremylt     PCMGType pcmg_type;
5366c5df90dSjeremylt     KSPConvergedReason reason;
5376c5df90dSjeremylt     PetscReal rnorm;
5386c5df90dSjeremylt     PetscInt its;
539*9b072555Sjeremylt     ierr = KSPGetType(ksp, &ksp_type); CHKERRQ(ierr);
5406c5df90dSjeremylt     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
5416c5df90dSjeremylt     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
5426c5df90dSjeremylt     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
543*9b072555Sjeremylt     ierr = PCMGGetType(pc, &pcmg_type); CHKERRQ(ierr);
5446c5df90dSjeremylt     if (!test_mode || reason < 0 || rnorm > 1e-8) {
5456c5df90dSjeremylt       ierr = PetscPrintf(comm,
5466c5df90dSjeremylt                          "  KSP:\n"
5476c5df90dSjeremylt                          "    KSP Type                           : %s\n"
5486c5df90dSjeremylt                          "    KSP Convergence                    : %s\n"
5496c5df90dSjeremylt                          "    Total KSP Iterations               : %D\n"
5506c5df90dSjeremylt                          "    Final rnorm                        : %e\n",
551*9b072555Sjeremylt                          ksp_type, KSPConvergedReasons[reason], its,
5526c5df90dSjeremylt                          (double)rnorm); CHKERRQ(ierr);
5536c5df90dSjeremylt       ierr = PetscPrintf(comm,
5546c5df90dSjeremylt                          "  PCMG:\n"
5556c5df90dSjeremylt                          "    PCMG Type                          : %s\n"
5566c5df90dSjeremylt                          "    PCMG Cycle Type                    : %s\n",
557*9b072555Sjeremylt                          PCMGTypes[pcmg_type],
558*9b072555Sjeremylt                          PCMGCycleTypes[pcmg_cycle_type]); CHKERRQ(ierr);
5596c5df90dSjeremylt     }
5606c5df90dSjeremylt     if (!test_mode) {
5616c5df90dSjeremylt       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
5626c5df90dSjeremylt     }
5636c5df90dSjeremylt     {
564*9b072555Sjeremylt       PetscReal max_error;
565*9b072555Sjeremylt       ierr = ComputeErrorMax(user_O[fine_level], op_error, X[fine_level], target,
566*9b072555Sjeremylt                              &max_error); CHKERRQ(ierr);
5676c5df90dSjeremylt       PetscReal tol = 5e-2;
568*9b072555Sjeremylt       if (!test_mode || max_error > tol) {
5696c5df90dSjeremylt         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
5706c5df90dSjeremylt         CHKERRQ(ierr);
5716c5df90dSjeremylt         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
5726c5df90dSjeremylt         CHKERRQ(ierr);
5736c5df90dSjeremylt         ierr = PetscPrintf(comm,
5746c5df90dSjeremylt                            "    Pointwise Error (max)              : %e\n"
5756c5df90dSjeremylt                            "    CG Solve Time                      : %g (%g) sec\n",
576*9b072555Sjeremylt                            (double)max_error, rt_max, rt_min); CHKERRQ(ierr);
5776c5df90dSjeremylt       }
5786c5df90dSjeremylt     }
5796c5df90dSjeremylt     if (benchmark_mode && (!test_mode)) {
5806c5df90dSjeremylt       ierr = PetscPrintf(comm,
5816c5df90dSjeremylt                          "    DoFs/Sec in CG                     : %g (%g) million\n",
582*9b072555Sjeremylt                          1e-6*g_size[fine_level]*its/rt_max,
583*9b072555Sjeremylt                          1e-6*g_size[fine_level]*its/rt_min);
5846c5df90dSjeremylt       CHKERRQ(ierr);
5856c5df90dSjeremylt     }
5866c5df90dSjeremylt   }
5876c5df90dSjeremylt 
5886c5df90dSjeremylt   if (write_solution) {
589*9b072555Sjeremylt     PetscViewer vtk_viewer_soln;
5906c5df90dSjeremylt 
591*9b072555Sjeremylt     ierr = PetscViewerCreate(comm, &vtk_viewer_soln); CHKERRQ(ierr);
592*9b072555Sjeremylt     ierr = PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK); CHKERRQ(ierr);
593*9b072555Sjeremylt     ierr = PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"); CHKERRQ(ierr);
594*9b072555Sjeremylt     ierr = VecView(X[fine_level], vtk_viewer_soln); CHKERRQ(ierr);
595*9b072555Sjeremylt     ierr = PetscViewerDestroy(&vtk_viewer_soln); CHKERRQ(ierr);
5966c5df90dSjeremylt   }
5976c5df90dSjeremylt 
5986c5df90dSjeremylt   // Cleanup
599*9b072555Sjeremylt   for (int i=0; i<num_levels; i++) {
6006c5df90dSjeremylt     ierr = VecDestroy(&X[i]); CHKERRQ(ierr);
601*9b072555Sjeremylt     ierr = VecDestroy(&X_loc[i]); CHKERRQ(ierr);
6026c5df90dSjeremylt     ierr = VecDestroy(&mult[i]); CHKERRQ(ierr);
603*9b072555Sjeremylt     ierr = VecDestroy(&user_O[i]->Y_loc); CHKERRQ(ierr);
604*9b072555Sjeremylt     ierr = MatDestroy(&mat_O[i]); CHKERRQ(ierr);
605*9b072555Sjeremylt     ierr = PetscFree(user_O[i]); CHKERRQ(ierr);
6066c5df90dSjeremylt     if (i > 0) {
607*9b072555Sjeremylt       ierr = MatDestroy(&mat_pr[i]); CHKERRQ(ierr);
608*9b072555Sjeremylt       ierr = PetscFree(user_pr[i]); CHKERRQ(ierr);
6096c5df90dSjeremylt     }
610*9b072555Sjeremylt     ierr = CeedDataDestroy(i, ceed_data[i]); CHKERRQ(ierr);
6116c5df90dSjeremylt     ierr = DMDestroy(&dm[i]); CHKERRQ(ierr);
6126c5df90dSjeremylt   }
613*9b072555Sjeremylt   ierr = PetscFree(level_degrees); CHKERRQ(ierr);
6146c5df90dSjeremylt   ierr = PetscFree(dm); CHKERRQ(ierr);
6156c5df90dSjeremylt   ierr = PetscFree(X); CHKERRQ(ierr);
616*9b072555Sjeremylt   ierr = PetscFree(X_loc); CHKERRQ(ierr);
6176c5df90dSjeremylt   ierr = PetscFree(mult); CHKERRQ(ierr);
618*9b072555Sjeremylt   ierr = PetscFree(mat_O); CHKERRQ(ierr);
619*9b072555Sjeremylt   ierr = PetscFree(mat_pr); CHKERRQ(ierr);
620*9b072555Sjeremylt   ierr = PetscFree(ceed_data); CHKERRQ(ierr);
621*9b072555Sjeremylt   ierr = PetscFree(user_O); CHKERRQ(ierr);
622*9b072555Sjeremylt   ierr = PetscFree(user_pr); CHKERRQ(ierr);
623*9b072555Sjeremylt   ierr = PetscFree(l_size); CHKERRQ(ierr);
624*9b072555Sjeremylt   ierr = PetscFree(xl_size); CHKERRQ(ierr);
625*9b072555Sjeremylt   ierr = PetscFree(g_size); CHKERRQ(ierr);
6266c5df90dSjeremylt   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
627*9b072555Sjeremylt   ierr = VecDestroy(&rhs_loc); CHKERRQ(ierr);
628*9b072555Sjeremylt   ierr = MatDestroy(&mat_coarse); CHKERRQ(ierr);
6296c5df90dSjeremylt   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
630*9b072555Sjeremylt   ierr = SNESDestroy(&snes_dummy); CHKERRQ(ierr);
631*9b072555Sjeremylt   ierr = DMDestroy(&dm_orig); CHKERRQ(ierr);
6326c5df90dSjeremylt   CeedVectorDestroy(&target);
633*9b072555Sjeremylt   CeedQFunctionDestroy(&qf_error);
634*9b072555Sjeremylt   CeedQFunctionDestroy(&qf_restrict);
635*9b072555Sjeremylt   CeedQFunctionDestroy(&qf_prolong);
636*9b072555Sjeremylt   CeedOperatorDestroy(&op_error);
6376c5df90dSjeremylt   CeedDestroy(&ceed);
6386c5df90dSjeremylt   return PetscFinalize();
6396c5df90dSjeremylt }
640