xref: /libCEED/examples/petsc/multigrid.c (revision f9342789ef8d4c153f494bcb8a6c8db3c5c40df0)
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
316c5df90dSjeremylt //     multigrid -problem bp4 -ceed /cpu/self
326c5df90dSjeremylt //     multigrid -problem bp5 -ceed /cpu/occa
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 
416c5df90dSjeremylt #define multigrid
426c5df90dSjeremylt #include "setup.h"
436c5df90dSjeremylt 
446c5df90dSjeremylt int main(int argc, char **argv) {
456c5df90dSjeremylt   PetscInt ierr;
466c5df90dSjeremylt   MPI_Comm comm;
47cb0b5415Sjeremylt   char filename[PETSC_MAX_PATH_LEN],
48cb0b5415Sjeremylt        ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self";
496c5df90dSjeremylt   double my_rt_start, my_rt, rt_min, rt_max;
5061608365Sjeremylt   PetscInt degree = 3, qextra, *lsize, *xlsize, *gsize, dim = 3, fineLevel,
516c5df90dSjeremylt            melem[3] = {3, 3, 3}, ncompu = 1, numlevels = degree, *leveldegrees;
526c5df90dSjeremylt   PetscScalar *r;
536c5df90dSjeremylt   PetscBool test_mode, benchmark_mode, read_mesh, write_solution;
5409a940d7Sjeremylt   PetscLogStage solvestage;
55199551f5Sjeremylt   DM  *dm, dmorig;
56c70bd2a0Sjeremylt   SNES snesdummy;
576c5df90dSjeremylt   KSP ksp;
586c5df90dSjeremylt   PC pc;
59199551f5Sjeremylt   Mat *matO, *matPR, matcoarse;
60cdb3667fSjeremylt   Vec *X, *Xloc, *mult, rhs, rhsloc;
616c5df90dSjeremylt   UserO *userO;
62a97643b0Sjeremylt   UserProlongRestr *userPR;
636c5df90dSjeremylt   Ceed ceed;
646c5df90dSjeremylt   CeedData *ceeddata;
659396343dSjeremylt   CeedMemType memtyperequested;
665e81177dSjeremylt   CeedVector rhsceed, target;
67c70bd2a0Sjeremylt   CeedQFunction qferror, qfrestrict, qfprolong;
68c70bd2a0Sjeremylt   CeedOperator operror;
69199551f5Sjeremylt   bpType bpchoice;
706c5df90dSjeremylt   coarsenType coarsen;
716c5df90dSjeremylt 
726c5df90dSjeremylt   ierr = PetscInitialize(&argc, &argv, NULL, help);
736c5df90dSjeremylt   if (ierr) return ierr;
746c5df90dSjeremylt   comm = PETSC_COMM_WORLD;
756c5df90dSjeremylt 
769396343dSjeremylt   // Check PETSc CUDA avaliability
779396343dSjeremylt   PetscBool petschavecuda, setmemtyperequest = PETSC_FALSE;
789396343dSjeremylt   // *INDENT-OFF*
799396343dSjeremylt   #ifdef PETSC_HAVE_CUDA
809396343dSjeremylt   petschavecuda = PETSC_TRUE;
819396343dSjeremylt   #else
829396343dSjeremylt   petschavecuda = PETSC_FALSE;
839396343dSjeremylt   #endif
849396343dSjeremylt   // *INDENT-ON*
859396343dSjeremylt 
866c5df90dSjeremylt   // Parse command line options
876c5df90dSjeremylt   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr);
88199551f5Sjeremylt   bpchoice = CEED_BP3;
896c5df90dSjeremylt   ierr = PetscOptionsEnum("-problem",
906c5df90dSjeremylt                           "CEED benchmark problem to solve", NULL,
91199551f5Sjeremylt                           bpTypes, (PetscEnum)bpchoice, (PetscEnum *)&bpchoice,
926c5df90dSjeremylt                           NULL); CHKERRQ(ierr);
93199551f5Sjeremylt   ncompu = bpOptions[bpchoice].ncompu;
946c5df90dSjeremylt   test_mode = PETSC_FALSE;
956c5df90dSjeremylt   ierr = PetscOptionsBool("-test",
966c5df90dSjeremylt                           "Testing mode (do not print unless error is large)",
976c5df90dSjeremylt                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
986c5df90dSjeremylt   benchmark_mode = PETSC_FALSE;
996c5df90dSjeremylt   ierr = PetscOptionsBool("-benchmark",
1006c5df90dSjeremylt                           "Benchmarking mode (prints benchmark statistics)",
1016c5df90dSjeremylt                           NULL, benchmark_mode, &benchmark_mode, NULL);
1026c5df90dSjeremylt   CHKERRQ(ierr);
1036c5df90dSjeremylt   write_solution = PETSC_FALSE;
1046c5df90dSjeremylt   ierr = PetscOptionsBool("-write_solution",
1056c5df90dSjeremylt                           "Write solution for visualization",
1066c5df90dSjeremylt                           NULL, write_solution, &write_solution, NULL);
1076c5df90dSjeremylt   CHKERRQ(ierr);
1086c5df90dSjeremylt   degree = test_mode ? 3 : 2;
1096c5df90dSjeremylt   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
1106c5df90dSjeremylt                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
1116c5df90dSjeremylt   if (degree < 1) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
1126c5df90dSjeremylt                              "-degree %D must be at least 1", degree);
113199551f5Sjeremylt   qextra = bpOptions[bpchoice].qextra;
1146c5df90dSjeremylt   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
1156c5df90dSjeremylt                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
1166c5df90dSjeremylt   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
1176c5df90dSjeremylt                             NULL, ceedresource, ceedresource,
1186c5df90dSjeremylt                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
1196c5df90dSjeremylt   coarsen = COARSEN_UNIFORM;
1206c5df90dSjeremylt   ierr = PetscOptionsEnum("-coarsen",
1216c5df90dSjeremylt                           "Coarsening strategy to use", NULL,
1226c5df90dSjeremylt                           coarsenTypes, (PetscEnum)coarsen,
1236c5df90dSjeremylt                           (PetscEnum *)&coarsen, NULL); CHKERRQ(ierr);
124cb32e2e7SValeria Barra   read_mesh = PETSC_FALSE;
1256c5df90dSjeremylt   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL,
1266c5df90dSjeremylt                             filename, filename, sizeof(filename), &read_mesh);
1276c5df90dSjeremylt   CHKERRQ(ierr);
1286c5df90dSjeremylt   if (!read_mesh) {
1296c5df90dSjeremylt     PetscInt tmp = dim;
1306c5df90dSjeremylt     ierr = PetscOptionsIntArray("-cells","Number of cells per dimension", NULL,
1316c5df90dSjeremylt                                 melem, &tmp, NULL); CHKERRQ(ierr);
1326c5df90dSjeremylt   }
1339396343dSjeremylt   memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST;
1349396343dSjeremylt   ierr = PetscOptionsEnum("-memtype",
1359396343dSjeremylt                           "CEED MemType requested", NULL,
1369396343dSjeremylt                           memTypes, (PetscEnum)memtyperequested,
1379396343dSjeremylt                           (PetscEnum *)&memtyperequested, &setmemtyperequest);
1389396343dSjeremylt   CHKERRQ(ierr);
1396c5df90dSjeremylt   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
1406c5df90dSjeremylt 
1419396343dSjeremylt   // Set up libCEED
1429396343dSjeremylt   CeedInit(ceedresource, &ceed);
1439396343dSjeremylt   CeedMemType memtypebackend;
1449396343dSjeremylt   CeedGetPreferredMemType(ceed, &memtypebackend);
1459396343dSjeremylt 
1469396343dSjeremylt   // Check memtype compatibility
1479396343dSjeremylt   if (!setmemtyperequest)
1489396343dSjeremylt     memtyperequested = memtypebackend;
1499396343dSjeremylt   else if (!petschavecuda && memtyperequested == CEED_MEM_DEVICE)
1509396343dSjeremylt     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS,
1519396343dSjeremylt              "PETSc was not built with CUDA. "
1529396343dSjeremylt              "Requested MemType CEED_MEM_DEVICE is not supported.", NULL);
1539396343dSjeremylt 
1546c5df90dSjeremylt   // Setup DM
1556c5df90dSjeremylt   if (read_mesh) {
156199551f5Sjeremylt     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dmorig);
1576c5df90dSjeremylt     CHKERRQ(ierr);
1586c5df90dSjeremylt   } else {
1596c5df90dSjeremylt     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_FALSE, melem, NULL,
160199551f5Sjeremylt                                NULL, NULL, PETSC_TRUE,&dmorig); CHKERRQ(ierr);
1616c5df90dSjeremylt   }
1626c5df90dSjeremylt 
1636c5df90dSjeremylt   {
1646c5df90dSjeremylt     DM dmDist = NULL;
1656c5df90dSjeremylt     PetscPartitioner part;
1666c5df90dSjeremylt 
167199551f5Sjeremylt     ierr = DMPlexGetPartitioner(dmorig, &part); CHKERRQ(ierr);
1686c5df90dSjeremylt     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
169199551f5Sjeremylt     ierr = DMPlexDistribute(dmorig, 0, NULL, &dmDist); CHKERRQ(ierr);
1706c5df90dSjeremylt     if (dmDist) {
171199551f5Sjeremylt       ierr = DMDestroy(&dmorig); CHKERRQ(ierr);
172199551f5Sjeremylt       dmorig = dmDist;
1736c5df90dSjeremylt     }
1746c5df90dSjeremylt   }
1756c5df90dSjeremylt 
1766c5df90dSjeremylt   // Allocate arrays for PETSc objects for each level
1776c5df90dSjeremylt   switch (coarsen) {
1786c5df90dSjeremylt   case COARSEN_UNIFORM:
1796c5df90dSjeremylt     numlevels = degree;
1806c5df90dSjeremylt     break;
181dc7d240cSValeria Barra   case COARSEN_LOGARITHMIC:
1826c5df90dSjeremylt     numlevels = ceil(log(degree)/log(2)) + 1;
1836c5df90dSjeremylt     break;
1846c5df90dSjeremylt   }
1856c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &leveldegrees); CHKERRQ(ierr);
18661608365Sjeremylt   fineLevel = numlevels - 1;
18761608365Sjeremylt 
1886c5df90dSjeremylt   switch (coarsen) {
1896c5df90dSjeremylt   case COARSEN_UNIFORM:
1906c5df90dSjeremylt     for (int i=0; i<numlevels; i++) leveldegrees[i] = i + 1;
1916c5df90dSjeremylt     break;
192dc7d240cSValeria Barra   case COARSEN_LOGARITHMIC:
1936c5df90dSjeremylt     for (int i=0; i<numlevels - 1; i++) leveldegrees[i] = pow(2,i);
19461608365Sjeremylt     leveldegrees[fineLevel] = degree;
1956c5df90dSjeremylt     break;
1966c5df90dSjeremylt   }
1976c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &dm); CHKERRQ(ierr);
1986c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &X); CHKERRQ(ierr);
1996c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &Xloc); CHKERRQ(ierr);
2006c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &mult); CHKERRQ(ierr);
2016c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &userO); CHKERRQ(ierr);
202a97643b0Sjeremylt   ierr = PetscMalloc1(numlevels, &userPR); CHKERRQ(ierr);
2036c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &matO); CHKERRQ(ierr);
204a97643b0Sjeremylt   ierr = PetscMalloc1(numlevels, &matPR); CHKERRQ(ierr);
2056c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &lsize); CHKERRQ(ierr);
2066c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &xlsize); CHKERRQ(ierr);
2076c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &gsize); CHKERRQ(ierr);
2086c5df90dSjeremylt 
2096c5df90dSjeremylt   // Setup DM and Operator Mat Shells for each level
2106c5df90dSjeremylt   for (CeedInt i=0; i<numlevels; i++) {
2116c5df90dSjeremylt     // Create DM
212199551f5Sjeremylt     ierr = DMClone(dmorig, &dm[i]); CHKERRQ(ierr);
213199551f5Sjeremylt     ierr = SetupDMByDegree(dm[i], leveldegrees[i], ncompu, bpchoice);
2146c5df90dSjeremylt     CHKERRQ(ierr);
2156c5df90dSjeremylt 
2166c5df90dSjeremylt     // Create vectors
2179396343dSjeremylt     if (memtyperequested == CEED_MEM_DEVICE) {
2189396343dSjeremylt       ierr = DMSetVecType(dm[i], VECCUDA); CHKERRQ(ierr);
2199396343dSjeremylt     }
2206c5df90dSjeremylt     ierr = DMCreateGlobalVector(dm[i], &X[i]); CHKERRQ(ierr);
2216c5df90dSjeremylt     ierr = VecGetLocalSize(X[i], &lsize[i]); CHKERRQ(ierr);
2226c5df90dSjeremylt     ierr = VecGetSize(X[i], &gsize[i]); CHKERRQ(ierr);
2236c5df90dSjeremylt     ierr = DMCreateLocalVector(dm[i], &Xloc[i]); CHKERRQ(ierr);
2246c5df90dSjeremylt     ierr = VecGetSize(Xloc[i], &xlsize[i]); CHKERRQ(ierr);
2256c5df90dSjeremylt 
2266c5df90dSjeremylt     // Operator
2276c5df90dSjeremylt     ierr = PetscMalloc1(1, &userO[i]); CHKERRQ(ierr);
2286c5df90dSjeremylt     ierr = MatCreateShell(comm, lsize[i], lsize[i], gsize[i], gsize[i],
2296c5df90dSjeremylt                           userO[i], &matO[i]); CHKERRQ(ierr);
2306c5df90dSjeremylt     ierr = MatShellSetOperation(matO[i], MATOP_MULT,
231ce74dcefSjeremylt                                 (void(*)(void))MatMult_Ceed); CHKERRQ(ierr);
2326c5df90dSjeremylt     ierr = MatShellSetOperation(matO[i], MATOP_GET_DIAGONAL,
233ce74dcefSjeremylt                                 (void(*)(void))MatGetDiag); CHKERRQ(ierr);
2349396343dSjeremylt     if (memtyperequested == CEED_MEM_DEVICE) {
2359396343dSjeremylt       ierr = MatShellSetVecType(matO[i], VECCUDA); CHKERRQ(ierr);
2369396343dSjeremylt     }
2376c5df90dSjeremylt 
2386c5df90dSjeremylt     // Level transfers
2396c5df90dSjeremylt     if (i > 0) {
2406c5df90dSjeremylt       // Interp
241a97643b0Sjeremylt       ierr = PetscMalloc1(1, &userPR[i]); CHKERRQ(ierr);
2426c5df90dSjeremylt       ierr = MatCreateShell(comm, lsize[i], lsize[i-1], gsize[i], gsize[i-1],
243a97643b0Sjeremylt                             userPR[i], &matPR[i]); CHKERRQ(ierr);
244a97643b0Sjeremylt       ierr = MatShellSetOperation(matPR[i], MATOP_MULT,
245a97643b0Sjeremylt                                   (void(*)(void))MatMult_Prolong);
2466c5df90dSjeremylt       CHKERRQ(ierr);
247a97643b0Sjeremylt       ierr = MatShellSetOperation(matPR[i], MATOP_MULT_TRANSPOSE,
2486c5df90dSjeremylt                                   (void(*)(void))MatMult_Restrict);
2496c5df90dSjeremylt       CHKERRQ(ierr);
2509396343dSjeremylt       if (memtyperequested == CEED_MEM_DEVICE) {
2519396343dSjeremylt         ierr = MatShellSetVecType(matPR[i], VECCUDA); CHKERRQ(ierr);
2529396343dSjeremylt       }
2536c5df90dSjeremylt     }
2546c5df90dSjeremylt   }
25561608365Sjeremylt   ierr = VecDuplicate(X[fineLevel], &rhs); CHKERRQ(ierr);
2566c5df90dSjeremylt 
2576c5df90dSjeremylt   // Print global grid information
2586c5df90dSjeremylt   if (!test_mode) {
2596c5df90dSjeremylt     PetscInt P = degree + 1, Q = P + qextra;
2609396343dSjeremylt 
2612d03409cSjeremylt     const char *usedresource;
2622d03409cSjeremylt     CeedGetResource(ceed, &usedresource);
2639396343dSjeremylt 
2649396343dSjeremylt     VecType vectype;
2659396343dSjeremylt     ierr = VecGetType(X[0], &vectype); CHKERRQ(ierr);
2669396343dSjeremylt 
2676c5df90dSjeremylt     ierr = PetscPrintf(comm,
2686c5df90dSjeremylt                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc + PCMG --\n"
2699396343dSjeremylt                        "  PETSc:\n"
2709396343dSjeremylt                        "    PETSc Vec Type                     : %s\n"
2716c5df90dSjeremylt                        "  libCEED:\n"
2726c5df90dSjeremylt                        "    libCEED Backend                    : %s\n"
2739396343dSjeremylt                        "    libCEED Backend MemType            : %s\n"
2749396343dSjeremylt                        "    libCEED User Requested MemType     : %s\n"
2756c5df90dSjeremylt                        "  Mesh:\n"
2766c5df90dSjeremylt                        "    Number of 1D Basis Nodes (p)       : %d\n"
2776c5df90dSjeremylt                        "    Number of 1D Quadrature Points (q) : %d\n"
2786c5df90dSjeremylt                        "    Global Nodes                       : %D\n"
2796c5df90dSjeremylt                        "    Owned Nodes                        : %D\n"
280db419314Sjeremylt                        "    DoF per node                       : %D\n"
2816c5df90dSjeremylt                        "  Multigrid:\n"
2826c5df90dSjeremylt                        "    Number of Levels                   : %d\n",
2839396343dSjeremylt                        bpchoice+1, vectype, usedresource,
2849396343dSjeremylt                        CeedMemTypes[memtypebackend],
2859396343dSjeremylt                        (setmemtyperequest) ?
2869396343dSjeremylt                        CeedMemTypes[memtyperequested] : "none",
2879396343dSjeremylt                        P, Q, gsize[fineLevel]/ncompu, lsize[fineLevel]/ncompu,
288db419314Sjeremylt                        ncompu, numlevels); CHKERRQ(ierr);
2896c5df90dSjeremylt   }
2906c5df90dSjeremylt 
2916c5df90dSjeremylt   // Create RHS vector
29261608365Sjeremylt   ierr = VecDuplicate(Xloc[fineLevel], &rhsloc); CHKERRQ(ierr);
2936c5df90dSjeremylt   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
2949396343dSjeremylt   if (memtyperequested == CEED_MEM_HOST) {
2956c5df90dSjeremylt     ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr);
2969396343dSjeremylt   } else {
2979396343dSjeremylt     ierr = VecCUDAGetArray(rhsloc, &r); CHKERRQ(ierr);
2989396343dSjeremylt   }
29961608365Sjeremylt   CeedVectorCreate(ceed, xlsize[fineLevel], &rhsceed);
300bcc01d93SJeremy L Thompson   CeedVectorSetArray(rhsceed, memtyperequested, CEED_USE_POINTER, r);
3016c5df90dSjeremylt 
3026c5df90dSjeremylt   // Set up libCEED operators on each level
3036c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &ceeddata); CHKERRQ(ierr);
3046c5df90dSjeremylt   for (int i=0; i<numlevels; i++) {
3056c5df90dSjeremylt     // Print level information
30661608365Sjeremylt     if (!test_mode && (i == 0 || i == fineLevel)) {
3076c5df90dSjeremylt       ierr = PetscPrintf(comm,"    Level %D (%s):\n"
3086c5df90dSjeremylt                          "      Number of 1D Basis Nodes (p)     : %d\n"
3096c5df90dSjeremylt                          "      Global Nodes                     : %D\n"
3106c5df90dSjeremylt                          "      Owned Nodes                      : %D\n",
3116c5df90dSjeremylt                          i, (i? "fine" : "coarse"), leveldegrees[i] + 1,
3126c5df90dSjeremylt                          gsize[i]/ncompu, lsize[i]/ncompu); CHKERRQ(ierr);
3136c5df90dSjeremylt     }
3146c5df90dSjeremylt     ierr = PetscMalloc1(1, &ceeddata[i]); CHKERRQ(ierr);
3156c5df90dSjeremylt     ierr = SetupLibceedByDegree(dm[i], ceed, leveldegrees[i], dim, qextra,
316199551f5Sjeremylt                                 ncompu, gsize[i], xlsize[i], bpchoice,
317*f9342789SJeremy L Thompson                                 ceeddata[i], i==(fineLevel), rhsceed, &target);
318*f9342789SJeremy L Thompson     CHKERRQ(ierr);
3196c5df90dSjeremylt   }
3206c5df90dSjeremylt 
3216c5df90dSjeremylt   // Gather RHS
3229396343dSjeremylt   CeedVectorSyncArray(rhsceed, memtyperequested);
3239396343dSjeremylt   if (memtyperequested == CEED_MEM_HOST) {
3246c5df90dSjeremylt     ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr);
3259396343dSjeremylt   } else {
3269396343dSjeremylt     ierr = VecCUDARestoreArray(rhsloc, &r); CHKERRQ(ierr);
3279396343dSjeremylt   }
3286c5df90dSjeremylt   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
3299396343dSjeremylt   ierr = DMLocalToGlobal(dm[fineLevel], rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr);
3306c5df90dSjeremylt   CeedVectorDestroy(&rhsceed);
3316c5df90dSjeremylt 
33243eb8658SJeremy L Thompson   // Create the restriction/interpolation QFunction
33360f77c51Sjeremylt   CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_NONE, CEED_EVAL_INTERP,
334c70bd2a0Sjeremylt                               &qfrestrict);
33560f77c51Sjeremylt   CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_INTERP, CEED_EVAL_NONE,
336c70bd2a0Sjeremylt                               &qfprolong);
3376c5df90dSjeremylt 
3386c5df90dSjeremylt   // Set up libCEED level transfer operators
339199551f5Sjeremylt   ierr = CeedLevelTransferSetup(ceed, numlevels, ncompu, bpchoice, ceeddata,
340c70bd2a0Sjeremylt                                 leveldegrees, qfrestrict, qfprolong);
3416c5df90dSjeremylt   CHKERRQ(ierr);
3426c5df90dSjeremylt 
34343eb8658SJeremy L Thompson   // Create the error QFunction
344199551f5Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].error,
345199551f5Sjeremylt                               bpOptions[bpchoice].errorfname, &qferror);
346c70bd2a0Sjeremylt   CeedQFunctionAddInput(qferror, "u", ncompu, CEED_EVAL_INTERP);
347c70bd2a0Sjeremylt   CeedQFunctionAddInput(qferror, "true_soln", ncompu, CEED_EVAL_NONE);
348c70bd2a0Sjeremylt   CeedQFunctionAddOutput(qferror, "error", ncompu, CEED_EVAL_NONE);
3496c5df90dSjeremylt 
3506c5df90dSjeremylt   // Create the error operator
351c70bd2a0Sjeremylt   CeedOperatorCreate(ceed, qferror, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
352c70bd2a0Sjeremylt                      &operror);
353c70bd2a0Sjeremylt   CeedOperatorSetField(operror, "u", ceeddata[fineLevel]->Erestrictu,
35461608365Sjeremylt                        ceeddata[fineLevel]->basisu, CEED_VECTOR_ACTIVE);
355c70bd2a0Sjeremylt   CeedOperatorSetField(operror, "true_soln", ceeddata[fineLevel]->Erestrictui,
356a8d32208Sjeremylt                        CEED_BASIS_COLLOCATED, target);
357c70bd2a0Sjeremylt   CeedOperatorSetField(operror, "error", ceeddata[fineLevel]->Erestrictui,
358a8d32208Sjeremylt                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
3596c5df90dSjeremylt 
3606c5df90dSjeremylt   // Calculate multiplicity
3616c5df90dSjeremylt   for (int i=0; i<numlevels; i++) {
3626c5df90dSjeremylt     PetscScalar *x;
3636c5df90dSjeremylt 
3646c5df90dSjeremylt     // CEED vector
365*f9342789SJeremy L Thompson     ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr);
3666c5df90dSjeremylt     ierr = VecGetArray(Xloc[i], &x); CHKERRQ(ierr);
3676c5df90dSjeremylt     CeedVectorSetArray(ceeddata[i]->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
3686c5df90dSjeremylt 
3696c5df90dSjeremylt     // Multiplicity
370a8d32208Sjeremylt     CeedElemRestrictionGetMultiplicity(ceeddata[i]->Erestrictu,
3716c5df90dSjeremylt                                        ceeddata[i]->xceed);
372*f9342789SJeremy L Thompson     CeedVectorSyncArray(ceeddata[i]->xceed, CEED_MEM_HOST);
3736c5df90dSjeremylt 
3746c5df90dSjeremylt     // Restore vector
3756c5df90dSjeremylt     ierr = VecRestoreArray(Xloc[i], &x); CHKERRQ(ierr);
3766c5df90dSjeremylt 
3776c5df90dSjeremylt     // Creat mult vector
3786c5df90dSjeremylt     ierr = VecDuplicate(Xloc[i], &mult[i]); CHKERRQ(ierr);
3796c5df90dSjeremylt 
3806c5df90dSjeremylt     // Local-to-global
3816c5df90dSjeremylt     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
382483f8b0dSjeremylt     ierr = DMLocalToGlobal(dm[i], Xloc[i], ADD_VALUES, X[i]);
3836c5df90dSjeremylt     CHKERRQ(ierr);
3846c5df90dSjeremylt     ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr);
3856c5df90dSjeremylt 
3866c5df90dSjeremylt     // Global-to-local
387483f8b0dSjeremylt     ierr = DMGlobalToLocal(dm[i], X[i], INSERT_VALUES, mult[i]);
3886c5df90dSjeremylt     CHKERRQ(ierr);
3896c5df90dSjeremylt     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
3906c5df90dSjeremylt 
3916c5df90dSjeremylt     // Multiplicity scaling
3926c5df90dSjeremylt     ierr = VecReciprocal(mult[i]);
3936c5df90dSjeremylt   }
3946c5df90dSjeremylt 
3956c5df90dSjeremylt   // Set up Mat
3966c5df90dSjeremylt   for (int i=0; i<numlevels; i++) {
397226c3a8fSjeremylt     // User Operator
3986c5df90dSjeremylt     userO[i]->comm = comm;
3996c5df90dSjeremylt     userO[i]->dm = dm[i];
4006c5df90dSjeremylt     userO[i]->Xloc = Xloc[i];
4016c5df90dSjeremylt     ierr = VecDuplicate(Xloc[i], &userO[i]->Yloc); CHKERRQ(ierr);
4026c5df90dSjeremylt     userO[i]->xceed = ceeddata[i]->xceed;
4036c5df90dSjeremylt     userO[i]->yceed = ceeddata[i]->yceed;
404c70bd2a0Sjeremylt     userO[i]->op = ceeddata[i]->opapply;
4056c5df90dSjeremylt     userO[i]->ceed = ceed;
4069396343dSjeremylt     userO[i]->memtype = memtyperequested;
4079396343dSjeremylt     if (memtyperequested == CEED_MEM_HOST) {
4089396343dSjeremylt       userO[i]->VecGetArray = VecGetArray;
4099396343dSjeremylt       userO[i]->VecGetArrayRead = VecGetArrayRead;
4109396343dSjeremylt       userO[i]->VecRestoreArray = VecRestoreArray;
4119396343dSjeremylt       userO[i]->VecRestoreArrayRead = VecRestoreArrayRead;
4129396343dSjeremylt     } else {
4139396343dSjeremylt       userO[i]->VecGetArray = VecCUDAGetArray;
4149396343dSjeremylt       userO[i]->VecGetArrayRead = VecCUDAGetArrayRead;
4159396343dSjeremylt       userO[i]->VecRestoreArray = VecCUDARestoreArray;
4169396343dSjeremylt       userO[i]->VecRestoreArrayRead = VecCUDARestoreArrayRead;
4179396343dSjeremylt     }
4186c5df90dSjeremylt 
4196c5df90dSjeremylt     if (i > 0) {
420a97643b0Sjeremylt       // Prolongation/Restriction Operator
421a97643b0Sjeremylt       userPR[i]->comm = comm;
422199551f5Sjeremylt       userPR[i]->dmf = dm[i];
423199551f5Sjeremylt       userPR[i]->dmc = dm[i-1];
424199551f5Sjeremylt       userPR[i]->locvecc = Xloc[i-1];
425199551f5Sjeremylt       userPR[i]->locvecf = userO[i]->Yloc;
426199551f5Sjeremylt       userPR[i]->multvec = mult[i];
427199551f5Sjeremylt       userPR[i]->ceedvecc = userO[i-1]->xceed;
428199551f5Sjeremylt       userPR[i]->ceedvecf = userO[i]->yceed;
429c70bd2a0Sjeremylt       userPR[i]->opprolong = ceeddata[i]->opprolong;
430c70bd2a0Sjeremylt       userPR[i]->oprestrict = ceeddata[i]->oprestrict;
431a97643b0Sjeremylt       userPR[i]->ceed = ceed;
4329396343dSjeremylt       userPR[i]->memtype = userO[i]->memtype;
4339396343dSjeremylt       userPR[i]->VecGetArray = userO[i]->VecGetArray;
4349396343dSjeremylt       userPR[i]->VecGetArrayRead = userO[i]->VecGetArrayRead;
4359396343dSjeremylt       userPR[i]->VecRestoreArray = userO[i]->VecRestoreArray;
4369396343dSjeremylt       userPR[i]->VecRestoreArrayRead = userO[i]->VecRestoreArrayRead;
4376c5df90dSjeremylt     }
4386c5df90dSjeremylt   }
4396c5df90dSjeremylt 
44015ce0ef0Sjeremylt   // Setup dummy SNES for AMG coarse solve
441c70bd2a0Sjeremylt   ierr = SNESCreate(comm, &snesdummy); CHKERRQ(ierr);
442c70bd2a0Sjeremylt   ierr = SNESSetDM(snesdummy, dm[0]); CHKERRQ(ierr);
443c70bd2a0Sjeremylt   ierr = SNESSetSolution(snesdummy, X[0]); CHKERRQ(ierr);
44415ce0ef0Sjeremylt 
44515ce0ef0Sjeremylt   // -- Jacobian matrix
44615ce0ef0Sjeremylt   ierr = DMSetMatType(dm[0], MATAIJ); CHKERRQ(ierr);
447199551f5Sjeremylt   ierr = DMCreateMatrix(dm[0], &matcoarse); CHKERRQ(ierr);
448199551f5Sjeremylt   ierr = SNESSetJacobian(snesdummy, matcoarse, matcoarse, NULL,
44915ce0ef0Sjeremylt                          NULL); CHKERRQ(ierr);
45015ce0ef0Sjeremylt 
45115ce0ef0Sjeremylt   // -- Residual evaluation function
452c70bd2a0Sjeremylt   ierr = SNESSetFunction(snesdummy, X[0], FormResidual_Ceed,
45315ce0ef0Sjeremylt                          userO[0]); CHKERRQ(ierr);
45415ce0ef0Sjeremylt 
45515ce0ef0Sjeremylt   // -- Form Jacobian
456c70bd2a0Sjeremylt   ierr = SNESComputeJacobianDefaultColor(snesdummy, X[0], matO[0],
457199551f5Sjeremylt                                          matcoarse, NULL); CHKERRQ(ierr);
45815ce0ef0Sjeremylt 
4596c5df90dSjeremylt   // Set up KSP
4606c5df90dSjeremylt   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
4616c5df90dSjeremylt   {
4626c5df90dSjeremylt     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
4636c5df90dSjeremylt     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
4646c5df90dSjeremylt     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
4656c5df90dSjeremylt                             PETSC_DEFAULT); CHKERRQ(ierr);
4666c5df90dSjeremylt   }
4676c5df90dSjeremylt   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
46861608365Sjeremylt   ierr = KSPSetOperators(ksp, matO[fineLevel], matO[fineLevel]);
4696c5df90dSjeremylt   CHKERRQ(ierr);
4706c5df90dSjeremylt 
4716c5df90dSjeremylt   // Set up PCMG
4726c5df90dSjeremylt   ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
4736c5df90dSjeremylt   PCMGCycleType pcgmcycletype = PC_MG_CYCLE_V;
4746c5df90dSjeremylt   {
4756c5df90dSjeremylt     ierr = PCSetType(pc, PCMG); CHKERRQ(ierr);
4766c5df90dSjeremylt 
4776c5df90dSjeremylt     // PCMG levels
4786c5df90dSjeremylt     ierr = PCMGSetLevels(pc, numlevels, NULL); CHKERRQ(ierr);
4796c5df90dSjeremylt     for (int i=0; i<numlevels; i++) {
4806c5df90dSjeremylt       // Smoother
4816c5df90dSjeremylt       KSP smoother;
4826c5df90dSjeremylt       PC smoother_pc;
4836c5df90dSjeremylt       ierr = PCMGGetSmoother(pc, i, &smoother); CHKERRQ(ierr);
4846c5df90dSjeremylt       ierr = KSPSetType(smoother, KSPCHEBYSHEV); CHKERRQ(ierr);
4856c5df90dSjeremylt       ierr = KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1); CHKERRQ(ierr);
4866c5df90dSjeremylt       ierr = KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE); CHKERRQ(ierr);
4876c5df90dSjeremylt       ierr = KSPSetOperators(smoother, matO[i], matO[i]); CHKERRQ(ierr);
4886c5df90dSjeremylt       ierr = KSPGetPC(smoother, &smoother_pc); CHKERRQ(ierr);
4896c5df90dSjeremylt       ierr = PCSetType(smoother_pc, PCJACOBI); CHKERRQ(ierr);
4906c5df90dSjeremylt       ierr = PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
4916c5df90dSjeremylt 
4926c5df90dSjeremylt       // Work vector
4936c5df90dSjeremylt       if (i < numlevels - 1) {
4946c5df90dSjeremylt         ierr = PCMGSetX(pc, i, X[i]); CHKERRQ(ierr);
4956c5df90dSjeremylt       }
4966c5df90dSjeremylt 
4976c5df90dSjeremylt       // Level transfers
4986c5df90dSjeremylt       if (i > 0) {
4996c5df90dSjeremylt         // Interpolation
500a97643b0Sjeremylt         ierr = PCMGSetInterpolation(pc, i, matPR[i]); CHKERRQ(ierr);
5016c5df90dSjeremylt       }
5026c5df90dSjeremylt 
5036c5df90dSjeremylt       // Coarse solve
5046c5df90dSjeremylt       KSP coarse;
5056c5df90dSjeremylt       PC coarse_pc;
5066c5df90dSjeremylt       ierr = PCMGGetCoarseSolve(pc, &coarse); CHKERRQ(ierr);
50715ce0ef0Sjeremylt       ierr = KSPSetType(coarse, KSPPREONLY); CHKERRQ(ierr);
508199551f5Sjeremylt       ierr = KSPSetOperators(coarse, matcoarse, matcoarse); CHKERRQ(ierr);
50915ce0ef0Sjeremylt 
5106c5df90dSjeremylt       ierr = KSPGetPC(coarse, &coarse_pc); CHKERRQ(ierr);
51115ce0ef0Sjeremylt       ierr = PCSetType(coarse_pc, PCGAMG); CHKERRQ(ierr);
51215ce0ef0Sjeremylt 
51315ce0ef0Sjeremylt       ierr = KSPSetOptionsPrefix(coarse, "coarse_"); CHKERRQ(ierr);
51415ce0ef0Sjeremylt       ierr = PCSetOptionsPrefix(coarse_pc, "coarse_"); CHKERRQ(ierr);
51515ce0ef0Sjeremylt       ierr = KSPSetFromOptions(coarse); CHKERRQ(ierr);
51615ce0ef0Sjeremylt       ierr = PCSetFromOptions(coarse_pc); CHKERRQ(ierr);
5176c5df90dSjeremylt     }
5186c5df90dSjeremylt 
5196c5df90dSjeremylt     // PCMG options
5206c5df90dSjeremylt     ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr);
5216c5df90dSjeremylt     ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr);
5226c5df90dSjeremylt     ierr = PCMGSetCycleType(pc, pcgmcycletype); CHKERRQ(ierr);
5236c5df90dSjeremylt   }
5246c5df90dSjeremylt 
5256c5df90dSjeremylt   // First run, if benchmarking
5266c5df90dSjeremylt   if (benchmark_mode) {
5276c5df90dSjeremylt     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
5286c5df90dSjeremylt     CHKERRQ(ierr);
52961608365Sjeremylt     ierr = VecZeroEntries(X[fineLevel]); CHKERRQ(ierr);
5306c5df90dSjeremylt     my_rt_start = MPI_Wtime();
53161608365Sjeremylt     ierr = KSPSolve(ksp, rhs, X[fineLevel]); CHKERRQ(ierr);
5326c5df90dSjeremylt     my_rt = MPI_Wtime() - my_rt_start;
5336c5df90dSjeremylt     ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
5346c5df90dSjeremylt     CHKERRQ(ierr);
5356c5df90dSjeremylt     // Set maxits based on first iteration timing
5366c5df90dSjeremylt     if (my_rt > 0.02) {
5376c5df90dSjeremylt       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
5386c5df90dSjeremylt       CHKERRQ(ierr);
5396c5df90dSjeremylt     } else {
5406c5df90dSjeremylt       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
5416c5df90dSjeremylt       CHKERRQ(ierr);
5426c5df90dSjeremylt     }
5436c5df90dSjeremylt   }
5446c5df90dSjeremylt 
5456c5df90dSjeremylt   // Timed solve
54661608365Sjeremylt   ierr = VecZeroEntries(X[fineLevel]); CHKERRQ(ierr);
5476c5df90dSjeremylt   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
54809a940d7Sjeremylt 
54909a940d7Sjeremylt   // -- Performance logging
55009a940d7Sjeremylt   ierr = PetscLogStageRegister("Solve Stage", &solvestage); CHKERRQ(ierr);
55109a940d7Sjeremylt   ierr = PetscLogStagePush(solvestage); CHKERRQ(ierr);
55209a940d7Sjeremylt 
55309a940d7Sjeremylt   // -- Solve
5546c5df90dSjeremylt   my_rt_start = MPI_Wtime();
55561608365Sjeremylt   ierr = KSPSolve(ksp, rhs, X[fineLevel]); CHKERRQ(ierr);
5566c5df90dSjeremylt   my_rt = MPI_Wtime() - my_rt_start;
5576c5df90dSjeremylt 
55809a940d7Sjeremylt 
55909a940d7Sjeremylt   // -- Performance logging
56009a940d7Sjeremylt   ierr = PetscLogStagePop();
56109a940d7Sjeremylt 
5626c5df90dSjeremylt   // Output results
5636c5df90dSjeremylt   {
5646c5df90dSjeremylt     KSPType ksptype;
5656c5df90dSjeremylt     PCMGType pcmgtype;
5666c5df90dSjeremylt     KSPConvergedReason reason;
5676c5df90dSjeremylt     PetscReal rnorm;
5686c5df90dSjeremylt     PetscInt its;
5696c5df90dSjeremylt     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
5706c5df90dSjeremylt     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
5716c5df90dSjeremylt     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
5726c5df90dSjeremylt     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
5736c5df90dSjeremylt     ierr = PCMGGetType(pc, &pcmgtype); CHKERRQ(ierr);
5746c5df90dSjeremylt     if (!test_mode || reason < 0 || rnorm > 1e-8) {
5756c5df90dSjeremylt       ierr = PetscPrintf(comm,
5766c5df90dSjeremylt                          "  KSP:\n"
5776c5df90dSjeremylt                          "    KSP Type                           : %s\n"
5786c5df90dSjeremylt                          "    KSP Convergence                    : %s\n"
5796c5df90dSjeremylt                          "    Total KSP Iterations               : %D\n"
5806c5df90dSjeremylt                          "    Final rnorm                        : %e\n",
5816c5df90dSjeremylt                          ksptype, KSPConvergedReasons[reason], its,
5826c5df90dSjeremylt                          (double)rnorm); CHKERRQ(ierr);
5836c5df90dSjeremylt       ierr = PetscPrintf(comm,
5846c5df90dSjeremylt                          "  PCMG:\n"
5856c5df90dSjeremylt                          "    PCMG Type                          : %s\n"
5866c5df90dSjeremylt                          "    PCMG Cycle Type                    : %s\n",
5876c5df90dSjeremylt                          PCMGTypes[pcmgtype],
5886c5df90dSjeremylt                          PCMGCycleTypes[pcgmcycletype]); CHKERRQ(ierr);
5896c5df90dSjeremylt     }
5906c5df90dSjeremylt     if (!test_mode) {
5916c5df90dSjeremylt       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
5926c5df90dSjeremylt     }
5936c5df90dSjeremylt     {
5946c5df90dSjeremylt       PetscReal maxerror;
595c70bd2a0Sjeremylt       ierr = ComputeErrorMax(userO[fineLevel], operror, X[fineLevel], target,
5966c5df90dSjeremylt                              &maxerror); CHKERRQ(ierr);
5976c5df90dSjeremylt       PetscReal tol = 5e-2;
5986c5df90dSjeremylt       if (!test_mode || maxerror > tol) {
5996c5df90dSjeremylt         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
6006c5df90dSjeremylt         CHKERRQ(ierr);
6016c5df90dSjeremylt         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
6026c5df90dSjeremylt         CHKERRQ(ierr);
6036c5df90dSjeremylt         ierr = PetscPrintf(comm,
6046c5df90dSjeremylt                            "    Pointwise Error (max)              : %e\n"
6056c5df90dSjeremylt                            "    CG Solve Time                      : %g (%g) sec\n",
6066c5df90dSjeremylt                            (double)maxerror, rt_max, rt_min); CHKERRQ(ierr);
6076c5df90dSjeremylt       }
6086c5df90dSjeremylt     }
6096c5df90dSjeremylt     if (benchmark_mode && (!test_mode)) {
6106c5df90dSjeremylt       ierr = PetscPrintf(comm,
6116c5df90dSjeremylt                          "    DoFs/Sec in CG                     : %g (%g) million\n",
61261608365Sjeremylt                          1e-6*gsize[fineLevel]*its/rt_max,
61361608365Sjeremylt                          1e-6*gsize[fineLevel]*its/rt_min);
6146c5df90dSjeremylt       CHKERRQ(ierr);
6156c5df90dSjeremylt     }
6166c5df90dSjeremylt   }
6176c5df90dSjeremylt 
6186c5df90dSjeremylt   if (write_solution) {
6196c5df90dSjeremylt     PetscViewer vtkviewersoln;
6206c5df90dSjeremylt 
6216c5df90dSjeremylt     ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr);
6226c5df90dSjeremylt     ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr);
6236c5df90dSjeremylt     ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr);
62461608365Sjeremylt     ierr = VecView(X[fineLevel], vtkviewersoln); CHKERRQ(ierr);
6256c5df90dSjeremylt     ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr);
6266c5df90dSjeremylt   }
6276c5df90dSjeremylt 
6286c5df90dSjeremylt   // Cleanup
6296c5df90dSjeremylt   for (int i=0; i<numlevels; i++) {
6306c5df90dSjeremylt     ierr = VecDestroy(&X[i]); CHKERRQ(ierr);
6316c5df90dSjeremylt     ierr = VecDestroy(&Xloc[i]); CHKERRQ(ierr);
6326c5df90dSjeremylt     ierr = VecDestroy(&mult[i]); CHKERRQ(ierr);
6336c5df90dSjeremylt     ierr = VecDestroy(&userO[i]->Yloc); CHKERRQ(ierr);
6346c5df90dSjeremylt     ierr = MatDestroy(&matO[i]); CHKERRQ(ierr);
6356c5df90dSjeremylt     ierr = PetscFree(userO[i]); CHKERRQ(ierr);
6366c5df90dSjeremylt     if (i > 0) {
637a97643b0Sjeremylt       ierr = MatDestroy(&matPR[i]); CHKERRQ(ierr);
638a97643b0Sjeremylt       ierr = PetscFree(userPR[i]); CHKERRQ(ierr);
6396c5df90dSjeremylt     }
6406c5df90dSjeremylt     ierr = CeedDataDestroy(i, ceeddata[i]); CHKERRQ(ierr);
6416c5df90dSjeremylt     ierr = DMDestroy(&dm[i]); CHKERRQ(ierr);
6426c5df90dSjeremylt   }
6436c5df90dSjeremylt   ierr = PetscFree(leveldegrees); CHKERRQ(ierr);
6446c5df90dSjeremylt   ierr = PetscFree(dm); CHKERRQ(ierr);
6456c5df90dSjeremylt   ierr = PetscFree(X); CHKERRQ(ierr);
6466c5df90dSjeremylt   ierr = PetscFree(Xloc); CHKERRQ(ierr);
6476c5df90dSjeremylt   ierr = PetscFree(mult); CHKERRQ(ierr);
6486c5df90dSjeremylt   ierr = PetscFree(matO); CHKERRQ(ierr);
649a97643b0Sjeremylt   ierr = PetscFree(matPR); CHKERRQ(ierr);
6506c5df90dSjeremylt   ierr = PetscFree(ceeddata); CHKERRQ(ierr);
6516c5df90dSjeremylt   ierr = PetscFree(userO); CHKERRQ(ierr);
652a97643b0Sjeremylt   ierr = PetscFree(userPR); CHKERRQ(ierr);
6536c5df90dSjeremylt   ierr = PetscFree(lsize); CHKERRQ(ierr);
6546c5df90dSjeremylt   ierr = PetscFree(xlsize); CHKERRQ(ierr);
6556c5df90dSjeremylt   ierr = PetscFree(gsize); CHKERRQ(ierr);
6566c5df90dSjeremylt   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
6576c5df90dSjeremylt   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
658199551f5Sjeremylt   ierr = MatDestroy(&matcoarse); CHKERRQ(ierr);
6596c5df90dSjeremylt   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
660c70bd2a0Sjeremylt   ierr = SNESDestroy(&snesdummy); CHKERRQ(ierr);
661199551f5Sjeremylt   ierr = DMDestroy(&dmorig); CHKERRQ(ierr);
6626c5df90dSjeremylt   CeedVectorDestroy(&target);
663c70bd2a0Sjeremylt   CeedQFunctionDestroy(&qferror);
664c70bd2a0Sjeremylt   CeedQFunctionDestroy(&qfrestrict);
665c70bd2a0Sjeremylt   CeedQFunctionDestroy(&qfprolong);
666c70bd2a0Sjeremylt   CeedOperatorDestroy(&operror);
6676c5df90dSjeremylt   CeedDestroy(&ceed);
6686c5df90dSjeremylt   return PetscFinalize();
6696c5df90dSjeremylt }
670