xref: /libCEED/examples/petsc/multigrid.c (revision e83e87a50f1b2e8810b376a735430565127e4d25)
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 
41*e83e87a5Sjeremylt #include "bps.h"
426c5df90dSjeremylt 
43cfa59c5bSRey // Transition from a value of "a" for x=0, to a value of "b" for x=1.  Optionally
44cfa59c5bSRey // smooth -- see the commented versions at the end.
45cfa59c5bSRey static double step(const double a, const double b, double x) {
46cfa59c5bSRey   if (x <= 0) return a;
47cfa59c5bSRey   if (x >= 1) return b;
48cfa59c5bSRey   return a + (b-a) * (x);
49cfa59c5bSRey }
50cfa59c5bSRey 
51cfa59c5bSRey // 1D transformation at the right boundary
52cfa59c5bSRey static double right(const double eps, const double x) {
53cfa59c5bSRey   return (x <= 0.5) ? (2-eps) * x : 1 + eps*(x-1);
54cfa59c5bSRey }
55cfa59c5bSRey 
56cfa59c5bSRey // 1D transformation at the left boundary
57cfa59c5bSRey static double left(const double eps, const double x) {
58cfa59c5bSRey   return 1-right(eps,1-x);
59cfa59c5bSRey }
60cfa59c5bSRey 
61cfa59c5bSRey // Apply 3D Kershaw mesh transformation
62cfa59c5bSRey // The eps parameters are in (0, 1]
63cfa59c5bSRey // Uniform mesh is recovered for eps=1
64cfa59c5bSRey static PetscErrorCode kershaw(DM dmorig, PetscScalar eps) {
65cfa59c5bSRey   PetscErrorCode ierr;
66cfa59c5bSRey   Vec coord;
67cfa59c5bSRey   PetscInt ncoord;
68cfa59c5bSRey   PetscScalar *c;
69cfa59c5bSRey 
70cfa59c5bSRey   PetscFunctionBeginUser;
71cfa59c5bSRey   ierr = DMGetCoordinatesLocal(dmorig, &coord); CHKERRQ(ierr);
72cfa59c5bSRey   ierr = VecGetLocalSize(coord, &ncoord); CHKERRQ(ierr);
73cfa59c5bSRey   ierr = VecGetArray(coord, &c); CHKERRQ(ierr);
74cfa59c5bSRey 
75cfa59c5bSRey   for (PetscInt i = 0; i < ncoord; i += 3) {
76cfa59c5bSRey     PetscScalar x = c[i], y = c[i+1], z = c[i+2];
77cfa59c5bSRey     PetscInt layer = x*6;
78cfa59c5bSRey     PetscScalar lambda = (x-layer/6.0)*6;
79cfa59c5bSRey     c[i] = x;
80cfa59c5bSRey 
81cfa59c5bSRey     switch (layer) {
82cfa59c5bSRey     case 0:
83cfa59c5bSRey       c[i+1] = left(eps, y);
84cfa59c5bSRey       c[i+2] = left(eps, z);
85cfa59c5bSRey       break;
86cfa59c5bSRey     case 1:
87cfa59c5bSRey     case 4:
88cfa59c5bSRey       c[i+1] = step(left(eps, y), right(eps, y), lambda);
89cfa59c5bSRey       c[i+2] = step(left(eps, z), right(eps, z), lambda);
90cfa59c5bSRey       break;
91cfa59c5bSRey     case 2:
92cfa59c5bSRey       c[i+1] = step(right(eps, y), left(eps, y), lambda/2);
93cfa59c5bSRey       c[i+2] = step(right(eps, z), left(eps, z), lambda/2);
94cfa59c5bSRey       break;
95cfa59c5bSRey     case 3:
96cfa59c5bSRey       c[i+1] = step(right(eps, y), left(eps, y), (1+lambda)/2);
97cfa59c5bSRey       c[i+2] = step(right(eps, z), left(eps, z), (1+lambda)/2);
98cfa59c5bSRey       break;
99cfa59c5bSRey     default:
100cfa59c5bSRey       c[i+1] = right(eps, y);
101cfa59c5bSRey       c[i+2] = right(eps, z);
102cfa59c5bSRey     }
103cfa59c5bSRey   }
104cfa59c5bSRey   ierr = VecRestoreArray(coord, &c); CHKERRQ(ierr);
105cfa59c5bSRey   PetscFunctionReturn(0);
106cfa59c5bSRey }
107cfa59c5bSRey 
1086c5df90dSjeremylt int main(int argc, char **argv) {
1096c5df90dSjeremylt   PetscInt ierr;
1106c5df90dSjeremylt   MPI_Comm comm;
111cb0b5415Sjeremylt   char filename[PETSC_MAX_PATH_LEN],
112cb0b5415Sjeremylt        ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self";
1136c5df90dSjeremylt   double my_rt_start, my_rt, rt_min, rt_max;
11461608365Sjeremylt   PetscInt degree = 3, qextra, *lsize, *xlsize, *gsize, dim = 3, fineLevel,
1156c5df90dSjeremylt            melem[3] = {3, 3, 3}, ncompu = 1, numlevels = degree, *leveldegrees;
1166c5df90dSjeremylt   PetscScalar *r;
117cfa59c5bSRey   PetscScalar eps = 1.0;
1186c5df90dSjeremylt   PetscBool test_mode, benchmark_mode, read_mesh, write_solution;
11909a940d7Sjeremylt   PetscLogStage solvestage;
120199551f5Sjeremylt   DM  *dm, dmorig;
121c70bd2a0Sjeremylt   SNES snesdummy;
1226c5df90dSjeremylt   KSP ksp;
1236c5df90dSjeremylt   PC pc;
124199551f5Sjeremylt   Mat *matO, *matPR, matcoarse;
125cdb3667fSjeremylt   Vec *X, *Xloc, *mult, rhs, rhsloc;
126b68a8d79SJed Brown   PetscMemType memtype;
1276c5df90dSjeremylt   UserO *userO;
128a97643b0Sjeremylt   UserProlongRestr *userPR;
1296c5df90dSjeremylt   Ceed ceed;
1306c5df90dSjeremylt   CeedData *ceeddata;
1315e81177dSjeremylt   CeedVector rhsceed, target;
132c70bd2a0Sjeremylt   CeedQFunction qferror, qfrestrict, qfprolong;
133c70bd2a0Sjeremylt   CeedOperator operror;
134199551f5Sjeremylt   bpType bpchoice;
1356c5df90dSjeremylt   coarsenType coarsen;
1366c5df90dSjeremylt 
1376c5df90dSjeremylt   ierr = PetscInitialize(&argc, &argv, NULL, help);
1386c5df90dSjeremylt   if (ierr) return ierr;
1396c5df90dSjeremylt   comm = PETSC_COMM_WORLD;
1406c5df90dSjeremylt 
1416c5df90dSjeremylt   // Parse command line options
1426c5df90dSjeremylt   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr);
143199551f5Sjeremylt   bpchoice = CEED_BP3;
1446c5df90dSjeremylt   ierr = PetscOptionsEnum("-problem",
1456c5df90dSjeremylt                           "CEED benchmark problem to solve", NULL,
146199551f5Sjeremylt                           bpTypes, (PetscEnum)bpchoice, (PetscEnum *)&bpchoice,
1476c5df90dSjeremylt                           NULL); CHKERRQ(ierr);
148199551f5Sjeremylt   ncompu = bpOptions[bpchoice].ncompu;
1496c5df90dSjeremylt   test_mode = PETSC_FALSE;
1506c5df90dSjeremylt   ierr = PetscOptionsBool("-test",
1516c5df90dSjeremylt                           "Testing mode (do not print unless error is large)",
1526c5df90dSjeremylt                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
1536c5df90dSjeremylt   benchmark_mode = PETSC_FALSE;
1546c5df90dSjeremylt   ierr = PetscOptionsBool("-benchmark",
1556c5df90dSjeremylt                           "Benchmarking mode (prints benchmark statistics)",
1566c5df90dSjeremylt                           NULL, benchmark_mode, &benchmark_mode, NULL);
1576c5df90dSjeremylt   CHKERRQ(ierr);
1586c5df90dSjeremylt   write_solution = PETSC_FALSE;
1596c5df90dSjeremylt   ierr = PetscOptionsBool("-write_solution",
1606c5df90dSjeremylt                           "Write solution for visualization",
1616c5df90dSjeremylt                           NULL, write_solution, &write_solution, NULL);
1626c5df90dSjeremylt   CHKERRQ(ierr);
163cfa59c5bSRey   ierr = PetscOptionsScalar("-eps",
164cfa59c5bSRey                             "Epsilon parameter for Kershaw mesh transformation",
165cfa59c5bSRey                             NULL, eps, &eps, NULL);
166cfa59c5bSRey   if (eps > 1 || eps <= 0) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
167cfa59c5bSRey                                       "-eps %D must be (0,1]", eps);
1686c5df90dSjeremylt   degree = test_mode ? 3 : 2;
1696c5df90dSjeremylt   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
1706c5df90dSjeremylt                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
1716c5df90dSjeremylt   if (degree < 1) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
1726c5df90dSjeremylt                              "-degree %D must be at least 1", degree);
173199551f5Sjeremylt   qextra = bpOptions[bpchoice].qextra;
1746c5df90dSjeremylt   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
1756c5df90dSjeremylt                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
1766c5df90dSjeremylt   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
1776c5df90dSjeremylt                             NULL, ceedresource, ceedresource,
1786c5df90dSjeremylt                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
1796c5df90dSjeremylt   coarsen = COARSEN_UNIFORM;
1806c5df90dSjeremylt   ierr = PetscOptionsEnum("-coarsen",
1816c5df90dSjeremylt                           "Coarsening strategy to use", NULL,
1826c5df90dSjeremylt                           coarsenTypes, (PetscEnum)coarsen,
1836c5df90dSjeremylt                           (PetscEnum *)&coarsen, NULL); CHKERRQ(ierr);
184cb32e2e7SValeria Barra   read_mesh = PETSC_FALSE;
1856c5df90dSjeremylt   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL,
1866c5df90dSjeremylt                             filename, filename, sizeof(filename), &read_mesh);
1876c5df90dSjeremylt   CHKERRQ(ierr);
1886c5df90dSjeremylt   if (!read_mesh) {
1896c5df90dSjeremylt     PetscInt tmp = dim;
1906c5df90dSjeremylt     ierr = PetscOptionsIntArray("-cells","Number of cells per dimension", NULL,
1916c5df90dSjeremylt                                 melem, &tmp, NULL); CHKERRQ(ierr);
1926c5df90dSjeremylt   }
1936c5df90dSjeremylt   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
1946c5df90dSjeremylt 
1959396343dSjeremylt   // Set up libCEED
1969396343dSjeremylt   CeedInit(ceedresource, &ceed);
1979396343dSjeremylt   CeedMemType memtypebackend;
1989396343dSjeremylt   CeedGetPreferredMemType(ceed, &memtypebackend);
1999396343dSjeremylt 
2006c5df90dSjeremylt   // Setup DM
2016c5df90dSjeremylt   if (read_mesh) {
202199551f5Sjeremylt     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dmorig);
2036c5df90dSjeremylt     CHKERRQ(ierr);
2046c5df90dSjeremylt   } else {
2056c5df90dSjeremylt     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_FALSE, melem, NULL,
206199551f5Sjeremylt                                NULL, NULL, PETSC_TRUE,&dmorig); CHKERRQ(ierr);
2076c5df90dSjeremylt   }
2086c5df90dSjeremylt 
2096c5df90dSjeremylt   {
2106c5df90dSjeremylt     DM dmDist = NULL;
2116c5df90dSjeremylt     PetscPartitioner part;
2126c5df90dSjeremylt 
213199551f5Sjeremylt     ierr = DMPlexGetPartitioner(dmorig, &part); CHKERRQ(ierr);
2146c5df90dSjeremylt     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
215199551f5Sjeremylt     ierr = DMPlexDistribute(dmorig, 0, NULL, &dmDist); CHKERRQ(ierr);
2166c5df90dSjeremylt     if (dmDist) {
217199551f5Sjeremylt       ierr = DMDestroy(&dmorig); CHKERRQ(ierr);
218199551f5Sjeremylt       dmorig = dmDist;
2196c5df90dSjeremylt     }
2206c5df90dSjeremylt   }
2216c5df90dSjeremylt 
222cfa59c5bSRey   // apply Kershaw mesh transformation
223cfa59c5bSRey   ierr = kershaw(dmorig, eps); CHKERRQ(ierr);
224cfa59c5bSRey 
225b68a8d79SJed Brown   VecType vectype;
226b68a8d79SJed Brown   switch (memtypebackend) {
227b68a8d79SJed Brown   case CEED_MEM_HOST: vectype = VECSTANDARD; break;
228b68a8d79SJed Brown   case CEED_MEM_DEVICE: {
229b68a8d79SJed Brown     const char *resolved;
230b68a8d79SJed Brown     CeedGetResource(ceed, &resolved);
231b68a8d79SJed Brown     if (strstr(resolved, "/gpu/cuda")) vectype = VECCUDA;
232ca2d516cSJed Brown     else if (strstr(resolved, "/gpu/hip/occa"))
233ca2d516cSJed Brown       vectype = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678
234b68a8d79SJed Brown     else if (strstr(resolved, "/gpu/hip")) vectype = VECHIP;
235b68a8d79SJed Brown     else vectype = VECSTANDARD;
236b68a8d79SJed Brown   }
237b68a8d79SJed Brown   }
238b68a8d79SJed Brown   ierr = DMSetVecType(dmorig, vectype); CHKERRQ(ierr);
239b68a8d79SJed Brown   ierr = DMSetFromOptions(dmorig); CHKERRQ(ierr);
240b68a8d79SJed Brown 
2416c5df90dSjeremylt   // Allocate arrays for PETSc objects for each level
2426c5df90dSjeremylt   switch (coarsen) {
2436c5df90dSjeremylt   case COARSEN_UNIFORM:
2446c5df90dSjeremylt     numlevels = degree;
2456c5df90dSjeremylt     break;
246dc7d240cSValeria Barra   case COARSEN_LOGARITHMIC:
2476c5df90dSjeremylt     numlevels = ceil(log(degree)/log(2)) + 1;
2486c5df90dSjeremylt     break;
2496c5df90dSjeremylt   }
2506c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &leveldegrees); CHKERRQ(ierr);
25161608365Sjeremylt   fineLevel = numlevels - 1;
25261608365Sjeremylt 
2536c5df90dSjeremylt   switch (coarsen) {
2546c5df90dSjeremylt   case COARSEN_UNIFORM:
2556c5df90dSjeremylt     for (int i=0; i<numlevels; i++) leveldegrees[i] = i + 1;
2566c5df90dSjeremylt     break;
257dc7d240cSValeria Barra   case COARSEN_LOGARITHMIC:
2586c5df90dSjeremylt     for (int i=0; i<numlevels - 1; i++) leveldegrees[i] = pow(2,i);
25961608365Sjeremylt     leveldegrees[fineLevel] = degree;
2606c5df90dSjeremylt     break;
2616c5df90dSjeremylt   }
2626c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &dm); CHKERRQ(ierr);
2636c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &X); CHKERRQ(ierr);
2646c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &Xloc); CHKERRQ(ierr);
2656c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &mult); CHKERRQ(ierr);
2666c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &userO); CHKERRQ(ierr);
267a97643b0Sjeremylt   ierr = PetscMalloc1(numlevels, &userPR); CHKERRQ(ierr);
2686c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &matO); CHKERRQ(ierr);
269a97643b0Sjeremylt   ierr = PetscMalloc1(numlevels, &matPR); CHKERRQ(ierr);
2706c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &lsize); CHKERRQ(ierr);
2716c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &xlsize); CHKERRQ(ierr);
2726c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &gsize); CHKERRQ(ierr);
2736c5df90dSjeremylt 
2746c5df90dSjeremylt   // Setup DM and Operator Mat Shells for each level
2756c5df90dSjeremylt   for (CeedInt i=0; i<numlevels; i++) {
2766c5df90dSjeremylt     // Create DM
277199551f5Sjeremylt     ierr = DMClone(dmorig, &dm[i]); CHKERRQ(ierr);
278b68a8d79SJed Brown     ierr = DMGetVecType(dmorig, &vectype); CHKERRQ(ierr);
279b68a8d79SJed Brown     ierr = DMSetVecType(dm[i], vectype); CHKERRQ(ierr);
280*e83e87a5Sjeremylt     PetscInt dim;
281*e83e87a5Sjeremylt     ierr = DMGetDimension(dm[i], &dim); CHKERRQ(ierr);
282*e83e87a5Sjeremylt     ierr = SetupDMByDegree(dm[i], leveldegrees[i], ncompu, dim,
283*e83e87a5Sjeremylt                            bpOptions[bpchoice].enforcebc, bpOptions[bpchoice].bcsfunc);
2846c5df90dSjeremylt     CHKERRQ(ierr);
2856c5df90dSjeremylt 
2866c5df90dSjeremylt     // Create vectors
2876c5df90dSjeremylt     ierr = DMCreateGlobalVector(dm[i], &X[i]); CHKERRQ(ierr);
2886c5df90dSjeremylt     ierr = VecGetLocalSize(X[i], &lsize[i]); CHKERRQ(ierr);
2896c5df90dSjeremylt     ierr = VecGetSize(X[i], &gsize[i]); CHKERRQ(ierr);
2906c5df90dSjeremylt     ierr = DMCreateLocalVector(dm[i], &Xloc[i]); CHKERRQ(ierr);
2916c5df90dSjeremylt     ierr = VecGetSize(Xloc[i], &xlsize[i]); CHKERRQ(ierr);
2926c5df90dSjeremylt 
2936c5df90dSjeremylt     // Operator
2946c5df90dSjeremylt     ierr = PetscMalloc1(1, &userO[i]); CHKERRQ(ierr);
2956c5df90dSjeremylt     ierr = MatCreateShell(comm, lsize[i], lsize[i], gsize[i], gsize[i],
2966c5df90dSjeremylt                           userO[i], &matO[i]); CHKERRQ(ierr);
2976c5df90dSjeremylt     ierr = MatShellSetOperation(matO[i], MATOP_MULT,
298ce74dcefSjeremylt                                 (void(*)(void))MatMult_Ceed); CHKERRQ(ierr);
2996c5df90dSjeremylt     ierr = MatShellSetOperation(matO[i], MATOP_GET_DIAGONAL,
300ce74dcefSjeremylt                                 (void(*)(void))MatGetDiag); CHKERRQ(ierr);
301b68a8d79SJed Brown     ierr = MatShellSetVecType(matO[i], vectype); CHKERRQ(ierr);
3026c5df90dSjeremylt 
3036c5df90dSjeremylt     // Level transfers
3046c5df90dSjeremylt     if (i > 0) {
3056c5df90dSjeremylt       // Interp
306a97643b0Sjeremylt       ierr = PetscMalloc1(1, &userPR[i]); CHKERRQ(ierr);
3076c5df90dSjeremylt       ierr = MatCreateShell(comm, lsize[i], lsize[i-1], gsize[i], gsize[i-1],
308a97643b0Sjeremylt                             userPR[i], &matPR[i]); CHKERRQ(ierr);
309a97643b0Sjeremylt       ierr = MatShellSetOperation(matPR[i], MATOP_MULT,
310a97643b0Sjeremylt                                   (void(*)(void))MatMult_Prolong);
3116c5df90dSjeremylt       CHKERRQ(ierr);
312a97643b0Sjeremylt       ierr = MatShellSetOperation(matPR[i], MATOP_MULT_TRANSPOSE,
3136c5df90dSjeremylt                                   (void(*)(void))MatMult_Restrict);
3146c5df90dSjeremylt       CHKERRQ(ierr);
315b68a8d79SJed Brown       ierr = MatShellSetVecType(matPR[i], vectype); CHKERRQ(ierr);
3166c5df90dSjeremylt     }
3176c5df90dSjeremylt   }
31861608365Sjeremylt   ierr = VecDuplicate(X[fineLevel], &rhs); CHKERRQ(ierr);
3196c5df90dSjeremylt 
3206c5df90dSjeremylt   // Print global grid information
3216c5df90dSjeremylt   if (!test_mode) {
3226c5df90dSjeremylt     PetscInt P = degree + 1, Q = P + qextra;
3239396343dSjeremylt 
3242d03409cSjeremylt     const char *usedresource;
3252d03409cSjeremylt     CeedGetResource(ceed, &usedresource);
3269396343dSjeremylt 
3279396343dSjeremylt     ierr = VecGetType(X[0], &vectype); CHKERRQ(ierr);
3289396343dSjeremylt 
3296c5df90dSjeremylt     ierr = PetscPrintf(comm,
3306c5df90dSjeremylt                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc + PCMG --\n"
3319396343dSjeremylt                        "  PETSc:\n"
3329396343dSjeremylt                        "    PETSc Vec Type                     : %s\n"
3336c5df90dSjeremylt                        "  libCEED:\n"
3346c5df90dSjeremylt                        "    libCEED Backend                    : %s\n"
3359396343dSjeremylt                        "    libCEED Backend MemType            : %s\n"
3366c5df90dSjeremylt                        "  Mesh:\n"
3376c5df90dSjeremylt                        "    Number of 1D Basis Nodes (p)       : %d\n"
3386c5df90dSjeremylt                        "    Number of 1D Quadrature Points (q) : %d\n"
3396c5df90dSjeremylt                        "    Global Nodes                       : %D\n"
3406c5df90dSjeremylt                        "    Owned Nodes                        : %D\n"
341db419314Sjeremylt                        "    DoF per node                       : %D\n"
3426c5df90dSjeremylt                        "  Multigrid:\n"
3436c5df90dSjeremylt                        "    Number of Levels                   : %d\n",
3449396343dSjeremylt                        bpchoice+1, vectype, usedresource,
3459396343dSjeremylt                        CeedMemTypes[memtypebackend],
3469396343dSjeremylt                        P, Q, gsize[fineLevel]/ncompu, lsize[fineLevel]/ncompu,
347db419314Sjeremylt                        ncompu, numlevels); CHKERRQ(ierr);
3486c5df90dSjeremylt   }
3496c5df90dSjeremylt 
3506c5df90dSjeremylt   // Create RHS vector
35161608365Sjeremylt   ierr = VecDuplicate(Xloc[fineLevel], &rhsloc); CHKERRQ(ierr);
3526c5df90dSjeremylt   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
353b68a8d79SJed Brown   ierr = VecGetArrayAndMemType(rhsloc, &r, &memtype); CHKERRQ(ierr);
35461608365Sjeremylt   CeedVectorCreate(ceed, xlsize[fineLevel], &rhsceed);
355b68a8d79SJed Brown   CeedVectorSetArray(rhsceed, MemTypeP2C(memtype), CEED_USE_POINTER, r);
3566c5df90dSjeremylt 
3576c5df90dSjeremylt   // Set up libCEED operators on each level
3586c5df90dSjeremylt   ierr = PetscMalloc1(numlevels, &ceeddata); CHKERRQ(ierr);
3596c5df90dSjeremylt   for (int i=0; i<numlevels; i++) {
3606c5df90dSjeremylt     // Print level information
36161608365Sjeremylt     if (!test_mode && (i == 0 || i == fineLevel)) {
3626c5df90dSjeremylt       ierr = PetscPrintf(comm,"    Level %D (%s):\n"
3636c5df90dSjeremylt                          "      Number of 1D Basis Nodes (p)     : %d\n"
3646c5df90dSjeremylt                          "      Global Nodes                     : %D\n"
3656c5df90dSjeremylt                          "      Owned Nodes                      : %D\n",
3666c5df90dSjeremylt                          i, (i? "fine" : "coarse"), leveldegrees[i] + 1,
3676c5df90dSjeremylt                          gsize[i]/ncompu, lsize[i]/ncompu); CHKERRQ(ierr);
3686c5df90dSjeremylt     }
3696c5df90dSjeremylt     ierr = PetscMalloc1(1, &ceeddata[i]); CHKERRQ(ierr);
3706c5df90dSjeremylt     ierr = SetupLibceedByDegree(dm[i], ceed, leveldegrees[i], dim, qextra,
371*e83e87a5Sjeremylt                                 dim, ncompu, gsize[i], xlsize[i], bpOptions[bpchoice],
372f9342789SJeremy L Thompson                                 ceeddata[i], i==(fineLevel), rhsceed, &target);
373f9342789SJeremy L Thompson     CHKERRQ(ierr);
3746c5df90dSjeremylt   }
3756c5df90dSjeremylt 
3766c5df90dSjeremylt   // Gather RHS
377b68a8d79SJed Brown   CeedVectorTakeArray(rhsceed, MemTypeP2C(memtype), NULL);
378b68a8d79SJed Brown   ierr = VecRestoreArrayAndMemType(rhsloc, &r); CHKERRQ(ierr);
3796c5df90dSjeremylt   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
3809396343dSjeremylt   ierr = DMLocalToGlobal(dm[fineLevel], rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr);
3816c5df90dSjeremylt   CeedVectorDestroy(&rhsceed);
3826c5df90dSjeremylt 
38343eb8658SJeremy L Thompson   // Create the restriction/interpolation QFunction
38460f77c51Sjeremylt   CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_NONE, CEED_EVAL_INTERP,
385c70bd2a0Sjeremylt                               &qfrestrict);
38660f77c51Sjeremylt   CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_INTERP, CEED_EVAL_NONE,
387c70bd2a0Sjeremylt                               &qfprolong);
3886c5df90dSjeremylt 
3896c5df90dSjeremylt   // Set up libCEED level transfer operators
390*e83e87a5Sjeremylt   ierr = CeedLevelTransferSetup(ceed, numlevels, ncompu, ceeddata, leveldegrees,
391*e83e87a5Sjeremylt                                 qfrestrict, qfprolong); CHKERRQ(ierr);
3926c5df90dSjeremylt 
39343eb8658SJeremy L Thompson   // Create the error QFunction
394199551f5Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].error,
395199551f5Sjeremylt                               bpOptions[bpchoice].errorfname, &qferror);
396c70bd2a0Sjeremylt   CeedQFunctionAddInput(qferror, "u", ncompu, CEED_EVAL_INTERP);
397c70bd2a0Sjeremylt   CeedQFunctionAddInput(qferror, "true_soln", ncompu, CEED_EVAL_NONE);
398c70bd2a0Sjeremylt   CeedQFunctionAddOutput(qferror, "error", ncompu, CEED_EVAL_NONE);
3996c5df90dSjeremylt 
4006c5df90dSjeremylt   // Create the error operator
401c70bd2a0Sjeremylt   CeedOperatorCreate(ceed, qferror, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
402c70bd2a0Sjeremylt                      &operror);
403c70bd2a0Sjeremylt   CeedOperatorSetField(operror, "u", ceeddata[fineLevel]->Erestrictu,
40461608365Sjeremylt                        ceeddata[fineLevel]->basisu, CEED_VECTOR_ACTIVE);
405c70bd2a0Sjeremylt   CeedOperatorSetField(operror, "true_soln", ceeddata[fineLevel]->Erestrictui,
406a8d32208Sjeremylt                        CEED_BASIS_COLLOCATED, target);
407c70bd2a0Sjeremylt   CeedOperatorSetField(operror, "error", ceeddata[fineLevel]->Erestrictui,
408a8d32208Sjeremylt                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
4096c5df90dSjeremylt 
4106c5df90dSjeremylt   // Calculate multiplicity
4116c5df90dSjeremylt   for (int i=0; i<numlevels; i++) {
4126c5df90dSjeremylt     PetscScalar *x;
4136c5df90dSjeremylt 
4146c5df90dSjeremylt     // CEED vector
415f9342789SJeremy L Thompson     ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr);
4166c5df90dSjeremylt     ierr = VecGetArray(Xloc[i], &x); CHKERRQ(ierr);
417*e83e87a5Sjeremylt     CeedVectorSetArray(ceeddata[i]->Xceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
4186c5df90dSjeremylt 
4196c5df90dSjeremylt     // Multiplicity
420a8d32208Sjeremylt     CeedElemRestrictionGetMultiplicity(ceeddata[i]->Erestrictu,
421*e83e87a5Sjeremylt                                        ceeddata[i]->Xceed);
422*e83e87a5Sjeremylt     CeedVectorSyncArray(ceeddata[i]->Xceed, CEED_MEM_HOST);
4236c5df90dSjeremylt 
4246c5df90dSjeremylt     // Restore vector
4256c5df90dSjeremylt     ierr = VecRestoreArray(Xloc[i], &x); CHKERRQ(ierr);
4266c5df90dSjeremylt 
4276c5df90dSjeremylt     // Creat mult vector
4286c5df90dSjeremylt     ierr = VecDuplicate(Xloc[i], &mult[i]); CHKERRQ(ierr);
4296c5df90dSjeremylt 
4306c5df90dSjeremylt     // Local-to-global
4316c5df90dSjeremylt     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
432483f8b0dSjeremylt     ierr = DMLocalToGlobal(dm[i], Xloc[i], ADD_VALUES, X[i]);
4336c5df90dSjeremylt     CHKERRQ(ierr);
4346c5df90dSjeremylt     ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr);
4356c5df90dSjeremylt 
4366c5df90dSjeremylt     // Global-to-local
437483f8b0dSjeremylt     ierr = DMGlobalToLocal(dm[i], X[i], INSERT_VALUES, mult[i]);
4386c5df90dSjeremylt     CHKERRQ(ierr);
4396c5df90dSjeremylt     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
4406c5df90dSjeremylt 
4416c5df90dSjeremylt     // Multiplicity scaling
4426c5df90dSjeremylt     ierr = VecReciprocal(mult[i]);
4436c5df90dSjeremylt   }
4446c5df90dSjeremylt 
4456c5df90dSjeremylt   // Set up Mat
4466c5df90dSjeremylt   for (int i=0; i<numlevels; i++) {
447226c3a8fSjeremylt     // User Operator
4486c5df90dSjeremylt     userO[i]->comm = comm;
4496c5df90dSjeremylt     userO[i]->dm = dm[i];
4506c5df90dSjeremylt     userO[i]->Xloc = Xloc[i];
4516c5df90dSjeremylt     ierr = VecDuplicate(Xloc[i], &userO[i]->Yloc); CHKERRQ(ierr);
452*e83e87a5Sjeremylt     userO[i]->Xceed = ceeddata[i]->Xceed;
453*e83e87a5Sjeremylt     userO[i]->Yceed = ceeddata[i]->Yceed;
454*e83e87a5Sjeremylt     userO[i]->op = ceeddata[i]->opApply;
4556c5df90dSjeremylt     userO[i]->ceed = ceed;
4566c5df90dSjeremylt 
4576c5df90dSjeremylt     if (i > 0) {
458a97643b0Sjeremylt       // Prolongation/Restriction Operator
459a97643b0Sjeremylt       userPR[i]->comm = comm;
460199551f5Sjeremylt       userPR[i]->dmf = dm[i];
461199551f5Sjeremylt       userPR[i]->dmc = dm[i-1];
462199551f5Sjeremylt       userPR[i]->locvecc = Xloc[i-1];
463199551f5Sjeremylt       userPR[i]->locvecf = userO[i]->Yloc;
464199551f5Sjeremylt       userPR[i]->multvec = mult[i];
465*e83e87a5Sjeremylt       userPR[i]->ceedvecc = userO[i-1]->Xceed;
466*e83e87a5Sjeremylt       userPR[i]->ceedvecf = userO[i]->Yceed;
467*e83e87a5Sjeremylt       userPR[i]->opProlong = ceeddata[i]->opProlong;
468*e83e87a5Sjeremylt       userPR[i]->opRestrict = ceeddata[i]->opRestrict;
469a97643b0Sjeremylt       userPR[i]->ceed = ceed;
4706c5df90dSjeremylt     }
4716c5df90dSjeremylt   }
4726c5df90dSjeremylt 
47315ce0ef0Sjeremylt   // Setup dummy SNES for AMG coarse solve
474c70bd2a0Sjeremylt   ierr = SNESCreate(comm, &snesdummy); CHKERRQ(ierr);
475c70bd2a0Sjeremylt   ierr = SNESSetDM(snesdummy, dm[0]); CHKERRQ(ierr);
476c70bd2a0Sjeremylt   ierr = SNESSetSolution(snesdummy, X[0]); CHKERRQ(ierr);
47715ce0ef0Sjeremylt 
47815ce0ef0Sjeremylt   // -- Jacobian matrix
47915ce0ef0Sjeremylt   ierr = DMSetMatType(dm[0], MATAIJ); CHKERRQ(ierr);
480199551f5Sjeremylt   ierr = DMCreateMatrix(dm[0], &matcoarse); CHKERRQ(ierr);
481199551f5Sjeremylt   ierr = SNESSetJacobian(snesdummy, matcoarse, matcoarse, NULL,
48215ce0ef0Sjeremylt                          NULL); CHKERRQ(ierr);
48315ce0ef0Sjeremylt 
48415ce0ef0Sjeremylt   // -- Residual evaluation function
485c70bd2a0Sjeremylt   ierr = SNESSetFunction(snesdummy, X[0], FormResidual_Ceed,
48615ce0ef0Sjeremylt                          userO[0]); CHKERRQ(ierr);
48715ce0ef0Sjeremylt 
48815ce0ef0Sjeremylt   // -- Form Jacobian
489c70bd2a0Sjeremylt   ierr = SNESComputeJacobianDefaultColor(snesdummy, X[0], matO[0],
490199551f5Sjeremylt                                          matcoarse, NULL); CHKERRQ(ierr);
49115ce0ef0Sjeremylt 
4926c5df90dSjeremylt   // Set up KSP
4936c5df90dSjeremylt   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
4946c5df90dSjeremylt   {
4956c5df90dSjeremylt     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
4966c5df90dSjeremylt     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
4976c5df90dSjeremylt     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
4986c5df90dSjeremylt                             PETSC_DEFAULT); CHKERRQ(ierr);
4996c5df90dSjeremylt   }
5006c5df90dSjeremylt   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
50161608365Sjeremylt   ierr = KSPSetOperators(ksp, matO[fineLevel], matO[fineLevel]);
5026c5df90dSjeremylt   CHKERRQ(ierr);
5036c5df90dSjeremylt 
5046c5df90dSjeremylt   // Set up PCMG
5056c5df90dSjeremylt   ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
5066c5df90dSjeremylt   PCMGCycleType pcgmcycletype = PC_MG_CYCLE_V;
5076c5df90dSjeremylt   {
5086c5df90dSjeremylt     ierr = PCSetType(pc, PCMG); CHKERRQ(ierr);
5096c5df90dSjeremylt 
5106c5df90dSjeremylt     // PCMG levels
5116c5df90dSjeremylt     ierr = PCMGSetLevels(pc, numlevels, NULL); CHKERRQ(ierr);
5126c5df90dSjeremylt     for (int i=0; i<numlevels; i++) {
5136c5df90dSjeremylt       // Smoother
5146c5df90dSjeremylt       KSP smoother;
5156c5df90dSjeremylt       PC smoother_pc;
5166c5df90dSjeremylt       ierr = PCMGGetSmoother(pc, i, &smoother); CHKERRQ(ierr);
5176c5df90dSjeremylt       ierr = KSPSetType(smoother, KSPCHEBYSHEV); CHKERRQ(ierr);
5186c5df90dSjeremylt       ierr = KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1); CHKERRQ(ierr);
5196c5df90dSjeremylt       ierr = KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE); CHKERRQ(ierr);
5206c5df90dSjeremylt       ierr = KSPSetOperators(smoother, matO[i], matO[i]); CHKERRQ(ierr);
5216c5df90dSjeremylt       ierr = KSPGetPC(smoother, &smoother_pc); CHKERRQ(ierr);
5226c5df90dSjeremylt       ierr = PCSetType(smoother_pc, PCJACOBI); CHKERRQ(ierr);
5236c5df90dSjeremylt       ierr = PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
5246c5df90dSjeremylt 
5256c5df90dSjeremylt       // Work vector
5266c5df90dSjeremylt       if (i < numlevels - 1) {
5276c5df90dSjeremylt         ierr = PCMGSetX(pc, i, X[i]); CHKERRQ(ierr);
5286c5df90dSjeremylt       }
5296c5df90dSjeremylt 
5306c5df90dSjeremylt       // Level transfers
5316c5df90dSjeremylt       if (i > 0) {
5326c5df90dSjeremylt         // Interpolation
533a97643b0Sjeremylt         ierr = PCMGSetInterpolation(pc, i, matPR[i]); CHKERRQ(ierr);
5346c5df90dSjeremylt       }
5356c5df90dSjeremylt 
5366c5df90dSjeremylt       // Coarse solve
5376c5df90dSjeremylt       KSP coarse;
5386c5df90dSjeremylt       PC coarse_pc;
5396c5df90dSjeremylt       ierr = PCMGGetCoarseSolve(pc, &coarse); CHKERRQ(ierr);
54015ce0ef0Sjeremylt       ierr = KSPSetType(coarse, KSPPREONLY); CHKERRQ(ierr);
541199551f5Sjeremylt       ierr = KSPSetOperators(coarse, matcoarse, matcoarse); CHKERRQ(ierr);
54215ce0ef0Sjeremylt 
5436c5df90dSjeremylt       ierr = KSPGetPC(coarse, &coarse_pc); CHKERRQ(ierr);
54415ce0ef0Sjeremylt       ierr = PCSetType(coarse_pc, PCGAMG); CHKERRQ(ierr);
54515ce0ef0Sjeremylt 
54615ce0ef0Sjeremylt       ierr = KSPSetOptionsPrefix(coarse, "coarse_"); CHKERRQ(ierr);
54715ce0ef0Sjeremylt       ierr = PCSetOptionsPrefix(coarse_pc, "coarse_"); CHKERRQ(ierr);
54815ce0ef0Sjeremylt       ierr = KSPSetFromOptions(coarse); CHKERRQ(ierr);
54915ce0ef0Sjeremylt       ierr = PCSetFromOptions(coarse_pc); CHKERRQ(ierr);
5506c5df90dSjeremylt     }
5516c5df90dSjeremylt 
5526c5df90dSjeremylt     // PCMG options
5536c5df90dSjeremylt     ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr);
5546c5df90dSjeremylt     ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr);
5556c5df90dSjeremylt     ierr = PCMGSetCycleType(pc, pcgmcycletype); CHKERRQ(ierr);
5566c5df90dSjeremylt   }
5576c5df90dSjeremylt 
5586c5df90dSjeremylt   // First run, if benchmarking
5596c5df90dSjeremylt   if (benchmark_mode) {
5606c5df90dSjeremylt     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
5616c5df90dSjeremylt     CHKERRQ(ierr);
56261608365Sjeremylt     ierr = VecZeroEntries(X[fineLevel]); CHKERRQ(ierr);
5636c5df90dSjeremylt     my_rt_start = MPI_Wtime();
56461608365Sjeremylt     ierr = KSPSolve(ksp, rhs, X[fineLevel]); CHKERRQ(ierr);
5656c5df90dSjeremylt     my_rt = MPI_Wtime() - my_rt_start;
5666c5df90dSjeremylt     ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
5676c5df90dSjeremylt     CHKERRQ(ierr);
5686c5df90dSjeremylt     // Set maxits based on first iteration timing
5696c5df90dSjeremylt     if (my_rt > 0.02) {
5706c5df90dSjeremylt       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
5716c5df90dSjeremylt       CHKERRQ(ierr);
5726c5df90dSjeremylt     } else {
5736c5df90dSjeremylt       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
5746c5df90dSjeremylt       CHKERRQ(ierr);
5756c5df90dSjeremylt     }
5766c5df90dSjeremylt   }
5776c5df90dSjeremylt 
5786c5df90dSjeremylt   // Timed solve
57961608365Sjeremylt   ierr = VecZeroEntries(X[fineLevel]); CHKERRQ(ierr);
5806c5df90dSjeremylt   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
58109a940d7Sjeremylt 
58209a940d7Sjeremylt   // -- Performance logging
58309a940d7Sjeremylt   ierr = PetscLogStageRegister("Solve Stage", &solvestage); CHKERRQ(ierr);
58409a940d7Sjeremylt   ierr = PetscLogStagePush(solvestage); CHKERRQ(ierr);
58509a940d7Sjeremylt 
58609a940d7Sjeremylt   // -- Solve
5876c5df90dSjeremylt   my_rt_start = MPI_Wtime();
58861608365Sjeremylt   ierr = KSPSolve(ksp, rhs, X[fineLevel]); CHKERRQ(ierr);
5896c5df90dSjeremylt   my_rt = MPI_Wtime() - my_rt_start;
5906c5df90dSjeremylt 
59109a940d7Sjeremylt 
59209a940d7Sjeremylt   // -- Performance logging
59309a940d7Sjeremylt   ierr = PetscLogStagePop();
59409a940d7Sjeremylt 
5956c5df90dSjeremylt   // Output results
5966c5df90dSjeremylt   {
5976c5df90dSjeremylt     KSPType ksptype;
5986c5df90dSjeremylt     PCMGType pcmgtype;
5996c5df90dSjeremylt     KSPConvergedReason reason;
6006c5df90dSjeremylt     PetscReal rnorm;
6016c5df90dSjeremylt     PetscInt its;
6026c5df90dSjeremylt     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
6036c5df90dSjeremylt     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
6046c5df90dSjeremylt     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
6056c5df90dSjeremylt     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
6066c5df90dSjeremylt     ierr = PCMGGetType(pc, &pcmgtype); CHKERRQ(ierr);
6076c5df90dSjeremylt     if (!test_mode || reason < 0 || rnorm > 1e-8) {
6086c5df90dSjeremylt       ierr = PetscPrintf(comm,
6096c5df90dSjeremylt                          "  KSP:\n"
6106c5df90dSjeremylt                          "    KSP Type                           : %s\n"
6116c5df90dSjeremylt                          "    KSP Convergence                    : %s\n"
6126c5df90dSjeremylt                          "    Total KSP Iterations               : %D\n"
6136c5df90dSjeremylt                          "    Final rnorm                        : %e\n",
6146c5df90dSjeremylt                          ksptype, KSPConvergedReasons[reason], its,
6156c5df90dSjeremylt                          (double)rnorm); CHKERRQ(ierr);
6166c5df90dSjeremylt       ierr = PetscPrintf(comm,
6176c5df90dSjeremylt                          "  PCMG:\n"
6186c5df90dSjeremylt                          "    PCMG Type                          : %s\n"
6196c5df90dSjeremylt                          "    PCMG Cycle Type                    : %s\n",
6206c5df90dSjeremylt                          PCMGTypes[pcmgtype],
6216c5df90dSjeremylt                          PCMGCycleTypes[pcgmcycletype]); CHKERRQ(ierr);
6226c5df90dSjeremylt     }
6236c5df90dSjeremylt     if (!test_mode) {
6246c5df90dSjeremylt       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
6256c5df90dSjeremylt     }
6266c5df90dSjeremylt     {
6276c5df90dSjeremylt       PetscReal maxerror;
628c70bd2a0Sjeremylt       ierr = ComputeErrorMax(userO[fineLevel], operror, X[fineLevel], target,
6296c5df90dSjeremylt                              &maxerror); CHKERRQ(ierr);
6306c5df90dSjeremylt       PetscReal tol = 5e-2;
6316c5df90dSjeremylt       if (!test_mode || maxerror > tol) {
6326c5df90dSjeremylt         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
6336c5df90dSjeremylt         CHKERRQ(ierr);
6346c5df90dSjeremylt         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
6356c5df90dSjeremylt         CHKERRQ(ierr);
6366c5df90dSjeremylt         ierr = PetscPrintf(comm,
6376c5df90dSjeremylt                            "    Pointwise Error (max)              : %e\n"
6386c5df90dSjeremylt                            "    CG Solve Time                      : %g (%g) sec\n",
6396c5df90dSjeremylt                            (double)maxerror, rt_max, rt_min); CHKERRQ(ierr);
6406c5df90dSjeremylt       }
6416c5df90dSjeremylt     }
6426c5df90dSjeremylt     if (benchmark_mode && (!test_mode)) {
6436c5df90dSjeremylt       ierr = PetscPrintf(comm,
6446c5df90dSjeremylt                          "    DoFs/Sec in CG                     : %g (%g) million\n",
64561608365Sjeremylt                          1e-6*gsize[fineLevel]*its/rt_max,
64661608365Sjeremylt                          1e-6*gsize[fineLevel]*its/rt_min);
6476c5df90dSjeremylt       CHKERRQ(ierr);
6486c5df90dSjeremylt     }
6496c5df90dSjeremylt   }
6506c5df90dSjeremylt 
6516c5df90dSjeremylt   if (write_solution) {
6526c5df90dSjeremylt     PetscViewer vtkviewersoln;
6536c5df90dSjeremylt 
6546c5df90dSjeremylt     ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr);
6556c5df90dSjeremylt     ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr);
656cfa59c5bSRey     ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtu"); CHKERRQ(ierr);
65761608365Sjeremylt     ierr = VecView(X[fineLevel], vtkviewersoln); CHKERRQ(ierr);
6586c5df90dSjeremylt     ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr);
6596c5df90dSjeremylt   }
6606c5df90dSjeremylt 
6616c5df90dSjeremylt   // Cleanup
6626c5df90dSjeremylt   for (int i=0; i<numlevels; i++) {
6636c5df90dSjeremylt     ierr = VecDestroy(&X[i]); CHKERRQ(ierr);
6646c5df90dSjeremylt     ierr = VecDestroy(&Xloc[i]); CHKERRQ(ierr);
6656c5df90dSjeremylt     ierr = VecDestroy(&mult[i]); CHKERRQ(ierr);
6666c5df90dSjeremylt     ierr = VecDestroy(&userO[i]->Yloc); CHKERRQ(ierr);
6676c5df90dSjeremylt     ierr = MatDestroy(&matO[i]); CHKERRQ(ierr);
6686c5df90dSjeremylt     ierr = PetscFree(userO[i]); CHKERRQ(ierr);
6696c5df90dSjeremylt     if (i > 0) {
670a97643b0Sjeremylt       ierr = MatDestroy(&matPR[i]); CHKERRQ(ierr);
671a97643b0Sjeremylt       ierr = PetscFree(userPR[i]); CHKERRQ(ierr);
6726c5df90dSjeremylt     }
6736c5df90dSjeremylt     ierr = CeedDataDestroy(i, ceeddata[i]); CHKERRQ(ierr);
6746c5df90dSjeremylt     ierr = DMDestroy(&dm[i]); CHKERRQ(ierr);
6756c5df90dSjeremylt   }
6766c5df90dSjeremylt   ierr = PetscFree(leveldegrees); CHKERRQ(ierr);
6776c5df90dSjeremylt   ierr = PetscFree(dm); CHKERRQ(ierr);
6786c5df90dSjeremylt   ierr = PetscFree(X); CHKERRQ(ierr);
6796c5df90dSjeremylt   ierr = PetscFree(Xloc); CHKERRQ(ierr);
6806c5df90dSjeremylt   ierr = PetscFree(mult); CHKERRQ(ierr);
6816c5df90dSjeremylt   ierr = PetscFree(matO); CHKERRQ(ierr);
682a97643b0Sjeremylt   ierr = PetscFree(matPR); CHKERRQ(ierr);
6836c5df90dSjeremylt   ierr = PetscFree(ceeddata); CHKERRQ(ierr);
6846c5df90dSjeremylt   ierr = PetscFree(userO); CHKERRQ(ierr);
685a97643b0Sjeremylt   ierr = PetscFree(userPR); CHKERRQ(ierr);
6866c5df90dSjeremylt   ierr = PetscFree(lsize); CHKERRQ(ierr);
6876c5df90dSjeremylt   ierr = PetscFree(xlsize); CHKERRQ(ierr);
6886c5df90dSjeremylt   ierr = PetscFree(gsize); CHKERRQ(ierr);
6896c5df90dSjeremylt   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
6906c5df90dSjeremylt   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
691199551f5Sjeremylt   ierr = MatDestroy(&matcoarse); CHKERRQ(ierr);
6926c5df90dSjeremylt   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
693c70bd2a0Sjeremylt   ierr = SNESDestroy(&snesdummy); CHKERRQ(ierr);
694199551f5Sjeremylt   ierr = DMDestroy(&dmorig); CHKERRQ(ierr);
6956c5df90dSjeremylt   CeedVectorDestroy(&target);
696c70bd2a0Sjeremylt   CeedQFunctionDestroy(&qferror);
697c70bd2a0Sjeremylt   CeedQFunctionDestroy(&qfrestrict);
698c70bd2a0Sjeremylt   CeedQFunctionDestroy(&qfprolong);
699c70bd2a0Sjeremylt   CeedOperatorDestroy(&operror);
7006c5df90dSjeremylt   CeedDestroy(&ceed);
7016c5df90dSjeremylt   return PetscFinalize();
7026c5df90dSjeremylt }
703