xref: /libCEED/examples/petsc/multigrid.c (revision 2d03409ca85173e0565dadbac5c9383b99bac8d1)
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;
476c5df90dSjeremylt   char ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self",
486c5df90dSjeremylt        filename[PETSC_MAX_PATH_LEN];
496c5df90dSjeremylt   double my_rt_start, my_rt, rt_min, rt_max;
506c5df90dSjeremylt   PetscInt degree = 3, qextra, *lsize, *xlsize, *gsize, dim = 3,
516c5df90dSjeremylt            melem[3] = {3, 3, 3}, ncompu = 1, numlevels = degree, *leveldegrees;
526c5df90dSjeremylt   PetscScalar *r;
536c5df90dSjeremylt   PetscBool test_mode, benchmark_mode, read_mesh, write_solution;
546c5df90dSjeremylt   DM  *dm, dmOrig;
556c5df90dSjeremylt   KSP ksp;
566c5df90dSjeremylt   PC pc;
576c5df90dSjeremylt   Mat *matO, *matI, *matR;
586c5df90dSjeremylt   Vec *X, *Xloc, *mult, rhs, rhsloc, diagloc;
596c5df90dSjeremylt   UserO *userO;
606c5df90dSjeremylt   UserIR *userI, *userR;
616c5df90dSjeremylt   Ceed ceed;
626c5df90dSjeremylt   CeedData *ceeddata;
636c5df90dSjeremylt   CeedVector rhsceed, diagceed, target;
646c5df90dSjeremylt   CeedQFunction qf_error, qf_restrict, qf_prolong;
656c5df90dSjeremylt   CeedOperator op_error;
666c5df90dSjeremylt   bpType bpChoice;
676c5df90dSjeremylt   coarsenType coarsen;
686c5df90dSjeremylt 
696c5df90dSjeremylt   ierr = PetscInitialize(&argc, &argv, NULL, help);
706c5df90dSjeremylt   if (ierr) return ierr;
716c5df90dSjeremylt   comm = PETSC_COMM_WORLD;
726c5df90dSjeremylt 
736c5df90dSjeremylt   // Parse command line options
746c5df90dSjeremylt   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr);
756c5df90dSjeremylt   bpChoice = CEED_BP3;
766c5df90dSjeremylt   ierr = PetscOptionsEnum("-problem",
776c5df90dSjeremylt                           "CEED benchmark problem to solve", NULL,
786c5df90dSjeremylt                           bpTypes, (PetscEnum)bpChoice, (PetscEnum *)&bpChoice,
796c5df90dSjeremylt                           NULL); CHKERRQ(ierr);
806c5df90dSjeremylt   ncompu = bpOptions[bpChoice].ncompu;
816c5df90dSjeremylt   test_mode = PETSC_FALSE;
826c5df90dSjeremylt   ierr = PetscOptionsBool("-test",
836c5df90dSjeremylt                           "Testing mode (do not print unless error is large)",
846c5df90dSjeremylt                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
856c5df90dSjeremylt   benchmark_mode = PETSC_FALSE;
866c5df90dSjeremylt   ierr = PetscOptionsBool("-benchmark",
876c5df90dSjeremylt                           "Benchmarking mode (prints benchmark statistics)",
886c5df90dSjeremylt                           NULL, benchmark_mode, &benchmark_mode, NULL);
896c5df90dSjeremylt   CHKERRQ(ierr);
906c5df90dSjeremylt   write_solution = PETSC_FALSE;
916c5df90dSjeremylt   ierr = PetscOptionsBool("-write_solution",
926c5df90dSjeremylt                           "Write solution for visualization",
936c5df90dSjeremylt                           NULL, write_solution, &write_solution, NULL);
946c5df90dSjeremylt   CHKERRQ(ierr);
956c5df90dSjeremylt   degree = test_mode ? 3 : 2;
966c5df90dSjeremylt   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
976c5df90dSjeremylt                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
986c5df90dSjeremylt   if (degree < 1) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
996c5df90dSjeremylt                              "-degree %D must be at least 1", degree);
1006c5df90dSjeremylt   qextra = bpOptions[bpChoice].qextra;
1016c5df90dSjeremylt   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
1026c5df90dSjeremylt                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
1036c5df90dSjeremylt   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
1046c5df90dSjeremylt                             NULL, ceedresource, ceedresource,
1056c5df90dSjeremylt                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
1066c5df90dSjeremylt   coarsen = COARSEN_UNIFORM;
1076c5df90dSjeremylt   ierr = PetscOptionsEnum("-coarsen",
1086c5df90dSjeremylt                           "Coarsening strategy to use", NULL,
1096c5df90dSjeremylt                           coarsenTypes, (PetscEnum)coarsen,
1106c5df90dSjeremylt                           (PetscEnum *)&coarsen, NULL); CHKERRQ(ierr);
1116c5df90dSjeremylt   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL,
1126c5df90dSjeremylt                             filename, filename, sizeof(filename), &read_mesh);
1136c5df90dSjeremylt   CHKERRQ(ierr);
1146c5df90dSjeremylt   if (!read_mesh) {
1156c5df90dSjeremylt     PetscInt tmp = dim;
1166c5df90dSjeremylt     ierr = PetscOptionsIntArray("-cells","Number of cells per dimension", NULL,
1176c5df90dSjeremylt                                 melem, &tmp, NULL); CHKERRQ(ierr);
1186c5df90dSjeremylt   }
1196c5df90dSjeremylt   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
1206c5df90dSjeremylt 
1216c5df90dSjeremylt   // Setup DM
1226c5df90dSjeremylt   if (read_mesh) {
1236c5df90dSjeremylt     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dmOrig);
1246c5df90dSjeremylt     CHKERRQ(ierr);
1256c5df90dSjeremylt   } else {
1266c5df90dSjeremylt     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_FALSE, melem, NULL,
1276c5df90dSjeremylt                                NULL, NULL, PETSC_TRUE,&dmOrig); CHKERRQ(ierr);
1286c5df90dSjeremylt   }
1296c5df90dSjeremylt 
1306c5df90dSjeremylt   {
1316c5df90dSjeremylt     DM dmDist = NULL;
1326c5df90dSjeremylt     PetscPartitioner part;
1336c5df90dSjeremylt 
1346c5df90dSjeremylt     ierr = DMPlexGetPartitioner(dmOrig, &part); CHKERRQ(ierr);
1356c5df90dSjeremylt     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
1366c5df90dSjeremylt     ierr = DMPlexDistribute(dmOrig, 0, NULL, &dmDist); CHKERRQ(ierr);
1376c5df90dSjeremylt     if (dmDist) {
1386c5df90dSjeremylt       ierr = DMDestroy(&dmOrig); CHKERRQ(ierr);
1396c5df90dSjeremylt       dmOrig = dmDist;
1406c5df90dSjeremylt     }
1416c5df90dSjeremylt   }
1426c5df90dSjeremylt 
1436c5df90dSjeremylt   // Allocate arrays for PETSc objects for each level
1446c5df90dSjeremylt   switch (coarsen) {
1456c5df90dSjeremylt   case COARSEN_UNIFORM:
1466c5df90dSjeremylt     numlevels = degree;
1476c5df90dSjeremylt     break;
1486c5df90dSjeremylt   case COARSEN_LOGRITHMIC:
1496c5df90dSjeremylt     numlevels = ceil(log(degree)/log(2)) + 1;
1506c5df90dSjeremylt     break;
1516c5df90dSjeremylt   }
1526c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &leveldegrees); CHKERRQ(ierr);
1536c5df90dSjeremylt   switch (coarsen) {
1546c5df90dSjeremylt   case COARSEN_UNIFORM:
1556c5df90dSjeremylt     for (int i=0; i<numlevels; i++) leveldegrees[i] = i + 1;
1566c5df90dSjeremylt     break;
1576c5df90dSjeremylt   case COARSEN_LOGRITHMIC:
1586c5df90dSjeremylt     for (int i=0; i<numlevels-1; i++) leveldegrees[i] = pow(2,i);
1596c5df90dSjeremylt     leveldegrees[numlevels-1] = degree;
1606c5df90dSjeremylt     break;
1616c5df90dSjeremylt   }
1626c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &dm); CHKERRQ(ierr);
1636c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &X); CHKERRQ(ierr);
1646c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &Xloc); CHKERRQ(ierr);
1656c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &mult); CHKERRQ(ierr);
1666c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &userO); CHKERRQ(ierr);
1676c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &userI); CHKERRQ(ierr);
1686c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &userR); CHKERRQ(ierr);
1696c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &matO); CHKERRQ(ierr);
1706c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &matI); CHKERRQ(ierr);
1716c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &matR); CHKERRQ(ierr);
1726c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &lsize); CHKERRQ(ierr);
1736c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &xlsize); CHKERRQ(ierr);
1746c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &gsize); CHKERRQ(ierr);
1756c5df90dSjeremylt 
1766c5df90dSjeremylt   // Setup DM and Operator Mat Shells for each level
1776c5df90dSjeremylt   for (CeedInt i=0; i<numlevels; i++) {
1786c5df90dSjeremylt     // Create DM
1796c5df90dSjeremylt     ierr = DMClone(dmOrig, &dm[i]); CHKERRQ(ierr);
1806c5df90dSjeremylt     ierr = SetupDMByDegree(dm[i], leveldegrees[i], ncompu, bpChoice);
1816c5df90dSjeremylt     CHKERRQ(ierr);
1826c5df90dSjeremylt 
1836c5df90dSjeremylt     // Create vectors
1846c5df90dSjeremylt     ierr = DMCreateGlobalVector(dm[i], &X[i]); CHKERRQ(ierr);
1856c5df90dSjeremylt     ierr = VecGetLocalSize(X[i], &lsize[i]); CHKERRQ(ierr);
1866c5df90dSjeremylt     ierr = VecGetSize(X[i], &gsize[i]); CHKERRQ(ierr);
1876c5df90dSjeremylt     ierr = DMCreateLocalVector(dm[i], &Xloc[i]); CHKERRQ(ierr);
1886c5df90dSjeremylt     ierr = VecGetSize(Xloc[i], &xlsize[i]); CHKERRQ(ierr);
1896c5df90dSjeremylt 
1906c5df90dSjeremylt     // Operator
1916c5df90dSjeremylt     ierr = PetscMalloc1(1, &userO[i]); CHKERRQ(ierr);
1926c5df90dSjeremylt     ierr = MatCreateShell(comm, lsize[i], lsize[i], gsize[i], gsize[i],
1936c5df90dSjeremylt                           userO[i], &matO[i]); CHKERRQ(ierr);
1946c5df90dSjeremylt     ierr = MatShellSetOperation(matO[i], MATOP_MULT,
1956c5df90dSjeremylt                                 (void(*)(void))MatMult_Ceed);
1966c5df90dSjeremylt     ierr = MatShellSetOperation(matO[i], MATOP_GET_DIAGONAL,
1976c5df90dSjeremylt                                 (void(*)(void))MatGetDiag);
1986c5df90dSjeremylt     CHKERRQ(ierr);
1996c5df90dSjeremylt 
2006c5df90dSjeremylt     // Level transfers
2016c5df90dSjeremylt     if (i > 0) {
2026c5df90dSjeremylt       // Interp
2036c5df90dSjeremylt       ierr = PetscMalloc1(1, &userI[i]); CHKERRQ(ierr);
2046c5df90dSjeremylt       ierr = MatCreateShell(comm, lsize[i], lsize[i-1], gsize[i], gsize[i-1],
2056c5df90dSjeremylt                             userI[i], &matI[i]); CHKERRQ(ierr);
2066c5df90dSjeremylt       ierr = MatShellSetOperation(matI[i], MATOP_MULT,
2076c5df90dSjeremylt                                   (void(*)(void))MatMult_Interp);
2086c5df90dSjeremylt       CHKERRQ(ierr);
2096c5df90dSjeremylt 
2106c5df90dSjeremylt       // Restrict
2116c5df90dSjeremylt       ierr = PetscMalloc1(1, &userR[i]); CHKERRQ(ierr);
2126c5df90dSjeremylt       ierr = MatCreateShell(comm, lsize[i-1], lsize[i], gsize[i-1], gsize[i],
2136c5df90dSjeremylt                             userR[i], &matR[i]); CHKERRQ(ierr);
2146c5df90dSjeremylt       ierr = MatShellSetOperation(matR[i], MATOP_MULT,
2156c5df90dSjeremylt                                   (void(*)(void))MatMult_Restrict);
2166c5df90dSjeremylt       CHKERRQ(ierr);
2176c5df90dSjeremylt     }
2186c5df90dSjeremylt   }
2196c5df90dSjeremylt   ierr = VecDuplicate(X[numlevels-1], &rhs); CHKERRQ(ierr);
2206c5df90dSjeremylt 
221*2d03409cSjeremylt   // Set up libCEED
222*2d03409cSjeremylt   CeedInit(ceedresource, &ceed);
223*2d03409cSjeremylt 
2246c5df90dSjeremylt   // Print global grid information
2256c5df90dSjeremylt   if (!test_mode) {
2266c5df90dSjeremylt     PetscInt P = degree + 1, Q = P + qextra;
227*2d03409cSjeremylt     const char *usedresource;
228*2d03409cSjeremylt     CeedGetResource(ceed, &usedresource);
2296c5df90dSjeremylt     ierr = PetscPrintf(comm,
2306c5df90dSjeremylt                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc + PCMG --\n"
2316c5df90dSjeremylt                        "  libCEED:\n"
2326c5df90dSjeremylt                        "    libCEED Backend                    : %s\n"
2336c5df90dSjeremylt                        "  Mesh:\n"
2346c5df90dSjeremylt                        "    Number of 1D Basis Nodes (p)       : %d\n"
2356c5df90dSjeremylt                        "    Number of 1D Quadrature Points (q) : %d\n"
2366c5df90dSjeremylt                        "    Global Nodes                       : %D\n"
2376c5df90dSjeremylt                        "    Owned Nodes                        : %D\n"
2386c5df90dSjeremylt                        "  Multigrid:\n"
2396c5df90dSjeremylt                        "    Number of Levels                   : %d\n",
240*2d03409cSjeremylt                        bpChoice+1, usedresource, P, Q,
2416c5df90dSjeremylt                        gsize[numlevels-1]/ncompu, lsize[numlevels-1]/ncompu,
2426c5df90dSjeremylt                        numlevels); CHKERRQ(ierr);
2436c5df90dSjeremylt   }
2446c5df90dSjeremylt 
2456c5df90dSjeremylt   // Create RHS vector
2466c5df90dSjeremylt   ierr = VecDuplicate(Xloc[numlevels-1], &rhsloc); CHKERRQ(ierr);
2476c5df90dSjeremylt   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
2486c5df90dSjeremylt   ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr);
2496c5df90dSjeremylt   CeedVectorCreate(ceed, xlsize[numlevels-1], &rhsceed);
2506c5df90dSjeremylt   CeedVectorSetArray(rhsceed, CEED_MEM_HOST, CEED_USE_POINTER, r);
2516c5df90dSjeremylt 
2526c5df90dSjeremylt   // Set up libCEED operators on each level
2536c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &ceeddata); CHKERRQ(ierr);
2546c5df90dSjeremylt   for (int i=0; i<numlevels; i++) {
2556c5df90dSjeremylt     // Print level information
2566c5df90dSjeremylt     if (!test_mode && (i == 0 || i == numlevels-1)) {
2576c5df90dSjeremylt       ierr = PetscPrintf(comm,"    Level %D (%s):\n"
2586c5df90dSjeremylt                          "      Number of 1D Basis Nodes (p)     : %d\n"
2596c5df90dSjeremylt                          "      Global Nodes                     : %D\n"
2606c5df90dSjeremylt                          "      Owned Nodes                      : %D\n",
2616c5df90dSjeremylt                          i, (i? "fine" : "coarse"), leveldegrees[i] + 1,
2626c5df90dSjeremylt                          gsize[i]/ncompu, lsize[i]/ncompu); CHKERRQ(ierr);
2636c5df90dSjeremylt     }
2646c5df90dSjeremylt     ierr = PetscMalloc1(1, &ceeddata[i]); CHKERRQ(ierr);
2656c5df90dSjeremylt     ierr = SetupLibceedByDegree(dm[i], ceed, leveldegrees[i], dim, qextra,
2667f823360Sjeremylt                                 ncompu, gsize[i], xlsize[i], bpChoice,
2677f823360Sjeremylt                                 ceeddata[i], i==(numlevels-1), rhsceed,
2687f823360Sjeremylt                                 &target); CHKERRQ(ierr);
2696c5df90dSjeremylt   }
2706c5df90dSjeremylt 
2716c5df90dSjeremylt   // Gather RHS
2726c5df90dSjeremylt   ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr);
2736c5df90dSjeremylt   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
2746c5df90dSjeremylt   ierr = DMLocalToGlobalBegin(dm[numlevels-1], rhsloc, ADD_VALUES, rhs);
2756c5df90dSjeremylt   CHKERRQ(ierr);
2766c5df90dSjeremylt   ierr = DMLocalToGlobalEnd(dm[numlevels-1], rhsloc, ADD_VALUES, rhs);
2776c5df90dSjeremylt   CHKERRQ(ierr);
2786c5df90dSjeremylt   CeedVectorDestroy(&rhsceed);
2796c5df90dSjeremylt 
2806c5df90dSjeremylt   // Create the restriction/interpolation Q-function
2816c5df90dSjeremylt   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].ident,
2826c5df90dSjeremylt                               bpOptions[bpChoice].identfname, &qf_restrict);
2836c5df90dSjeremylt   CeedQFunctionAddInput(qf_restrict, "uin", ncompu, CEED_EVAL_NONE);
2846c5df90dSjeremylt   CeedQFunctionAddOutput(qf_restrict, "uout", ncompu, CEED_EVAL_INTERP);
2856c5df90dSjeremylt   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].ident,
2866c5df90dSjeremylt                               bpOptions[bpChoice].identfname, &qf_prolong);
2876c5df90dSjeremylt   CeedQFunctionAddInput(qf_prolong, "uin", ncompu, CEED_EVAL_INTERP);
2886c5df90dSjeremylt   CeedQFunctionAddOutput(qf_prolong, "uout", ncompu, CEED_EVAL_NONE);
2896c5df90dSjeremylt 
2906c5df90dSjeremylt   // Set up libCEED level transfer operators
2916c5df90dSjeremylt   ierr = CeedLevelTransferSetup(ceed, numlevels, ncompu, bpChoice, ceeddata,
2926c5df90dSjeremylt                                 leveldegrees, qf_restrict, qf_prolong);
2936c5df90dSjeremylt   CHKERRQ(ierr);
2946c5df90dSjeremylt 
2956c5df90dSjeremylt   // Create the error Q-function
2966c5df90dSjeremylt   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].error,
2976c5df90dSjeremylt                               bpOptions[bpChoice].errorfname, &qf_error);
2986c5df90dSjeremylt   CeedQFunctionAddInput(qf_error, "u", ncompu, CEED_EVAL_INTERP);
2996c5df90dSjeremylt   CeedQFunctionAddInput(qf_error, "true_soln", ncompu, CEED_EVAL_NONE);
3006c5df90dSjeremylt   CeedQFunctionAddOutput(qf_error, "error", ncompu, CEED_EVAL_NONE);
3016c5df90dSjeremylt 
3026c5df90dSjeremylt   // Create the error operator
3036c5df90dSjeremylt   CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error);
3046c5df90dSjeremylt   CeedOperatorSetField(op_error, "u", ceeddata[numlevels-1]->Erestrictu,
3056c5df90dSjeremylt                        CEED_TRANSPOSE, ceeddata[numlevels-1]->basisu,
3066c5df90dSjeremylt                        CEED_VECTOR_ACTIVE);
3076c5df90dSjeremylt   CeedOperatorSetField(op_error, "true_soln", ceeddata[numlevels-1]->Erestrictui,
3086c5df90dSjeremylt                        CEED_NOTRANSPOSE, CEED_BASIS_COLLOCATED, target);
3096c5df90dSjeremylt   CeedOperatorSetField(op_error, "error", ceeddata[numlevels-1]->Erestrictui,
3106c5df90dSjeremylt                        CEED_NOTRANSPOSE, CEED_BASIS_COLLOCATED,
3116c5df90dSjeremylt                        CEED_VECTOR_ACTIVE);
3126c5df90dSjeremylt 
3136c5df90dSjeremylt   // Calculate multiplicity
3146c5df90dSjeremylt   for (int i=0; i<numlevels; i++) {
3156c5df90dSjeremylt     PetscScalar *x;
3166c5df90dSjeremylt 
3176c5df90dSjeremylt     // CEED vector
3186c5df90dSjeremylt     ierr = VecGetArray(Xloc[i], &x); CHKERRQ(ierr);
3196c5df90dSjeremylt     CeedVectorSetArray(ceeddata[i]->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
3206c5df90dSjeremylt 
3216c5df90dSjeremylt     // Multiplicity
3226c5df90dSjeremylt     CeedElemRestrictionGetMultiplicity(ceeddata[i]->Erestrictu,
3236c5df90dSjeremylt                                        ceeddata[i]->xceed);
3246c5df90dSjeremylt 
3256c5df90dSjeremylt     // Restore vector
3266c5df90dSjeremylt     ierr = VecRestoreArray(Xloc[i], &x); CHKERRQ(ierr);
3276c5df90dSjeremylt 
3286c5df90dSjeremylt     // Creat mult vector
3296c5df90dSjeremylt     ierr = VecDuplicate(Xloc[i], &mult[i]); CHKERRQ(ierr);
3306c5df90dSjeremylt 
3316c5df90dSjeremylt     // Local-to-global
3326c5df90dSjeremylt     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
3336c5df90dSjeremylt     ierr = DMLocalToGlobalBegin(dm[i], Xloc[i], ADD_VALUES, X[i]);
3346c5df90dSjeremylt     CHKERRQ(ierr);
3356c5df90dSjeremylt     ierr = DMLocalToGlobalEnd(dm[i], Xloc[i], ADD_VALUES, X[i]);
3366c5df90dSjeremylt     CHKERRQ(ierr);
3376c5df90dSjeremylt     ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr);
3386c5df90dSjeremylt 
3396c5df90dSjeremylt     // Global-to-local
3406c5df90dSjeremylt     ierr = DMGlobalToLocalBegin(dm[i], X[i], INSERT_VALUES, mult[i]);
3416c5df90dSjeremylt     CHKERRQ(ierr);
3426c5df90dSjeremylt     ierr = DMGlobalToLocalEnd(dm[i], X[i], INSERT_VALUES, mult[i]);
3436c5df90dSjeremylt     CHKERRQ(ierr);
3446c5df90dSjeremylt     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
3456c5df90dSjeremylt 
3466c5df90dSjeremylt     // Multiplicity scaling
3476c5df90dSjeremylt     ierr = VecReciprocal(mult[i]);
3486c5df90dSjeremylt   }
3496c5df90dSjeremylt 
3506c5df90dSjeremylt   // Set up Mat
3516c5df90dSjeremylt   for (int i=0; i<numlevels; i++) {
352226c3a8fSjeremylt     // User Operator
3536c5df90dSjeremylt     userO[i]->comm = comm;
3546c5df90dSjeremylt     userO[i]->dm = dm[i];
3556c5df90dSjeremylt     userO[i]->Xloc = Xloc[i];
3566c5df90dSjeremylt     ierr = VecDuplicate(Xloc[i], &userO[i]->Yloc); CHKERRQ(ierr);
3576c5df90dSjeremylt     userO[i]->xceed = ceeddata[i]->xceed;
3586c5df90dSjeremylt     userO[i]->yceed = ceeddata[i]->yceed;
3596c5df90dSjeremylt     userO[i]->op = ceeddata[i]->op_apply;
3606c5df90dSjeremylt     userO[i]->ceed = ceed;
3616c5df90dSjeremylt 
3626c5df90dSjeremylt     // Set up diagonal
3636c5df90dSjeremylt     const CeedScalar *ceedarray;
3646c5df90dSjeremylt     PetscScalar *petscarray;
3656c5df90dSjeremylt     CeedInt length;
3666c5df90dSjeremylt 
3676c5df90dSjeremylt     ierr = VecDuplicate(X[i], &userO[i]->diag); CHKERRQ(ierr);
3686c5df90dSjeremylt     ierr = VecDuplicate(Xloc[i], &diagloc); CHKERRQ(ierr);
3696c5df90dSjeremylt 
3706c5df90dSjeremylt     // -- Local diagonal
3716c5df90dSjeremylt     CeedOperatorAssembleLinearDiagonal(userO[i]->op, &diagceed,
3726c5df90dSjeremylt                                        CEED_REQUEST_IMMEDIATE);
3736c5df90dSjeremylt 
3746c5df90dSjeremylt     // -- Copy values
3756c5df90dSjeremylt     CeedVectorGetArrayRead(diagceed, CEED_MEM_HOST, &ceedarray);
3766c5df90dSjeremylt     ierr = VecGetArray(diagloc, &petscarray); CHKERRQ(ierr);
3776c5df90dSjeremylt     CeedVectorGetLength(diagceed, &length);
3786c5df90dSjeremylt     for (CeedInt i=0; i<length; i++)
3796c5df90dSjeremylt       petscarray[i] = ceedarray[i];
3806c5df90dSjeremylt     CeedVectorRestoreArrayRead(diagceed, &ceedarray);
3816c5df90dSjeremylt     ierr = VecRestoreArray(diagloc, &petscarray); CHKERRQ(ierr);
3826c5df90dSjeremylt 
3836c5df90dSjeremylt     // -- Global diagonal
3846c5df90dSjeremylt     ierr = VecZeroEntries(userO[i]->diag); CHKERRQ(ierr);
3856c5df90dSjeremylt     ierr = DMLocalToGlobalBegin(userO[i]->dm, diagloc, ADD_VALUES,
3866c5df90dSjeremylt                                 userO[i]->diag); CHKERRQ(ierr);
3876c5df90dSjeremylt     ierr = DMLocalToGlobalEnd(userO[i]->dm, diagloc, ADD_VALUES,
3886c5df90dSjeremylt                               userO[i]->diag); CHKERRQ(ierr);
3896c5df90dSjeremylt 
3906c5df90dSjeremylt     // -- Cleanup
3916c5df90dSjeremylt     ierr = VecDestroy(&diagloc); CHKERRQ(ierr);
3926c5df90dSjeremylt     CeedVectorDestroy(&diagceed);
3936c5df90dSjeremylt 
3946c5df90dSjeremylt     if (i > 0) {
3956c5df90dSjeremylt       // Interp Operator
3966c5df90dSjeremylt       userI[i]->comm = comm;
3976c5df90dSjeremylt       userI[i]->dmc = dm[i-1];
3986c5df90dSjeremylt       userI[i]->dmf = dm[i];
3996c5df90dSjeremylt       userI[i]->Xloc = Xloc[i-1];
4006c5df90dSjeremylt       userI[i]->Yloc = userO[i]->Yloc;
4016c5df90dSjeremylt       userI[i]->mult = mult[i];
4026c5df90dSjeremylt       userI[i]->ceedvecc = userO[i-1]->xceed;
4036c5df90dSjeremylt       userI[i]->ceedvecf = userO[i]->yceed;
4046c5df90dSjeremylt       userI[i]->op = ceeddata[i]->op_interp;
4056c5df90dSjeremylt       userI[i]->ceed = ceed;
4066c5df90dSjeremylt 
4076c5df90dSjeremylt       // Restrict Operator
4086c5df90dSjeremylt       userR[i]->comm = comm;
4096c5df90dSjeremylt       userR[i]->dmc = dm[i-1];
4106c5df90dSjeremylt       userR[i]->dmf = dm[i];
4116c5df90dSjeremylt       userR[i]->Xloc = Xloc[i];
4126c5df90dSjeremylt       userR[i]->Yloc = userO[i-1]->Yloc;
4136c5df90dSjeremylt       userR[i]->mult = mult[i];
4146c5df90dSjeremylt       userR[i]->ceedvecf = userO[i]->xceed;
4156c5df90dSjeremylt       userR[i]->ceedvecc = userO[i-1]->yceed;
4166c5df90dSjeremylt       userR[i]->op = ceeddata[i]->op_restrict;
4176c5df90dSjeremylt       userR[i]->ceed = ceed;
4186c5df90dSjeremylt     }
4196c5df90dSjeremylt   }
4206c5df90dSjeremylt 
4216c5df90dSjeremylt   // Set up KSP
4226c5df90dSjeremylt   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
4236c5df90dSjeremylt   {
4246c5df90dSjeremylt     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
4256c5df90dSjeremylt     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
4266c5df90dSjeremylt     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
4276c5df90dSjeremylt                             PETSC_DEFAULT); CHKERRQ(ierr);
4286c5df90dSjeremylt   }
4296c5df90dSjeremylt   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
4306c5df90dSjeremylt   ierr = KSPSetOperators(ksp, matO[numlevels-1], matO[numlevels-1]);
4316c5df90dSjeremylt   CHKERRQ(ierr);
4326c5df90dSjeremylt 
4336c5df90dSjeremylt   // Set up PCMG
4346c5df90dSjeremylt   ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
4356c5df90dSjeremylt   PCMGCycleType pcgmcycletype = PC_MG_CYCLE_V;
4366c5df90dSjeremylt   {
4376c5df90dSjeremylt     ierr = PCSetType(pc, PCMG); CHKERRQ(ierr);
4386c5df90dSjeremylt 
4396c5df90dSjeremylt     // PCMG levels
4406c5df90dSjeremylt     ierr = PCMGSetLevels(pc, numlevels, NULL); CHKERRQ(ierr);
4416c5df90dSjeremylt     for (int i=0; i<numlevels; i++) {
4426c5df90dSjeremylt       // Smoother
4436c5df90dSjeremylt       KSP smoother;
4446c5df90dSjeremylt       PC smoother_pc;
4456c5df90dSjeremylt       ierr = PCMGGetSmoother(pc, i, &smoother); CHKERRQ(ierr);
4466c5df90dSjeremylt       ierr = KSPSetType(smoother, KSPCHEBYSHEV); CHKERRQ(ierr);
4476c5df90dSjeremylt       ierr = KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1); CHKERRQ(ierr);
4486c5df90dSjeremylt       ierr = KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE); CHKERRQ(ierr);
4496c5df90dSjeremylt       ierr = KSPSetOperators(smoother, matO[i], matO[i]); CHKERRQ(ierr);
4506c5df90dSjeremylt       ierr = KSPGetPC(smoother, &smoother_pc); CHKERRQ(ierr);
4516c5df90dSjeremylt       ierr = PCSetType(smoother_pc, PCJACOBI); CHKERRQ(ierr);
4526c5df90dSjeremylt       ierr = PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
4536c5df90dSjeremylt 
4546c5df90dSjeremylt       // Work vector
4556c5df90dSjeremylt       if (i < numlevels-1) {
4566c5df90dSjeremylt         ierr = PCMGSetX(pc, i, X[i]); CHKERRQ(ierr);
4576c5df90dSjeremylt       }
4586c5df90dSjeremylt 
4596c5df90dSjeremylt       // Level transfers
4606c5df90dSjeremylt       if (i > 0) {
4616c5df90dSjeremylt         // Interpolation
4626c5df90dSjeremylt         ierr = PCMGSetInterpolation(pc, i, matI[i]); CHKERRQ(ierr);
4636c5df90dSjeremylt 
4646c5df90dSjeremylt         // Restriction
4656c5df90dSjeremylt         ierr = PCMGSetRestriction(pc, i, matR[i]); CHKERRQ(ierr);
4666c5df90dSjeremylt       }
4676c5df90dSjeremylt 
4686c5df90dSjeremylt       // Coarse solve
4696c5df90dSjeremylt       KSP coarse;
4706c5df90dSjeremylt       PC coarse_pc;
4716c5df90dSjeremylt       ierr = PCMGGetCoarseSolve(pc, &coarse); CHKERRQ(ierr);
4726c5df90dSjeremylt       ierr = KSPSetType(coarse, KSPCG); CHKERRQ(ierr);
4736c5df90dSjeremylt       ierr = KSPSetOperators(coarse, matO[0], matO[0]); CHKERRQ(ierr);
4746c5df90dSjeremylt       ierr = KSPSetTolerances(coarse, 1e-10, 1e-10, PETSC_DEFAULT,
4756c5df90dSjeremylt                               PETSC_DEFAULT); CHKERRQ(ierr);
4766c5df90dSjeremylt       ierr = KSPGetPC(coarse, &coarse_pc); CHKERRQ(ierr);
4776c5df90dSjeremylt       ierr = PCSetType(coarse_pc, PCJACOBI); CHKERRQ(ierr);
4786c5df90dSjeremylt       ierr = PCJacobiSetType(coarse_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
4796c5df90dSjeremylt     }
4806c5df90dSjeremylt 
4816c5df90dSjeremylt     // PCMG options
4826c5df90dSjeremylt     ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr);
4836c5df90dSjeremylt     ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr);
4846c5df90dSjeremylt     ierr = PCMGSetCycleType(pc, pcgmcycletype); CHKERRQ(ierr);
4856c5df90dSjeremylt   }
4866c5df90dSjeremylt 
4876c5df90dSjeremylt   // First run, if benchmarking
4886c5df90dSjeremylt   if (benchmark_mode) {
4896c5df90dSjeremylt     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
4906c5df90dSjeremylt     CHKERRQ(ierr);
4916c5df90dSjeremylt     ierr = VecZeroEntries(X[numlevels-1]); CHKERRQ(ierr);
4926c5df90dSjeremylt     my_rt_start = MPI_Wtime();
4936c5df90dSjeremylt     ierr = KSPSolve(ksp, rhs, X[numlevels-1]); CHKERRQ(ierr);
4946c5df90dSjeremylt     my_rt = MPI_Wtime() - my_rt_start;
4956c5df90dSjeremylt     ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
4966c5df90dSjeremylt     CHKERRQ(ierr);
4976c5df90dSjeremylt     // Set maxits based on first iteration timing
4986c5df90dSjeremylt     if (my_rt > 0.02) {
4996c5df90dSjeremylt       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
5006c5df90dSjeremylt       CHKERRQ(ierr);
5016c5df90dSjeremylt     } else {
5026c5df90dSjeremylt       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
5036c5df90dSjeremylt       CHKERRQ(ierr);
5046c5df90dSjeremylt     }
5056c5df90dSjeremylt   }
5066c5df90dSjeremylt 
5076c5df90dSjeremylt   // Timed solve
5086c5df90dSjeremylt   ierr = VecZeroEntries(X[numlevels-1]); CHKERRQ(ierr);
5096c5df90dSjeremylt   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
5106c5df90dSjeremylt   my_rt_start = MPI_Wtime();
5116c5df90dSjeremylt   ierr = KSPSolve(ksp, rhs, X[numlevels-1]); CHKERRQ(ierr);
5126c5df90dSjeremylt   my_rt = MPI_Wtime() - my_rt_start;
5136c5df90dSjeremylt 
5146c5df90dSjeremylt   // Output results
5156c5df90dSjeremylt   {
5166c5df90dSjeremylt     KSPType ksptype;
5176c5df90dSjeremylt     PCMGType pcmgtype;
5186c5df90dSjeremylt     KSPConvergedReason reason;
5196c5df90dSjeremylt     PetscReal rnorm;
5206c5df90dSjeremylt     PetscInt its;
5216c5df90dSjeremylt     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
5226c5df90dSjeremylt     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
5236c5df90dSjeremylt     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
5246c5df90dSjeremylt     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
5256c5df90dSjeremylt     ierr = PCMGGetType(pc, &pcmgtype); CHKERRQ(ierr);
5266c5df90dSjeremylt     if (!test_mode || reason < 0 || rnorm > 1e-8) {
5276c5df90dSjeremylt       ierr = PetscPrintf(comm,
5286c5df90dSjeremylt                          "  KSP:\n"
5296c5df90dSjeremylt                          "    KSP Type                           : %s\n"
5306c5df90dSjeremylt                          "    KSP Convergence                    : %s\n"
5316c5df90dSjeremylt                          "    Total KSP Iterations               : %D\n"
5326c5df90dSjeremylt                          "    Final rnorm                        : %e\n",
5336c5df90dSjeremylt                          ksptype, KSPConvergedReasons[reason], its,
5346c5df90dSjeremylt                          (double)rnorm); CHKERRQ(ierr);
5356c5df90dSjeremylt       ierr = PetscPrintf(comm,
5366c5df90dSjeremylt                          "  PCMG:\n"
5376c5df90dSjeremylt                          "    PCMG Type                          : %s\n"
5386c5df90dSjeremylt                          "    PCMG Cycle Type                    : %s\n",
5396c5df90dSjeremylt                          PCMGTypes[pcmgtype],
5406c5df90dSjeremylt                          PCMGCycleTypes[pcgmcycletype]); CHKERRQ(ierr);
5416c5df90dSjeremylt     }
5426c5df90dSjeremylt     if (!test_mode) {
5436c5df90dSjeremylt       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
5446c5df90dSjeremylt     }
5456c5df90dSjeremylt     {
5466c5df90dSjeremylt       PetscReal maxerror;
5476c5df90dSjeremylt       ierr = ComputeErrorMax(userO[numlevels-1], op_error, X[numlevels-1], target,
5486c5df90dSjeremylt                              &maxerror); CHKERRQ(ierr);
5496c5df90dSjeremylt       PetscReal tol = 5e-2;
5506c5df90dSjeremylt       if (!test_mode || maxerror > tol) {
5516c5df90dSjeremylt         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
5526c5df90dSjeremylt         CHKERRQ(ierr);
5536c5df90dSjeremylt         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
5546c5df90dSjeremylt         CHKERRQ(ierr);
5556c5df90dSjeremylt         ierr = PetscPrintf(comm,
5566c5df90dSjeremylt                            "    Pointwise Error (max)              : %e\n"
5576c5df90dSjeremylt                            "    CG Solve Time                      : %g (%g) sec\n",
5586c5df90dSjeremylt                            (double)maxerror, rt_max, rt_min); CHKERRQ(ierr);
5596c5df90dSjeremylt       }
5606c5df90dSjeremylt     }
5616c5df90dSjeremylt     if (benchmark_mode && (!test_mode)) {
5626c5df90dSjeremylt       ierr = PetscPrintf(comm,
5636c5df90dSjeremylt                          "    DoFs/Sec in CG                     : %g (%g) million\n",
5646c5df90dSjeremylt                          1e-6*gsize[numlevels-1]*its/rt_max,
5656c5df90dSjeremylt                          1e-6*gsize[numlevels-1]*its/rt_min);
5666c5df90dSjeremylt       CHKERRQ(ierr);
5676c5df90dSjeremylt     }
5686c5df90dSjeremylt   }
5696c5df90dSjeremylt 
5706c5df90dSjeremylt   if (write_solution) {
5716c5df90dSjeremylt     PetscViewer vtkviewersoln;
5726c5df90dSjeremylt 
5736c5df90dSjeremylt     ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr);
5746c5df90dSjeremylt     ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr);
5756c5df90dSjeremylt     ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr);
5766c5df90dSjeremylt     ierr = VecView(X[numlevels-1], vtkviewersoln); CHKERRQ(ierr);
5776c5df90dSjeremylt     ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr);
5786c5df90dSjeremylt   }
5796c5df90dSjeremylt 
5806c5df90dSjeremylt   // Cleanup
5816c5df90dSjeremylt   for (int i=0; i<numlevels; i++) {
5826c5df90dSjeremylt     ierr = VecDestroy(&X[i]); CHKERRQ(ierr);
5836c5df90dSjeremylt     ierr = VecDestroy(&Xloc[i]); CHKERRQ(ierr);
5846c5df90dSjeremylt     ierr = VecDestroy(&mult[i]); CHKERRQ(ierr);
5856c5df90dSjeremylt     ierr = VecDestroy(&userO[i]->Yloc); CHKERRQ(ierr);
5866c5df90dSjeremylt     ierr = VecDestroy(&userO[i]->diag); CHKERRQ(ierr);
5876c5df90dSjeremylt     ierr = MatDestroy(&matO[i]); CHKERRQ(ierr);
5886c5df90dSjeremylt     ierr = PetscFree(userO[i]); CHKERRQ(ierr);
5896c5df90dSjeremylt     if (i > 0) {
5906c5df90dSjeremylt       ierr = MatDestroy(&matI[i]); CHKERRQ(ierr);
5916c5df90dSjeremylt       ierr = PetscFree(userI[i]); CHKERRQ(ierr);
5926c5df90dSjeremylt       ierr = MatDestroy(&matR[i]); CHKERRQ(ierr);
5936c5df90dSjeremylt       ierr = PetscFree(userR[i]); CHKERRQ(ierr);
5946c5df90dSjeremylt     }
5956c5df90dSjeremylt     ierr = CeedDataDestroy(i, ceeddata[i]); CHKERRQ(ierr);
5966c5df90dSjeremylt     ierr = DMDestroy(&dm[i]); CHKERRQ(ierr);
5976c5df90dSjeremylt   }
5986c5df90dSjeremylt   ierr = PetscFree(leveldegrees); CHKERRQ(ierr);
5996c5df90dSjeremylt   ierr = PetscFree(dm); CHKERRQ(ierr);
6006c5df90dSjeremylt   ierr = PetscFree(X); CHKERRQ(ierr);
6016c5df90dSjeremylt   ierr = PetscFree(Xloc); CHKERRQ(ierr);
6026c5df90dSjeremylt   ierr = PetscFree(mult); CHKERRQ(ierr);
6036c5df90dSjeremylt   ierr = PetscFree(matO); CHKERRQ(ierr);
6046c5df90dSjeremylt   ierr = PetscFree(matI); CHKERRQ(ierr);
6056c5df90dSjeremylt   ierr = PetscFree(matR); CHKERRQ(ierr);
6066c5df90dSjeremylt   ierr = PetscFree(ceeddata); CHKERRQ(ierr);
6076c5df90dSjeremylt   ierr = PetscFree(userO); CHKERRQ(ierr);
6086c5df90dSjeremylt   ierr = PetscFree(userI); CHKERRQ(ierr);
6096c5df90dSjeremylt   ierr = PetscFree(userR); CHKERRQ(ierr);
6106c5df90dSjeremylt   ierr = PetscFree(lsize); CHKERRQ(ierr);
6116c5df90dSjeremylt   ierr = PetscFree(xlsize); CHKERRQ(ierr);
6126c5df90dSjeremylt   ierr = PetscFree(gsize); CHKERRQ(ierr);
6136c5df90dSjeremylt   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
6146c5df90dSjeremylt   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
6156c5df90dSjeremylt   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
6166c5df90dSjeremylt   ierr = DMDestroy(&dmOrig); CHKERRQ(ierr);
6176c5df90dSjeremylt   CeedVectorDestroy(&target);
6186c5df90dSjeremylt   CeedQFunctionDestroy(&qf_error);
6196c5df90dSjeremylt   CeedQFunctionDestroy(&qf_restrict);
6206c5df90dSjeremylt   CeedQFunctionDestroy(&qf_prolong);
6216c5df90dSjeremylt   CeedOperatorDestroy(&op_error);
6226c5df90dSjeremylt   CeedDestroy(&ceed);
6236c5df90dSjeremylt   return PetscFinalize();
6246c5df90dSjeremylt }
625