xref: /libCEED/examples/petsc/multigrid.c (revision cfa59c5b0c41ee3c7ca5285b7ce2579e871e5f29)
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 
416c5df90dSjeremylt #define multigrid
426c5df90dSjeremylt #include "setup.h"
436c5df90dSjeremylt 
44*cfa59c5bSRey // Transition from a value of "a" for x=0, to a value of "b" for x=1.  Optionally
45*cfa59c5bSRey // smooth -- see the commented versions at the end.
46*cfa59c5bSRey static double step(const double a, const double b, double x) {
47*cfa59c5bSRey   if (x <= 0) return a;
48*cfa59c5bSRey   if (x >= 1) return b;
49*cfa59c5bSRey   return a + (b-a) * (x);
50*cfa59c5bSRey }
51*cfa59c5bSRey 
52*cfa59c5bSRey // 1D transformation at the right boundary
53*cfa59c5bSRey static double right(const double eps, const double x) {
54*cfa59c5bSRey   return (x <= 0.5) ? (2-eps) * x : 1 + eps*(x-1);
55*cfa59c5bSRey }
56*cfa59c5bSRey 
57*cfa59c5bSRey // 1D transformation at the left boundary
58*cfa59c5bSRey static double left(const double eps, const double x) {
59*cfa59c5bSRey   return 1-right(eps,1-x);
60*cfa59c5bSRey }
61*cfa59c5bSRey 
62*cfa59c5bSRey // Apply 3D Kershaw mesh transformation
63*cfa59c5bSRey // The eps parameters are in (0, 1]
64*cfa59c5bSRey // Uniform mesh is recovered for eps=1
65*cfa59c5bSRey static PetscErrorCode kershaw(DM dmorig, PetscScalar eps) {
66*cfa59c5bSRey   PetscErrorCode ierr;
67*cfa59c5bSRey   Vec coord;
68*cfa59c5bSRey   PetscInt ncoord;
69*cfa59c5bSRey   PetscScalar *c;
70*cfa59c5bSRey 
71*cfa59c5bSRey   PetscFunctionBeginUser;
72*cfa59c5bSRey   ierr = DMGetCoordinatesLocal(dmorig, &coord); CHKERRQ(ierr);
73*cfa59c5bSRey   ierr = VecGetLocalSize(coord, &ncoord); CHKERRQ(ierr);
74*cfa59c5bSRey   ierr = VecGetArray(coord, &c); CHKERRQ(ierr);
75*cfa59c5bSRey 
76*cfa59c5bSRey   for (PetscInt i = 0; i < ncoord; i += 3) {
77*cfa59c5bSRey     PetscScalar x = c[i], y = c[i+1], z = c[i+2];
78*cfa59c5bSRey     PetscInt layer = x*6;
79*cfa59c5bSRey     PetscScalar lambda = (x-layer/6.0)*6;
80*cfa59c5bSRey     c[i] = x;
81*cfa59c5bSRey 
82*cfa59c5bSRey     switch (layer) {
83*cfa59c5bSRey     case 0:
84*cfa59c5bSRey       c[i+1] = left(eps, y);
85*cfa59c5bSRey       c[i+2] = left(eps, z);
86*cfa59c5bSRey       break;
87*cfa59c5bSRey     case 1:
88*cfa59c5bSRey     case 4:
89*cfa59c5bSRey       c[i+1] = step(left(eps, y), right(eps, y), lambda);
90*cfa59c5bSRey       c[i+2] = step(left(eps, z), right(eps, z), lambda);
91*cfa59c5bSRey       break;
92*cfa59c5bSRey     case 2:
93*cfa59c5bSRey       c[i+1] = step(right(eps, y), left(eps, y), lambda/2);
94*cfa59c5bSRey       c[i+2] = step(right(eps, z), left(eps, z), lambda/2);
95*cfa59c5bSRey       break;
96*cfa59c5bSRey     case 3:
97*cfa59c5bSRey       c[i+1] = step(right(eps, y), left(eps, y), (1+lambda)/2);
98*cfa59c5bSRey       c[i+2] = step(right(eps, z), left(eps, z), (1+lambda)/2);
99*cfa59c5bSRey       break;
100*cfa59c5bSRey     default:
101*cfa59c5bSRey       c[i+1] = right(eps, y);
102*cfa59c5bSRey       c[i+2] = right(eps, z);
103*cfa59c5bSRey     }
104*cfa59c5bSRey   }
105*cfa59c5bSRey   ierr = VecRestoreArray(coord, &c); CHKERRQ(ierr);
106*cfa59c5bSRey   PetscFunctionReturn(0);
107*cfa59c5bSRey }
108*cfa59c5bSRey 
1096c5df90dSjeremylt int main(int argc, char **argv) {
1106c5df90dSjeremylt   PetscInt ierr;
1116c5df90dSjeremylt   MPI_Comm comm;
112cb0b5415Sjeremylt   char filename[PETSC_MAX_PATH_LEN],
113cb0b5415Sjeremylt        ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self";
1146c5df90dSjeremylt   double my_rt_start, my_rt, rt_min, rt_max;
11561608365Sjeremylt   PetscInt degree = 3, qextra, *lsize, *xlsize, *gsize, dim = 3, fineLevel,
1166c5df90dSjeremylt            melem[3] = {3, 3, 3}, ncompu = 1, numlevels = degree, *leveldegrees;
1176c5df90dSjeremylt   PetscScalar *r;
118*cfa59c5bSRey   PetscScalar eps = 1.0;
1196c5df90dSjeremylt   PetscBool test_mode, benchmark_mode, read_mesh, write_solution;
12009a940d7Sjeremylt   PetscLogStage solvestage;
121199551f5Sjeremylt   DM  *dm, dmorig;
122c70bd2a0Sjeremylt   SNES snesdummy;
1236c5df90dSjeremylt   KSP ksp;
1246c5df90dSjeremylt   PC pc;
125199551f5Sjeremylt   Mat *matO, *matPR, matcoarse;
126cdb3667fSjeremylt   Vec *X, *Xloc, *mult, rhs, rhsloc;
127b68a8d79SJed Brown   PetscMemType memtype;
1286c5df90dSjeremylt   UserO *userO;
129a97643b0Sjeremylt   UserProlongRestr *userPR;
1306c5df90dSjeremylt   Ceed ceed;
1316c5df90dSjeremylt   CeedData *ceeddata;
1325e81177dSjeremylt   CeedVector rhsceed, target;
133c70bd2a0Sjeremylt   CeedQFunction qferror, qfrestrict, qfprolong;
134c70bd2a0Sjeremylt   CeedOperator operror;
135199551f5Sjeremylt   bpType bpchoice;
1366c5df90dSjeremylt   coarsenType coarsen;
1376c5df90dSjeremylt 
1386c5df90dSjeremylt   ierr = PetscInitialize(&argc, &argv, NULL, help);
1396c5df90dSjeremylt   if (ierr) return ierr;
1406c5df90dSjeremylt   comm = PETSC_COMM_WORLD;
1416c5df90dSjeremylt 
1426c5df90dSjeremylt   // Parse command line options
1436c5df90dSjeremylt   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr);
144199551f5Sjeremylt   bpchoice = CEED_BP3;
1456c5df90dSjeremylt   ierr = PetscOptionsEnum("-problem",
1466c5df90dSjeremylt                           "CEED benchmark problem to solve", NULL,
147199551f5Sjeremylt                           bpTypes, (PetscEnum)bpchoice, (PetscEnum *)&bpchoice,
1486c5df90dSjeremylt                           NULL); CHKERRQ(ierr);
149199551f5Sjeremylt   ncompu = bpOptions[bpchoice].ncompu;
1506c5df90dSjeremylt   test_mode = PETSC_FALSE;
1516c5df90dSjeremylt   ierr = PetscOptionsBool("-test",
1526c5df90dSjeremylt                           "Testing mode (do not print unless error is large)",
1536c5df90dSjeremylt                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
1546c5df90dSjeremylt   benchmark_mode = PETSC_FALSE;
1556c5df90dSjeremylt   ierr = PetscOptionsBool("-benchmark",
1566c5df90dSjeremylt                           "Benchmarking mode (prints benchmark statistics)",
1576c5df90dSjeremylt                           NULL, benchmark_mode, &benchmark_mode, NULL);
1586c5df90dSjeremylt   CHKERRQ(ierr);
1596c5df90dSjeremylt   write_solution = PETSC_FALSE;
1606c5df90dSjeremylt   ierr = PetscOptionsBool("-write_solution",
1616c5df90dSjeremylt                           "Write solution for visualization",
1626c5df90dSjeremylt                           NULL, write_solution, &write_solution, NULL);
1636c5df90dSjeremylt   CHKERRQ(ierr);
164*cfa59c5bSRey   ierr = PetscOptionsScalar("-eps",
165*cfa59c5bSRey                             "Epsilon parameter for Kershaw mesh transformation",
166*cfa59c5bSRey                             NULL, eps, &eps, NULL);
167*cfa59c5bSRey   if (eps > 1 || eps <= 0) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
168*cfa59c5bSRey                                       "-eps %D must be (0,1]", eps);
1696c5df90dSjeremylt   degree = test_mode ? 3 : 2;
1706c5df90dSjeremylt   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
1716c5df90dSjeremylt                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
1726c5df90dSjeremylt   if (degree < 1) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
1736c5df90dSjeremylt                              "-degree %D must be at least 1", degree);
174199551f5Sjeremylt   qextra = bpOptions[bpchoice].qextra;
1756c5df90dSjeremylt   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
1766c5df90dSjeremylt                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
1776c5df90dSjeremylt   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
1786c5df90dSjeremylt                             NULL, ceedresource, ceedresource,
1796c5df90dSjeremylt                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
1806c5df90dSjeremylt   coarsen = COARSEN_UNIFORM;
1816c5df90dSjeremylt   ierr = PetscOptionsEnum("-coarsen",
1826c5df90dSjeremylt                           "Coarsening strategy to use", NULL,
1836c5df90dSjeremylt                           coarsenTypes, (PetscEnum)coarsen,
1846c5df90dSjeremylt                           (PetscEnum *)&coarsen, NULL); CHKERRQ(ierr);
185cb32e2e7SValeria Barra   read_mesh = PETSC_FALSE;
1866c5df90dSjeremylt   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL,
1876c5df90dSjeremylt                             filename, filename, sizeof(filename), &read_mesh);
1886c5df90dSjeremylt   CHKERRQ(ierr);
1896c5df90dSjeremylt   if (!read_mesh) {
1906c5df90dSjeremylt     PetscInt tmp = dim;
1916c5df90dSjeremylt     ierr = PetscOptionsIntArray("-cells","Number of cells per dimension", NULL,
1926c5df90dSjeremylt                                 melem, &tmp, NULL); CHKERRQ(ierr);
1936c5df90dSjeremylt   }
1946c5df90dSjeremylt   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
1956c5df90dSjeremylt 
1969396343dSjeremylt   // Set up libCEED
1979396343dSjeremylt   CeedInit(ceedresource, &ceed);
1989396343dSjeremylt   CeedMemType memtypebackend;
1999396343dSjeremylt   CeedGetPreferredMemType(ceed, &memtypebackend);
2009396343dSjeremylt 
2016c5df90dSjeremylt   // Setup DM
2026c5df90dSjeremylt   if (read_mesh) {
203199551f5Sjeremylt     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dmorig);
2046c5df90dSjeremylt     CHKERRQ(ierr);
2056c5df90dSjeremylt   } else {
2066c5df90dSjeremylt     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_FALSE, melem, NULL,
207199551f5Sjeremylt                                NULL, NULL, PETSC_TRUE,&dmorig); CHKERRQ(ierr);
2086c5df90dSjeremylt   }
2096c5df90dSjeremylt 
2106c5df90dSjeremylt   {
2116c5df90dSjeremylt     DM dmDist = NULL;
2126c5df90dSjeremylt     PetscPartitioner part;
2136c5df90dSjeremylt 
214199551f5Sjeremylt     ierr = DMPlexGetPartitioner(dmorig, &part); CHKERRQ(ierr);
2156c5df90dSjeremylt     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
216199551f5Sjeremylt     ierr = DMPlexDistribute(dmorig, 0, NULL, &dmDist); CHKERRQ(ierr);
2176c5df90dSjeremylt     if (dmDist) {
218199551f5Sjeremylt       ierr = DMDestroy(&dmorig); CHKERRQ(ierr);
219199551f5Sjeremylt       dmorig = dmDist;
2206c5df90dSjeremylt     }
2216c5df90dSjeremylt   }
2226c5df90dSjeremylt 
223*cfa59c5bSRey   // apply Kershaw mesh transformation
224*cfa59c5bSRey   ierr = kershaw(dmorig, eps); CHKERRQ(ierr);
225*cfa59c5bSRey 
226b68a8d79SJed Brown   VecType vectype;
227b68a8d79SJed Brown   switch (memtypebackend) {
228b68a8d79SJed Brown   case CEED_MEM_HOST: vectype = VECSTANDARD; break;
229b68a8d79SJed Brown   case CEED_MEM_DEVICE: {
230b68a8d79SJed Brown     const char *resolved;
231b68a8d79SJed Brown     CeedGetResource(ceed, &resolved);
232b68a8d79SJed Brown     if (strstr(resolved, "/gpu/cuda")) vectype = VECCUDA;
233ca2d516cSJed Brown     else if (strstr(resolved, "/gpu/hip/occa"))
234ca2d516cSJed Brown       vectype = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678
235b68a8d79SJed Brown     else if (strstr(resolved, "/gpu/hip")) vectype = VECHIP;
236b68a8d79SJed Brown     else vectype = VECSTANDARD;
237b68a8d79SJed Brown   }
238b68a8d79SJed Brown   }
239b68a8d79SJed Brown   ierr = DMSetVecType(dmorig, vectype); CHKERRQ(ierr);
240b68a8d79SJed Brown   ierr = DMSetFromOptions(dmorig); CHKERRQ(ierr);
241b68a8d79SJed Brown 
2426c5df90dSjeremylt   // Allocate arrays for PETSc objects for each level
2436c5df90dSjeremylt   switch (coarsen) {
2446c5df90dSjeremylt   case COARSEN_UNIFORM:
2456c5df90dSjeremylt     numlevels = degree;
2466c5df90dSjeremylt     break;
247dc7d240cSValeria Barra   case COARSEN_LOGARITHMIC:
2486c5df90dSjeremylt     numlevels = ceil(log(degree)/log(2)) + 1;
2496c5df90dSjeremylt     break;
2506c5df90dSjeremylt   }
2516c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &leveldegrees); CHKERRQ(ierr);
25261608365Sjeremylt   fineLevel = numlevels - 1;
25361608365Sjeremylt 
2546c5df90dSjeremylt   switch (coarsen) {
2556c5df90dSjeremylt   case COARSEN_UNIFORM:
2566c5df90dSjeremylt     for (int i=0; i<numlevels; i++) leveldegrees[i] = i + 1;
2576c5df90dSjeremylt     break;
258dc7d240cSValeria Barra   case COARSEN_LOGARITHMIC:
2596c5df90dSjeremylt     for (int i=0; i<numlevels - 1; i++) leveldegrees[i] = pow(2,i);
26061608365Sjeremylt     leveldegrees[fineLevel] = degree;
2616c5df90dSjeremylt     break;
2626c5df90dSjeremylt   }
2636c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &dm); CHKERRQ(ierr);
2646c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &X); CHKERRQ(ierr);
2656c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &Xloc); CHKERRQ(ierr);
2666c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &mult); CHKERRQ(ierr);
2676c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &userO); CHKERRQ(ierr);
268a97643b0Sjeremylt   ierr = PetscMalloc1(numlevels, &userPR); CHKERRQ(ierr);
2696c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &matO); CHKERRQ(ierr);
270a97643b0Sjeremylt   ierr = PetscMalloc1(numlevels, &matPR); CHKERRQ(ierr);
2716c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &lsize); CHKERRQ(ierr);
2726c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &xlsize); CHKERRQ(ierr);
2736c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &gsize); CHKERRQ(ierr);
2746c5df90dSjeremylt 
2756c5df90dSjeremylt   // Setup DM and Operator Mat Shells for each level
2766c5df90dSjeremylt   for (CeedInt i=0; i<numlevels; i++) {
2776c5df90dSjeremylt     // Create DM
278199551f5Sjeremylt     ierr = DMClone(dmorig, &dm[i]); CHKERRQ(ierr);
279b68a8d79SJed Brown     ierr = DMGetVecType(dmorig, &vectype); CHKERRQ(ierr);
280b68a8d79SJed Brown     ierr = DMSetVecType(dm[i], vectype); CHKERRQ(ierr);
281199551f5Sjeremylt     ierr = SetupDMByDegree(dm[i], leveldegrees[i], ncompu, bpchoice);
2826c5df90dSjeremylt     CHKERRQ(ierr);
2836c5df90dSjeremylt 
2846c5df90dSjeremylt     // Create vectors
2856c5df90dSjeremylt     ierr = DMCreateGlobalVector(dm[i], &X[i]); CHKERRQ(ierr);
2866c5df90dSjeremylt     ierr = VecGetLocalSize(X[i], &lsize[i]); CHKERRQ(ierr);
2876c5df90dSjeremylt     ierr = VecGetSize(X[i], &gsize[i]); CHKERRQ(ierr);
2886c5df90dSjeremylt     ierr = DMCreateLocalVector(dm[i], &Xloc[i]); CHKERRQ(ierr);
2896c5df90dSjeremylt     ierr = VecGetSize(Xloc[i], &xlsize[i]); CHKERRQ(ierr);
2906c5df90dSjeremylt 
2916c5df90dSjeremylt     // Operator
2926c5df90dSjeremylt     ierr = PetscMalloc1(1, &userO[i]); CHKERRQ(ierr);
2936c5df90dSjeremylt     ierr = MatCreateShell(comm, lsize[i], lsize[i], gsize[i], gsize[i],
2946c5df90dSjeremylt                           userO[i], &matO[i]); CHKERRQ(ierr);
2956c5df90dSjeremylt     ierr = MatShellSetOperation(matO[i], MATOP_MULT,
296ce74dcefSjeremylt                                 (void(*)(void))MatMult_Ceed); CHKERRQ(ierr);
2976c5df90dSjeremylt     ierr = MatShellSetOperation(matO[i], MATOP_GET_DIAGONAL,
298ce74dcefSjeremylt                                 (void(*)(void))MatGetDiag); CHKERRQ(ierr);
299b68a8d79SJed Brown     ierr = MatShellSetVecType(matO[i], vectype); CHKERRQ(ierr);
3006c5df90dSjeremylt 
3016c5df90dSjeremylt     // Level transfers
3026c5df90dSjeremylt     if (i > 0) {
3036c5df90dSjeremylt       // Interp
304a97643b0Sjeremylt       ierr = PetscMalloc1(1, &userPR[i]); CHKERRQ(ierr);
3056c5df90dSjeremylt       ierr = MatCreateShell(comm, lsize[i], lsize[i-1], gsize[i], gsize[i-1],
306a97643b0Sjeremylt                             userPR[i], &matPR[i]); CHKERRQ(ierr);
307a97643b0Sjeremylt       ierr = MatShellSetOperation(matPR[i], MATOP_MULT,
308a97643b0Sjeremylt                                   (void(*)(void))MatMult_Prolong);
3096c5df90dSjeremylt       CHKERRQ(ierr);
310a97643b0Sjeremylt       ierr = MatShellSetOperation(matPR[i], MATOP_MULT_TRANSPOSE,
3116c5df90dSjeremylt                                   (void(*)(void))MatMult_Restrict);
3126c5df90dSjeremylt       CHKERRQ(ierr);
313b68a8d79SJed Brown       ierr = MatShellSetVecType(matPR[i], vectype); CHKERRQ(ierr);
3146c5df90dSjeremylt     }
3156c5df90dSjeremylt   }
31661608365Sjeremylt   ierr = VecDuplicate(X[fineLevel], &rhs); CHKERRQ(ierr);
3176c5df90dSjeremylt 
3186c5df90dSjeremylt   // Print global grid information
3196c5df90dSjeremylt   if (!test_mode) {
3206c5df90dSjeremylt     PetscInt P = degree + 1, Q = P + qextra;
3219396343dSjeremylt 
3222d03409cSjeremylt     const char *usedresource;
3232d03409cSjeremylt     CeedGetResource(ceed, &usedresource);
3249396343dSjeremylt 
3259396343dSjeremylt     ierr = VecGetType(X[0], &vectype); CHKERRQ(ierr);
3269396343dSjeremylt 
3276c5df90dSjeremylt     ierr = PetscPrintf(comm,
3286c5df90dSjeremylt                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc + PCMG --\n"
3299396343dSjeremylt                        "  PETSc:\n"
3309396343dSjeremylt                        "    PETSc Vec Type                     : %s\n"
3316c5df90dSjeremylt                        "  libCEED:\n"
3326c5df90dSjeremylt                        "    libCEED Backend                    : %s\n"
3339396343dSjeremylt                        "    libCEED Backend MemType            : %s\n"
3346c5df90dSjeremylt                        "  Mesh:\n"
3356c5df90dSjeremylt                        "    Number of 1D Basis Nodes (p)       : %d\n"
3366c5df90dSjeremylt                        "    Number of 1D Quadrature Points (q) : %d\n"
3376c5df90dSjeremylt                        "    Global Nodes                       : %D\n"
3386c5df90dSjeremylt                        "    Owned Nodes                        : %D\n"
339db419314Sjeremylt                        "    DoF per node                       : %D\n"
3406c5df90dSjeremylt                        "  Multigrid:\n"
3416c5df90dSjeremylt                        "    Number of Levels                   : %d\n",
3429396343dSjeremylt                        bpchoice+1, vectype, usedresource,
3439396343dSjeremylt                        CeedMemTypes[memtypebackend],
3449396343dSjeremylt                        P, Q, gsize[fineLevel]/ncompu, lsize[fineLevel]/ncompu,
345db419314Sjeremylt                        ncompu, numlevels); CHKERRQ(ierr);
3466c5df90dSjeremylt   }
3476c5df90dSjeremylt 
3486c5df90dSjeremylt   // Create RHS vector
34961608365Sjeremylt   ierr = VecDuplicate(Xloc[fineLevel], &rhsloc); CHKERRQ(ierr);
3506c5df90dSjeremylt   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
351b68a8d79SJed Brown   ierr = VecGetArrayAndMemType(rhsloc, &r, &memtype); CHKERRQ(ierr);
35261608365Sjeremylt   CeedVectorCreate(ceed, xlsize[fineLevel], &rhsceed);
353b68a8d79SJed Brown   CeedVectorSetArray(rhsceed, MemTypeP2C(memtype), CEED_USE_POINTER, r);
3546c5df90dSjeremylt 
3556c5df90dSjeremylt   // Set up libCEED operators on each level
3566c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &ceeddata); CHKERRQ(ierr);
3576c5df90dSjeremylt   for (int i=0; i<numlevels; i++) {
3586c5df90dSjeremylt     // Print level information
35961608365Sjeremylt     if (!test_mode && (i == 0 || i == fineLevel)) {
3606c5df90dSjeremylt       ierr = PetscPrintf(comm,"    Level %D (%s):\n"
3616c5df90dSjeremylt                          "      Number of 1D Basis Nodes (p)     : %d\n"
3626c5df90dSjeremylt                          "      Global Nodes                     : %D\n"
3636c5df90dSjeremylt                          "      Owned Nodes                      : %D\n",
3646c5df90dSjeremylt                          i, (i? "fine" : "coarse"), leveldegrees[i] + 1,
3656c5df90dSjeremylt                          gsize[i]/ncompu, lsize[i]/ncompu); CHKERRQ(ierr);
3666c5df90dSjeremylt     }
3676c5df90dSjeremylt     ierr = PetscMalloc1(1, &ceeddata[i]); CHKERRQ(ierr);
3686c5df90dSjeremylt     ierr = SetupLibceedByDegree(dm[i], ceed, leveldegrees[i], dim, qextra,
369199551f5Sjeremylt                                 ncompu, gsize[i], xlsize[i], bpchoice,
370f9342789SJeremy L Thompson                                 ceeddata[i], i==(fineLevel), rhsceed, &target);
371f9342789SJeremy L Thompson     CHKERRQ(ierr);
3726c5df90dSjeremylt   }
3736c5df90dSjeremylt 
3746c5df90dSjeremylt   // Gather RHS
375b68a8d79SJed Brown   CeedVectorTakeArray(rhsceed, MemTypeP2C(memtype), NULL);
376b68a8d79SJed Brown   ierr = VecRestoreArrayAndMemType(rhsloc, &r); CHKERRQ(ierr);
3776c5df90dSjeremylt   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
3789396343dSjeremylt   ierr = DMLocalToGlobal(dm[fineLevel], rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr);
3796c5df90dSjeremylt   CeedVectorDestroy(&rhsceed);
3806c5df90dSjeremylt 
38143eb8658SJeremy L Thompson   // Create the restriction/interpolation QFunction
38260f77c51Sjeremylt   CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_NONE, CEED_EVAL_INTERP,
383c70bd2a0Sjeremylt                               &qfrestrict);
38460f77c51Sjeremylt   CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_INTERP, CEED_EVAL_NONE,
385c70bd2a0Sjeremylt                               &qfprolong);
3866c5df90dSjeremylt 
3876c5df90dSjeremylt   // Set up libCEED level transfer operators
388199551f5Sjeremylt   ierr = CeedLevelTransferSetup(ceed, numlevels, ncompu, bpchoice, ceeddata,
389c70bd2a0Sjeremylt                                 leveldegrees, qfrestrict, qfprolong);
3906c5df90dSjeremylt   CHKERRQ(ierr);
3916c5df90dSjeremylt 
39243eb8658SJeremy L Thompson   // Create the error QFunction
393199551f5Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].error,
394199551f5Sjeremylt                               bpOptions[bpchoice].errorfname, &qferror);
395c70bd2a0Sjeremylt   CeedQFunctionAddInput(qferror, "u", ncompu, CEED_EVAL_INTERP);
396c70bd2a0Sjeremylt   CeedQFunctionAddInput(qferror, "true_soln", ncompu, CEED_EVAL_NONE);
397c70bd2a0Sjeremylt   CeedQFunctionAddOutput(qferror, "error", ncompu, CEED_EVAL_NONE);
3986c5df90dSjeremylt 
3996c5df90dSjeremylt   // Create the error operator
400c70bd2a0Sjeremylt   CeedOperatorCreate(ceed, qferror, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
401c70bd2a0Sjeremylt                      &operror);
402c70bd2a0Sjeremylt   CeedOperatorSetField(operror, "u", ceeddata[fineLevel]->Erestrictu,
40361608365Sjeremylt                        ceeddata[fineLevel]->basisu, CEED_VECTOR_ACTIVE);
404c70bd2a0Sjeremylt   CeedOperatorSetField(operror, "true_soln", ceeddata[fineLevel]->Erestrictui,
405a8d32208Sjeremylt                        CEED_BASIS_COLLOCATED, target);
406c70bd2a0Sjeremylt   CeedOperatorSetField(operror, "error", ceeddata[fineLevel]->Erestrictui,
407a8d32208Sjeremylt                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
4086c5df90dSjeremylt 
4096c5df90dSjeremylt   // Calculate multiplicity
4106c5df90dSjeremylt   for (int i=0; i<numlevels; i++) {
4116c5df90dSjeremylt     PetscScalar *x;
4126c5df90dSjeremylt 
4136c5df90dSjeremylt     // CEED vector
414f9342789SJeremy L Thompson     ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr);
4156c5df90dSjeremylt     ierr = VecGetArray(Xloc[i], &x); CHKERRQ(ierr);
4166c5df90dSjeremylt     CeedVectorSetArray(ceeddata[i]->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
4176c5df90dSjeremylt 
4186c5df90dSjeremylt     // Multiplicity
419a8d32208Sjeremylt     CeedElemRestrictionGetMultiplicity(ceeddata[i]->Erestrictu,
4206c5df90dSjeremylt                                        ceeddata[i]->xceed);
421f9342789SJeremy L Thompson     CeedVectorSyncArray(ceeddata[i]->xceed, CEED_MEM_HOST);
4226c5df90dSjeremylt 
4236c5df90dSjeremylt     // Restore vector
4246c5df90dSjeremylt     ierr = VecRestoreArray(Xloc[i], &x); CHKERRQ(ierr);
4256c5df90dSjeremylt 
4266c5df90dSjeremylt     // Creat mult vector
4276c5df90dSjeremylt     ierr = VecDuplicate(Xloc[i], &mult[i]); CHKERRQ(ierr);
4286c5df90dSjeremylt 
4296c5df90dSjeremylt     // Local-to-global
4306c5df90dSjeremylt     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
431483f8b0dSjeremylt     ierr = DMLocalToGlobal(dm[i], Xloc[i], ADD_VALUES, X[i]);
4326c5df90dSjeremylt     CHKERRQ(ierr);
4336c5df90dSjeremylt     ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr);
4346c5df90dSjeremylt 
4356c5df90dSjeremylt     // Global-to-local
436483f8b0dSjeremylt     ierr = DMGlobalToLocal(dm[i], X[i], INSERT_VALUES, mult[i]);
4376c5df90dSjeremylt     CHKERRQ(ierr);
4386c5df90dSjeremylt     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
4396c5df90dSjeremylt 
4406c5df90dSjeremylt     // Multiplicity scaling
4416c5df90dSjeremylt     ierr = VecReciprocal(mult[i]);
4426c5df90dSjeremylt   }
4436c5df90dSjeremylt 
4446c5df90dSjeremylt   // Set up Mat
4456c5df90dSjeremylt   for (int i=0; i<numlevels; i++) {
446226c3a8fSjeremylt     // User Operator
4476c5df90dSjeremylt     userO[i]->comm = comm;
4486c5df90dSjeremylt     userO[i]->dm = dm[i];
4496c5df90dSjeremylt     userO[i]->Xloc = Xloc[i];
4506c5df90dSjeremylt     ierr = VecDuplicate(Xloc[i], &userO[i]->Yloc); CHKERRQ(ierr);
4516c5df90dSjeremylt     userO[i]->xceed = ceeddata[i]->xceed;
4526c5df90dSjeremylt     userO[i]->yceed = ceeddata[i]->yceed;
453c70bd2a0Sjeremylt     userO[i]->op = ceeddata[i]->opapply;
4546c5df90dSjeremylt     userO[i]->ceed = ceed;
4556c5df90dSjeremylt 
4566c5df90dSjeremylt     if (i > 0) {
457a97643b0Sjeremylt       // Prolongation/Restriction Operator
458a97643b0Sjeremylt       userPR[i]->comm = comm;
459199551f5Sjeremylt       userPR[i]->dmf = dm[i];
460199551f5Sjeremylt       userPR[i]->dmc = dm[i-1];
461199551f5Sjeremylt       userPR[i]->locvecc = Xloc[i-1];
462199551f5Sjeremylt       userPR[i]->locvecf = userO[i]->Yloc;
463199551f5Sjeremylt       userPR[i]->multvec = mult[i];
464199551f5Sjeremylt       userPR[i]->ceedvecc = userO[i-1]->xceed;
465199551f5Sjeremylt       userPR[i]->ceedvecf = userO[i]->yceed;
466c70bd2a0Sjeremylt       userPR[i]->opprolong = ceeddata[i]->opprolong;
467c70bd2a0Sjeremylt       userPR[i]->oprestrict = ceeddata[i]->oprestrict;
468a97643b0Sjeremylt       userPR[i]->ceed = ceed;
4696c5df90dSjeremylt     }
4706c5df90dSjeremylt   }
4716c5df90dSjeremylt 
47215ce0ef0Sjeremylt   // Setup dummy SNES for AMG coarse solve
473c70bd2a0Sjeremylt   ierr = SNESCreate(comm, &snesdummy); CHKERRQ(ierr);
474c70bd2a0Sjeremylt   ierr = SNESSetDM(snesdummy, dm[0]); CHKERRQ(ierr);
475c70bd2a0Sjeremylt   ierr = SNESSetSolution(snesdummy, X[0]); CHKERRQ(ierr);
47615ce0ef0Sjeremylt 
47715ce0ef0Sjeremylt   // -- Jacobian matrix
47815ce0ef0Sjeremylt   ierr = DMSetMatType(dm[0], MATAIJ); CHKERRQ(ierr);
479199551f5Sjeremylt   ierr = DMCreateMatrix(dm[0], &matcoarse); CHKERRQ(ierr);
480199551f5Sjeremylt   ierr = SNESSetJacobian(snesdummy, matcoarse, matcoarse, NULL,
48115ce0ef0Sjeremylt                          NULL); CHKERRQ(ierr);
48215ce0ef0Sjeremylt 
48315ce0ef0Sjeremylt   // -- Residual evaluation function
484c70bd2a0Sjeremylt   ierr = SNESSetFunction(snesdummy, X[0], FormResidual_Ceed,
48515ce0ef0Sjeremylt                          userO[0]); CHKERRQ(ierr);
48615ce0ef0Sjeremylt 
48715ce0ef0Sjeremylt   // -- Form Jacobian
488c70bd2a0Sjeremylt   ierr = SNESComputeJacobianDefaultColor(snesdummy, X[0], matO[0],
489199551f5Sjeremylt                                          matcoarse, NULL); CHKERRQ(ierr);
49015ce0ef0Sjeremylt 
4916c5df90dSjeremylt   // Set up KSP
4926c5df90dSjeremylt   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
4936c5df90dSjeremylt   {
4946c5df90dSjeremylt     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
4956c5df90dSjeremylt     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
4966c5df90dSjeremylt     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
4976c5df90dSjeremylt                             PETSC_DEFAULT); CHKERRQ(ierr);
4986c5df90dSjeremylt   }
4996c5df90dSjeremylt   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
50061608365Sjeremylt   ierr = KSPSetOperators(ksp, matO[fineLevel], matO[fineLevel]);
5016c5df90dSjeremylt   CHKERRQ(ierr);
5026c5df90dSjeremylt 
5036c5df90dSjeremylt   // Set up PCMG
5046c5df90dSjeremylt   ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
5056c5df90dSjeremylt   PCMGCycleType pcgmcycletype = PC_MG_CYCLE_V;
5066c5df90dSjeremylt   {
5076c5df90dSjeremylt     ierr = PCSetType(pc, PCMG); CHKERRQ(ierr);
5086c5df90dSjeremylt 
5096c5df90dSjeremylt     // PCMG levels
5106c5df90dSjeremylt     ierr = PCMGSetLevels(pc, numlevels, NULL); CHKERRQ(ierr);
5116c5df90dSjeremylt     for (int i=0; i<numlevels; i++) {
5126c5df90dSjeremylt       // Smoother
5136c5df90dSjeremylt       KSP smoother;
5146c5df90dSjeremylt       PC smoother_pc;
5156c5df90dSjeremylt       ierr = PCMGGetSmoother(pc, i, &smoother); CHKERRQ(ierr);
5166c5df90dSjeremylt       ierr = KSPSetType(smoother, KSPCHEBYSHEV); CHKERRQ(ierr);
5176c5df90dSjeremylt       ierr = KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1); CHKERRQ(ierr);
5186c5df90dSjeremylt       ierr = KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE); CHKERRQ(ierr);
5196c5df90dSjeremylt       ierr = KSPSetOperators(smoother, matO[i], matO[i]); CHKERRQ(ierr);
5206c5df90dSjeremylt       ierr = KSPGetPC(smoother, &smoother_pc); CHKERRQ(ierr);
5216c5df90dSjeremylt       ierr = PCSetType(smoother_pc, PCJACOBI); CHKERRQ(ierr);
5226c5df90dSjeremylt       ierr = PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
5236c5df90dSjeremylt 
5246c5df90dSjeremylt       // Work vector
5256c5df90dSjeremylt       if (i < numlevels - 1) {
5266c5df90dSjeremylt         ierr = PCMGSetX(pc, i, X[i]); CHKERRQ(ierr);
5276c5df90dSjeremylt       }
5286c5df90dSjeremylt 
5296c5df90dSjeremylt       // Level transfers
5306c5df90dSjeremylt       if (i > 0) {
5316c5df90dSjeremylt         // Interpolation
532a97643b0Sjeremylt         ierr = PCMGSetInterpolation(pc, i, matPR[i]); CHKERRQ(ierr);
5336c5df90dSjeremylt       }
5346c5df90dSjeremylt 
5356c5df90dSjeremylt       // Coarse solve
5366c5df90dSjeremylt       KSP coarse;
5376c5df90dSjeremylt       PC coarse_pc;
5386c5df90dSjeremylt       ierr = PCMGGetCoarseSolve(pc, &coarse); CHKERRQ(ierr);
53915ce0ef0Sjeremylt       ierr = KSPSetType(coarse, KSPPREONLY); CHKERRQ(ierr);
540199551f5Sjeremylt       ierr = KSPSetOperators(coarse, matcoarse, matcoarse); CHKERRQ(ierr);
54115ce0ef0Sjeremylt 
5426c5df90dSjeremylt       ierr = KSPGetPC(coarse, &coarse_pc); CHKERRQ(ierr);
54315ce0ef0Sjeremylt       ierr = PCSetType(coarse_pc, PCGAMG); CHKERRQ(ierr);
54415ce0ef0Sjeremylt 
54515ce0ef0Sjeremylt       ierr = KSPSetOptionsPrefix(coarse, "coarse_"); CHKERRQ(ierr);
54615ce0ef0Sjeremylt       ierr = PCSetOptionsPrefix(coarse_pc, "coarse_"); CHKERRQ(ierr);
54715ce0ef0Sjeremylt       ierr = KSPSetFromOptions(coarse); CHKERRQ(ierr);
54815ce0ef0Sjeremylt       ierr = PCSetFromOptions(coarse_pc); CHKERRQ(ierr);
5496c5df90dSjeremylt     }
5506c5df90dSjeremylt 
5516c5df90dSjeremylt     // PCMG options
5526c5df90dSjeremylt     ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr);
5536c5df90dSjeremylt     ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr);
5546c5df90dSjeremylt     ierr = PCMGSetCycleType(pc, pcgmcycletype); CHKERRQ(ierr);
5556c5df90dSjeremylt   }
5566c5df90dSjeremylt 
5576c5df90dSjeremylt   // First run, if benchmarking
5586c5df90dSjeremylt   if (benchmark_mode) {
5596c5df90dSjeremylt     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
5606c5df90dSjeremylt     CHKERRQ(ierr);
56161608365Sjeremylt     ierr = VecZeroEntries(X[fineLevel]); CHKERRQ(ierr);
5626c5df90dSjeremylt     my_rt_start = MPI_Wtime();
56361608365Sjeremylt     ierr = KSPSolve(ksp, rhs, X[fineLevel]); CHKERRQ(ierr);
5646c5df90dSjeremylt     my_rt = MPI_Wtime() - my_rt_start;
5656c5df90dSjeremylt     ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
5666c5df90dSjeremylt     CHKERRQ(ierr);
5676c5df90dSjeremylt     // Set maxits based on first iteration timing
5686c5df90dSjeremylt     if (my_rt > 0.02) {
5696c5df90dSjeremylt       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
5706c5df90dSjeremylt       CHKERRQ(ierr);
5716c5df90dSjeremylt     } else {
5726c5df90dSjeremylt       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
5736c5df90dSjeremylt       CHKERRQ(ierr);
5746c5df90dSjeremylt     }
5756c5df90dSjeremylt   }
5766c5df90dSjeremylt 
5776c5df90dSjeremylt   // Timed solve
57861608365Sjeremylt   ierr = VecZeroEntries(X[fineLevel]); CHKERRQ(ierr);
5796c5df90dSjeremylt   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
58009a940d7Sjeremylt 
58109a940d7Sjeremylt   // -- Performance logging
58209a940d7Sjeremylt   ierr = PetscLogStageRegister("Solve Stage", &solvestage); CHKERRQ(ierr);
58309a940d7Sjeremylt   ierr = PetscLogStagePush(solvestage); CHKERRQ(ierr);
58409a940d7Sjeremylt 
58509a940d7Sjeremylt   // -- Solve
5866c5df90dSjeremylt   my_rt_start = MPI_Wtime();
58761608365Sjeremylt   ierr = KSPSolve(ksp, rhs, X[fineLevel]); CHKERRQ(ierr);
5886c5df90dSjeremylt   my_rt = MPI_Wtime() - my_rt_start;
5896c5df90dSjeremylt 
59009a940d7Sjeremylt 
59109a940d7Sjeremylt   // -- Performance logging
59209a940d7Sjeremylt   ierr = PetscLogStagePop();
59309a940d7Sjeremylt 
5946c5df90dSjeremylt   // Output results
5956c5df90dSjeremylt   {
5966c5df90dSjeremylt     KSPType ksptype;
5976c5df90dSjeremylt     PCMGType pcmgtype;
5986c5df90dSjeremylt     KSPConvergedReason reason;
5996c5df90dSjeremylt     PetscReal rnorm;
6006c5df90dSjeremylt     PetscInt its;
6016c5df90dSjeremylt     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
6026c5df90dSjeremylt     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
6036c5df90dSjeremylt     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
6046c5df90dSjeremylt     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
6056c5df90dSjeremylt     ierr = PCMGGetType(pc, &pcmgtype); CHKERRQ(ierr);
6066c5df90dSjeremylt     if (!test_mode || reason < 0 || rnorm > 1e-8) {
6076c5df90dSjeremylt       ierr = PetscPrintf(comm,
6086c5df90dSjeremylt                          "  KSP:\n"
6096c5df90dSjeremylt                          "    KSP Type                           : %s\n"
6106c5df90dSjeremylt                          "    KSP Convergence                    : %s\n"
6116c5df90dSjeremylt                          "    Total KSP Iterations               : %D\n"
6126c5df90dSjeremylt                          "    Final rnorm                        : %e\n",
6136c5df90dSjeremylt                          ksptype, KSPConvergedReasons[reason], its,
6146c5df90dSjeremylt                          (double)rnorm); CHKERRQ(ierr);
6156c5df90dSjeremylt       ierr = PetscPrintf(comm,
6166c5df90dSjeremylt                          "  PCMG:\n"
6176c5df90dSjeremylt                          "    PCMG Type                          : %s\n"
6186c5df90dSjeremylt                          "    PCMG Cycle Type                    : %s\n",
6196c5df90dSjeremylt                          PCMGTypes[pcmgtype],
6206c5df90dSjeremylt                          PCMGCycleTypes[pcgmcycletype]); CHKERRQ(ierr);
6216c5df90dSjeremylt     }
6226c5df90dSjeremylt     if (!test_mode) {
6236c5df90dSjeremylt       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
6246c5df90dSjeremylt     }
6256c5df90dSjeremylt     {
6266c5df90dSjeremylt       PetscReal maxerror;
627c70bd2a0Sjeremylt       ierr = ComputeErrorMax(userO[fineLevel], operror, X[fineLevel], target,
6286c5df90dSjeremylt                              &maxerror); CHKERRQ(ierr);
6296c5df90dSjeremylt       PetscReal tol = 5e-2;
6306c5df90dSjeremylt       if (!test_mode || maxerror > tol) {
6316c5df90dSjeremylt         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
6326c5df90dSjeremylt         CHKERRQ(ierr);
6336c5df90dSjeremylt         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
6346c5df90dSjeremylt         CHKERRQ(ierr);
6356c5df90dSjeremylt         ierr = PetscPrintf(comm,
6366c5df90dSjeremylt                            "    Pointwise Error (max)              : %e\n"
6376c5df90dSjeremylt                            "    CG Solve Time                      : %g (%g) sec\n",
6386c5df90dSjeremylt                            (double)maxerror, rt_max, rt_min); CHKERRQ(ierr);
6396c5df90dSjeremylt       }
6406c5df90dSjeremylt     }
6416c5df90dSjeremylt     if (benchmark_mode && (!test_mode)) {
6426c5df90dSjeremylt       ierr = PetscPrintf(comm,
6436c5df90dSjeremylt                          "    DoFs/Sec in CG                     : %g (%g) million\n",
64461608365Sjeremylt                          1e-6*gsize[fineLevel]*its/rt_max,
64561608365Sjeremylt                          1e-6*gsize[fineLevel]*its/rt_min);
6466c5df90dSjeremylt       CHKERRQ(ierr);
6476c5df90dSjeremylt     }
6486c5df90dSjeremylt   }
6496c5df90dSjeremylt 
6506c5df90dSjeremylt   if (write_solution) {
6516c5df90dSjeremylt     PetscViewer vtkviewersoln;
6526c5df90dSjeremylt 
6536c5df90dSjeremylt     ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr);
6546c5df90dSjeremylt     ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr);
655*cfa59c5bSRey     ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtu"); CHKERRQ(ierr);
65661608365Sjeremylt     ierr = VecView(X[fineLevel], vtkviewersoln); CHKERRQ(ierr);
6576c5df90dSjeremylt     ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr);
6586c5df90dSjeremylt   }
6596c5df90dSjeremylt 
6606c5df90dSjeremylt   // Cleanup
6616c5df90dSjeremylt   for (int i=0; i<numlevels; i++) {
6626c5df90dSjeremylt     ierr = VecDestroy(&X[i]); CHKERRQ(ierr);
6636c5df90dSjeremylt     ierr = VecDestroy(&Xloc[i]); CHKERRQ(ierr);
6646c5df90dSjeremylt     ierr = VecDestroy(&mult[i]); CHKERRQ(ierr);
6656c5df90dSjeremylt     ierr = VecDestroy(&userO[i]->Yloc); CHKERRQ(ierr);
6666c5df90dSjeremylt     ierr = MatDestroy(&matO[i]); CHKERRQ(ierr);
6676c5df90dSjeremylt     ierr = PetscFree(userO[i]); CHKERRQ(ierr);
6686c5df90dSjeremylt     if (i > 0) {
669a97643b0Sjeremylt       ierr = MatDestroy(&matPR[i]); CHKERRQ(ierr);
670a97643b0Sjeremylt       ierr = PetscFree(userPR[i]); CHKERRQ(ierr);
6716c5df90dSjeremylt     }
6726c5df90dSjeremylt     ierr = CeedDataDestroy(i, ceeddata[i]); CHKERRQ(ierr);
6736c5df90dSjeremylt     ierr = DMDestroy(&dm[i]); CHKERRQ(ierr);
6746c5df90dSjeremylt   }
6756c5df90dSjeremylt   ierr = PetscFree(leveldegrees); CHKERRQ(ierr);
6766c5df90dSjeremylt   ierr = PetscFree(dm); CHKERRQ(ierr);
6776c5df90dSjeremylt   ierr = PetscFree(X); CHKERRQ(ierr);
6786c5df90dSjeremylt   ierr = PetscFree(Xloc); CHKERRQ(ierr);
6796c5df90dSjeremylt   ierr = PetscFree(mult); CHKERRQ(ierr);
6806c5df90dSjeremylt   ierr = PetscFree(matO); CHKERRQ(ierr);
681a97643b0Sjeremylt   ierr = PetscFree(matPR); CHKERRQ(ierr);
6826c5df90dSjeremylt   ierr = PetscFree(ceeddata); CHKERRQ(ierr);
6836c5df90dSjeremylt   ierr = PetscFree(userO); CHKERRQ(ierr);
684a97643b0Sjeremylt   ierr = PetscFree(userPR); CHKERRQ(ierr);
6856c5df90dSjeremylt   ierr = PetscFree(lsize); CHKERRQ(ierr);
6866c5df90dSjeremylt   ierr = PetscFree(xlsize); CHKERRQ(ierr);
6876c5df90dSjeremylt   ierr = PetscFree(gsize); CHKERRQ(ierr);
6886c5df90dSjeremylt   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
6896c5df90dSjeremylt   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
690199551f5Sjeremylt   ierr = MatDestroy(&matcoarse); CHKERRQ(ierr);
6916c5df90dSjeremylt   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
692c70bd2a0Sjeremylt   ierr = SNESDestroy(&snesdummy); CHKERRQ(ierr);
693199551f5Sjeremylt   ierr = DMDestroy(&dmorig); CHKERRQ(ierr);
6946c5df90dSjeremylt   CeedVectorDestroy(&target);
695c70bd2a0Sjeremylt   CeedQFunctionDestroy(&qferror);
696c70bd2a0Sjeremylt   CeedQFunctionDestroy(&qfrestrict);
697c70bd2a0Sjeremylt   CeedQFunctionDestroy(&qfprolong);
698c70bd2a0Sjeremylt   CeedOperatorDestroy(&operror);
6996c5df90dSjeremylt   CeedDestroy(&ceed);
7006c5df90dSjeremylt   return PetscFinalize();
7016c5df90dSjeremylt }
702