xref: /libCEED/examples/petsc/multigrid.c (revision 483f8b0d18b49624d6de41b1f2037a7fe5b66ef5)
1 // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3 // reserved. See files LICENSE and NOTICE for details.
4 //
5 // This file is part of CEED, a collection of benchmarks, miniapps, software
6 // libraries and APIs for efficient high-order finite element and spectral
7 // element discretizations for exascale applications. For more information and
8 // source code availability see http://github.com/ceed.
9 //
10 // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11 // a collaborative effort of two U.S. Department of Energy organizations (Office
12 // of Science and the National Nuclear Security Administration) responsible for
13 // the planning and preparation of a capable exascale ecosystem, including
14 // software, applications, hardware, advanced system engineering and early
15 // testbed platforms, in support of the nation's exascale computing imperative.
16 
17 //                        libCEED + PETSc Example: CEED BPs 3-6 with Multigrid
18 //
19 // This example demonstrates a simple usage of libCEED with PETSc to solve the
20 // CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps.
21 //
22 // The code uses higher level communication protocols in DMPlex.
23 //
24 // Build with:
25 //
26 //     make multigrid [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
27 //
28 // Sample runs:
29 //
30 //     multigrid -problem bp3
31 //     multigrid -problem bp4 -ceed /cpu/self
32 //     multigrid -problem bp5 -ceed /cpu/occa
33 //     multigrid -problem bp6 -ceed /gpu/cuda
34 //
35 //TESTARGS -ceed {ceed_resource} -test -problem bp3 -degree 3
36 
37 /// @file
38 /// CEED BPs 1-6 multigrid example using PETSc
39 const char help[] = "Solve CEED BPs using p-multigrid with PETSc and DMPlex\n";
40 
41 #define multigrid
42 #include "setup.h"
43 
44 int main(int argc, char **argv) {
45   PetscInt ierr;
46   MPI_Comm comm;
47   char filename[PETSC_MAX_PATH_LEN],
48        ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self";
49   double my_rt_start, my_rt, rt_min, rt_max;
50   PetscInt degree = 3, qextra, *lsize, *xlsize, *gsize, dim = 3,
51            melem[3] = {3, 3, 3}, ncompu = 1, numlevels = degree, *leveldegrees;
52   PetscScalar *r;
53   PetscBool test_mode, benchmark_mode, read_mesh, write_solution;
54   DM  *dm, dmOrig;
55   SNES snes_dummy;
56   KSP ksp;
57   PC pc;
58   Mat *matO, *matI, *matR, matCoarse;
59   Vec *X, *Xloc, *mult, rhs, rhsloc;
60   UserO *userO;
61   UserIR *userI, *userR;
62   Ceed ceed;
63   CeedData *ceeddata;
64   CeedVector rhsceed, diagceed, target;
65   CeedQFunction qf_error, qf_restrict, qf_prolong;
66   CeedOperator op_error;
67   bpType bpChoice;
68   coarsenType coarsen;
69 
70   ierr = PetscInitialize(&argc, &argv, NULL, help);
71   if (ierr) return ierr;
72   comm = PETSC_COMM_WORLD;
73 
74   // Parse command line options
75   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr);
76   bpChoice = CEED_BP3;
77   ierr = PetscOptionsEnum("-problem",
78                           "CEED benchmark problem to solve", NULL,
79                           bpTypes, (PetscEnum)bpChoice, (PetscEnum *)&bpChoice,
80                           NULL); CHKERRQ(ierr);
81   ncompu = bpOptions[bpChoice].ncompu;
82   test_mode = PETSC_FALSE;
83   ierr = PetscOptionsBool("-test",
84                           "Testing mode (do not print unless error is large)",
85                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
86   benchmark_mode = PETSC_FALSE;
87   ierr = PetscOptionsBool("-benchmark",
88                           "Benchmarking mode (prints benchmark statistics)",
89                           NULL, benchmark_mode, &benchmark_mode, NULL);
90   CHKERRQ(ierr);
91   write_solution = PETSC_FALSE;
92   ierr = PetscOptionsBool("-write_solution",
93                           "Write solution for visualization",
94                           NULL, write_solution, &write_solution, NULL);
95   CHKERRQ(ierr);
96   degree = test_mode ? 3 : 2;
97   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
98                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
99   if (degree < 1) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
100                              "-degree %D must be at least 1", degree);
101   qextra = bpOptions[bpChoice].qextra;
102   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
103                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
104   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
105                             NULL, ceedresource, ceedresource,
106                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
107   coarsen = COARSEN_UNIFORM;
108   ierr = PetscOptionsEnum("-coarsen",
109                           "Coarsening strategy to use", NULL,
110                           coarsenTypes, (PetscEnum)coarsen,
111                           (PetscEnum *)&coarsen, NULL); CHKERRQ(ierr);
112   read_mesh = PETSC_FALSE;
113   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL,
114                             filename, filename, sizeof(filename), &read_mesh);
115   CHKERRQ(ierr);
116   if (!read_mesh) {
117     PetscInt tmp = dim;
118     ierr = PetscOptionsIntArray("-cells","Number of cells per dimension", NULL,
119                                 melem, &tmp, NULL); CHKERRQ(ierr);
120   }
121   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
122 
123   // Setup DM
124   if (read_mesh) {
125     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dmOrig);
126     CHKERRQ(ierr);
127   } else {
128     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_FALSE, melem, NULL,
129                                NULL, NULL, PETSC_TRUE,&dmOrig); CHKERRQ(ierr);
130   }
131 
132   {
133     DM dmDist = NULL;
134     PetscPartitioner part;
135 
136     ierr = DMPlexGetPartitioner(dmOrig, &part); CHKERRQ(ierr);
137     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
138     ierr = DMPlexDistribute(dmOrig, 0, NULL, &dmDist); CHKERRQ(ierr);
139     if (dmDist) {
140       ierr = DMDestroy(&dmOrig); CHKERRQ(ierr);
141       dmOrig = dmDist;
142     }
143   }
144 
145   // Allocate arrays for PETSc objects for each level
146   switch (coarsen) {
147   case COARSEN_UNIFORM:
148     numlevels = degree;
149     break;
150   case COARSEN_LOGARITHMIC:
151     numlevels = ceil(log(degree)/log(2)) + 1;
152     break;
153   }
154   ierr = PetscMalloc1(numlevels, &leveldegrees); CHKERRQ(ierr);
155   switch (coarsen) {
156   case COARSEN_UNIFORM:
157     for (int i=0; i<numlevels; i++) leveldegrees[i] = i + 1;
158     break;
159   case COARSEN_LOGARITHMIC:
160     for (int i=0; i<numlevels-1; i++) leveldegrees[i] = pow(2,i);
161     leveldegrees[numlevels-1] = degree;
162     break;
163   }
164   ierr = PetscMalloc1(numlevels, &dm); CHKERRQ(ierr);
165   ierr = PetscMalloc1(numlevels, &X); CHKERRQ(ierr);
166   ierr = PetscMalloc1(numlevels, &Xloc); CHKERRQ(ierr);
167   ierr = PetscMalloc1(numlevels, &mult); CHKERRQ(ierr);
168   ierr = PetscMalloc1(numlevels, &userO); CHKERRQ(ierr);
169   ierr = PetscMalloc1(numlevels, &userI); CHKERRQ(ierr);
170   ierr = PetscMalloc1(numlevels, &userR); CHKERRQ(ierr);
171   ierr = PetscMalloc1(numlevels, &matO); CHKERRQ(ierr);
172   ierr = PetscMalloc1(numlevels, &matI); CHKERRQ(ierr);
173   ierr = PetscMalloc1(numlevels, &matR); CHKERRQ(ierr);
174   ierr = PetscMalloc1(numlevels, &lsize); CHKERRQ(ierr);
175   ierr = PetscMalloc1(numlevels, &xlsize); CHKERRQ(ierr);
176   ierr = PetscMalloc1(numlevels, &gsize); CHKERRQ(ierr);
177 
178   // Setup DM and Operator Mat Shells for each level
179   for (CeedInt i=0; i<numlevels; i++) {
180     // Create DM
181     ierr = DMClone(dmOrig, &dm[i]); CHKERRQ(ierr);
182     ierr = SetupDMByDegree(dm[i], leveldegrees[i], ncompu, bpChoice);
183     CHKERRQ(ierr);
184 
185     // Create vectors
186     ierr = DMCreateGlobalVector(dm[i], &X[i]); CHKERRQ(ierr);
187     ierr = VecGetLocalSize(X[i], &lsize[i]); CHKERRQ(ierr);
188     ierr = VecGetSize(X[i], &gsize[i]); CHKERRQ(ierr);
189     ierr = DMCreateLocalVector(dm[i], &Xloc[i]); CHKERRQ(ierr);
190     ierr = VecGetSize(Xloc[i], &xlsize[i]); CHKERRQ(ierr);
191 
192     // Operator
193     ierr = PetscMalloc1(1, &userO[i]); CHKERRQ(ierr);
194     ierr = MatCreateShell(comm, lsize[i], lsize[i], gsize[i], gsize[i],
195                           userO[i], &matO[i]); CHKERRQ(ierr);
196     ierr = MatShellSetOperation(matO[i], MATOP_MULT,
197                                 (void(*)(void))MatMult_Ceed);
198     ierr = MatShellSetOperation(matO[i], MATOP_GET_DIAGONAL,
199                                 (void(*)(void))MatGetDiag);
200     CHKERRQ(ierr);
201 
202     // Level transfers
203     if (i > 0) {
204       // Interp
205       ierr = PetscMalloc1(1, &userI[i]); CHKERRQ(ierr);
206       ierr = MatCreateShell(comm, lsize[i], lsize[i-1], gsize[i], gsize[i-1],
207                             userI[i], &matI[i]); CHKERRQ(ierr);
208       ierr = MatShellSetOperation(matI[i], MATOP_MULT,
209                                   (void(*)(void))MatMult_Interp);
210       CHKERRQ(ierr);
211 
212       // Restrict
213       ierr = PetscMalloc1(1, &userR[i]); CHKERRQ(ierr);
214       ierr = MatCreateShell(comm, lsize[i-1], lsize[i], gsize[i-1], gsize[i],
215                             userR[i], &matR[i]); CHKERRQ(ierr);
216       ierr = MatShellSetOperation(matR[i], MATOP_MULT,
217                                   (void(*)(void))MatMult_Restrict);
218       CHKERRQ(ierr);
219     }
220   }
221   ierr = VecDuplicate(X[numlevels-1], &rhs); CHKERRQ(ierr);
222 
223   // Set up libCEED
224   CeedInit(ceedresource, &ceed);
225 
226   // Print global grid information
227   if (!test_mode) {
228     PetscInt P = degree + 1, Q = P + qextra;
229     const char *usedresource;
230     CeedGetResource(ceed, &usedresource);
231     ierr = PetscPrintf(comm,
232                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc + PCMG --\n"
233                        "  libCEED:\n"
234                        "    libCEED Backend                    : %s\n"
235                        "  Mesh:\n"
236                        "    Number of 1D Basis Nodes (p)       : %d\n"
237                        "    Number of 1D Quadrature Points (q) : %d\n"
238                        "    Global Nodes                       : %D\n"
239                        "    Owned Nodes                        : %D\n"
240                        "    DoF per node                       : %D\n"
241                        "  Multigrid:\n"
242                        "    Number of Levels                   : %d\n",
243                        bpChoice+1, usedresource, P, Q,
244                        gsize[numlevels-1]/ncompu, lsize[numlevels-1]/ncompu,
245                        ncompu, numlevels); CHKERRQ(ierr);
246   }
247 
248   // Create RHS vector
249   ierr = VecDuplicate(Xloc[numlevels-1], &rhsloc); CHKERRQ(ierr);
250   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
251   ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr);
252   CeedVectorCreate(ceed, xlsize[numlevels-1], &rhsceed);
253   CeedVectorSetArray(rhsceed, CEED_MEM_HOST, CEED_USE_POINTER, r);
254 
255   // Set up libCEED operators on each level
256   ierr = PetscMalloc1(numlevels, &ceeddata); CHKERRQ(ierr);
257   for (int i=0; i<numlevels; i++) {
258     // Print level information
259     if (!test_mode && (i == 0 || i == numlevels-1)) {
260       ierr = PetscPrintf(comm,"    Level %D (%s):\n"
261                          "      Number of 1D Basis Nodes (p)     : %d\n"
262                          "      Global Nodes                     : %D\n"
263                          "      Owned Nodes                      : %D\n",
264                          i, (i? "fine" : "coarse"), leveldegrees[i] + 1,
265                          gsize[i]/ncompu, lsize[i]/ncompu); CHKERRQ(ierr);
266     }
267     ierr = PetscMalloc1(1, &ceeddata[i]); CHKERRQ(ierr);
268     ierr = SetupLibceedByDegree(dm[i], ceed, leveldegrees[i], dim, qextra,
269                                 ncompu, gsize[i], xlsize[i], bpChoice,
270                                 ceeddata[i], i==(numlevels-1), rhsceed,
271                                 &target); CHKERRQ(ierr);
272   }
273 
274   // Gather RHS
275   ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr);
276   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
277   ierr = DMLocalToGlobal(dm[numlevels-1], rhsloc, ADD_VALUES, rhs);
278   CHKERRQ(ierr);
279   CeedVectorDestroy(&rhsceed);
280 
281   // Create the restriction/interpolation Q-function
282   CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_NONE, CEED_EVAL_INTERP,
283                               &qf_restrict);
284   CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_INTERP, CEED_EVAL_NONE,
285                               &qf_prolong);
286 
287   // Set up libCEED level transfer operators
288   ierr = CeedLevelTransferSetup(ceed, numlevels, ncompu, bpChoice, ceeddata,
289                                 leveldegrees, qf_restrict, qf_prolong);
290   CHKERRQ(ierr);
291 
292   // Create the error Q-function
293   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].error,
294                               bpOptions[bpChoice].errorfname, &qf_error);
295   CeedQFunctionAddInput(qf_error, "u", ncompu, CEED_EVAL_INTERP);
296   CeedQFunctionAddInput(qf_error, "true_soln", ncompu, CEED_EVAL_NONE);
297   CeedQFunctionAddOutput(qf_error, "error", ncompu, CEED_EVAL_NONE);
298 
299   // Create the error operator
300   CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
301                      &op_error);
302   CeedOperatorSetField(op_error, "u", ceeddata[numlevels-1]->Erestrictu,
303                        ceeddata[numlevels-1]->basisu, CEED_VECTOR_ACTIVE);
304   CeedOperatorSetField(op_error, "true_soln", ceeddata[numlevels-1]->Erestrictui,
305                        CEED_BASIS_COLLOCATED, target);
306   CeedOperatorSetField(op_error, "error", ceeddata[numlevels-1]->Erestrictui,
307                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
308 
309   // Calculate multiplicity
310   for (int i=0; i<numlevels; i++) {
311     PetscScalar *x;
312 
313     // CEED vector
314     ierr = VecGetArray(Xloc[i], &x); CHKERRQ(ierr);
315     CeedVectorSetArray(ceeddata[i]->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
316 
317     // Multiplicity
318     CeedElemRestrictionGetMultiplicity(ceeddata[i]->Erestrictu,
319                                        ceeddata[i]->xceed);
320 
321     // Restore vector
322     ierr = VecRestoreArray(Xloc[i], &x); CHKERRQ(ierr);
323 
324     // Creat mult vector
325     ierr = VecDuplicate(Xloc[i], &mult[i]); CHKERRQ(ierr);
326 
327     // Local-to-global
328     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
329     ierr = DMLocalToGlobal(dm[i], Xloc[i], ADD_VALUES, X[i]);
330     CHKERRQ(ierr);
331     ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr);
332 
333     // Global-to-local
334     ierr = DMGlobalToLocal(dm[i], X[i], INSERT_VALUES, mult[i]);
335     CHKERRQ(ierr);
336     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
337 
338     // Multiplicity scaling
339     ierr = VecReciprocal(mult[i]);
340   }
341 
342   // Set up Mat
343   for (int i=0; i<numlevels; i++) {
344     // User Operator
345     userO[i]->comm = comm;
346     userO[i]->dm = dm[i];
347     userO[i]->Xloc = Xloc[i];
348     ierr = VecDuplicate(Xloc[i], &userO[i]->Yloc); CHKERRQ(ierr);
349     userO[i]->xceed = ceeddata[i]->xceed;
350     userO[i]->yceed = ceeddata[i]->yceed;
351     userO[i]->op = ceeddata[i]->op_apply;
352     userO[i]->ceed = ceed;
353 
354     // Set up diagonal
355     const CeedScalar *ceedarray;
356     ierr = VecDuplicate(X[i], &userO[i]->diag); CHKERRQ(ierr);
357 
358     // -- Local diagonal
359     CeedOperatorAssembleLinearDiagonal(userO[i]->op, &diagceed,
360                                        CEED_REQUEST_IMMEDIATE);
361 
362     // -- Set PETSc array
363     CeedVectorGetArrayRead(diagceed, CEED_MEM_HOST, &ceedarray);
364     ierr = VecPlaceArray(Xloc[i], ceedarray); CHKERRQ(ierr);
365     CeedVectorRestoreArrayRead(diagceed, &ceedarray);
366 
367     // -- Global diagonal
368     ierr = VecZeroEntries(userO[i]->diag); CHKERRQ(ierr);
369     ierr = DMLocalToGlobal(userO[i]->dm, Xloc[i], ADD_VALUES,
370                                 userO[i]->diag); CHKERRQ(ierr);
371 
372     // -- Cleanup
373     ierr = VecResetArray(Xloc[i]); CHKERRQ(ierr);
374     CeedVectorDestroy(&diagceed);
375 
376     if (i > 0) {
377       // Interp Operator
378       userI[i]->comm = comm;
379       userI[i]->dmc = dm[i-1];
380       userI[i]->dmf = dm[i];
381       userI[i]->Xloc = Xloc[i-1];
382       userI[i]->Yloc = userO[i]->Yloc;
383       userI[i]->mult = mult[i];
384       userI[i]->ceedvecc = userO[i-1]->xceed;
385       userI[i]->ceedvecf = userO[i]->yceed;
386       userI[i]->op = ceeddata[i]->op_interp;
387       userI[i]->ceed = ceed;
388 
389       // Restrict Operator
390       userR[i]->comm = comm;
391       userR[i]->dmc = dm[i-1];
392       userR[i]->dmf = dm[i];
393       userR[i]->Xloc = Xloc[i];
394       userR[i]->Yloc = userO[i-1]->Yloc;
395       userR[i]->mult = mult[i];
396       userR[i]->ceedvecf = userO[i]->xceed;
397       userR[i]->ceedvecc = userO[i-1]->yceed;
398       userR[i]->op = ceeddata[i]->op_restrict;
399       userR[i]->ceed = ceed;
400     }
401   }
402 
403   // Setup dummy SNES for AMG coarse solve
404   ierr = SNESCreate(comm, &snes_dummy); CHKERRQ(ierr);
405   ierr = SNESSetDM(snes_dummy, dm[0]); CHKERRQ(ierr);
406   ierr = SNESSetSolution(snes_dummy, X[0]); CHKERRQ(ierr);
407 
408   // -- Jacobian matrix
409   ierr = DMSetMatType(dm[0], MATAIJ); CHKERRQ(ierr);
410   ierr = DMCreateMatrix(dm[0], &matCoarse); CHKERRQ(ierr);
411   ierr = SNESSetJacobian(snes_dummy, matCoarse, matCoarse, NULL,
412                          NULL); CHKERRQ(ierr);
413 
414   // -- Residual evaluation function
415   ierr = SNESSetFunction(snes_dummy, X[0], FormResidual_Ceed,
416                          userO[0]); CHKERRQ(ierr);
417 
418   // -- Form Jacobian
419   ierr = SNESComputeJacobianDefaultColor(snes_dummy, X[0], matO[0],
420                                          matCoarse, NULL); CHKERRQ(ierr);
421 
422   // Set up KSP
423   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
424   {
425     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
426     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
427     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
428                             PETSC_DEFAULT); CHKERRQ(ierr);
429   }
430   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
431   ierr = KSPSetOperators(ksp, matO[numlevels-1], matO[numlevels-1]);
432   CHKERRQ(ierr);
433 
434   // Set up PCMG
435   ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
436   PCMGCycleType pcgmcycletype = PC_MG_CYCLE_V;
437   {
438     ierr = PCSetType(pc, PCMG); CHKERRQ(ierr);
439 
440     // PCMG levels
441     ierr = PCMGSetLevels(pc, numlevels, NULL); CHKERRQ(ierr);
442     for (int i=0; i<numlevels; i++) {
443       // Smoother
444       KSP smoother;
445       PC smoother_pc;
446       ierr = PCMGGetSmoother(pc, i, &smoother); CHKERRQ(ierr);
447       ierr = KSPSetType(smoother, KSPCHEBYSHEV); CHKERRQ(ierr);
448       ierr = KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1); CHKERRQ(ierr);
449       ierr = KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE); CHKERRQ(ierr);
450       ierr = KSPSetOperators(smoother, matO[i], matO[i]); CHKERRQ(ierr);
451       ierr = KSPGetPC(smoother, &smoother_pc); CHKERRQ(ierr);
452       ierr = PCSetType(smoother_pc, PCJACOBI); CHKERRQ(ierr);
453       ierr = PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
454 
455       // Work vector
456       if (i < numlevels-1) {
457         ierr = PCMGSetX(pc, i, X[i]); CHKERRQ(ierr);
458       }
459 
460       // Level transfers
461       if (i > 0) {
462         // Interpolation
463         ierr = PCMGSetInterpolation(pc, i, matI[i]); CHKERRQ(ierr);
464 
465         // Restriction
466         ierr = PCMGSetRestriction(pc, i, matR[i]); CHKERRQ(ierr);
467       }
468 
469       // Coarse solve
470       KSP coarse;
471       PC coarse_pc;
472       ierr = PCMGGetCoarseSolve(pc, &coarse); CHKERRQ(ierr);
473       ierr = KSPSetType(coarse, KSPPREONLY); CHKERRQ(ierr);
474       ierr = KSPSetOperators(coarse, matCoarse, matCoarse); CHKERRQ(ierr);
475 
476       ierr = KSPGetPC(coarse, &coarse_pc); CHKERRQ(ierr);
477       ierr = PCSetType(coarse_pc, PCGAMG); CHKERRQ(ierr);
478 
479       ierr = KSPSetOptionsPrefix(coarse, "coarse_"); CHKERRQ(ierr);
480       ierr = PCSetOptionsPrefix(coarse_pc, "coarse_"); CHKERRQ(ierr);
481       ierr = KSPSetFromOptions(coarse); CHKERRQ(ierr);
482       ierr = PCSetFromOptions(coarse_pc); CHKERRQ(ierr);
483     }
484 
485     // PCMG options
486     ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr);
487     ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr);
488     ierr = PCMGSetCycleType(pc, pcgmcycletype); CHKERRQ(ierr);
489   }
490 
491   // First run, if benchmarking
492   if (benchmark_mode) {
493     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
494     CHKERRQ(ierr);
495     ierr = VecZeroEntries(X[numlevels-1]); CHKERRQ(ierr);
496     my_rt_start = MPI_Wtime();
497     ierr = KSPSolve(ksp, rhs, X[numlevels-1]); CHKERRQ(ierr);
498     my_rt = MPI_Wtime() - my_rt_start;
499     ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
500     CHKERRQ(ierr);
501     // Set maxits based on first iteration timing
502     if (my_rt > 0.02) {
503       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
504       CHKERRQ(ierr);
505     } else {
506       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
507       CHKERRQ(ierr);
508     }
509   }
510 
511   // Timed solve
512   ierr = VecZeroEntries(X[numlevels-1]); CHKERRQ(ierr);
513   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
514   my_rt_start = MPI_Wtime();
515   ierr = KSPSolve(ksp, rhs, X[numlevels-1]); CHKERRQ(ierr);
516   my_rt = MPI_Wtime() - my_rt_start;
517 
518   // Output results
519   {
520     KSPType ksptype;
521     PCMGType pcmgtype;
522     KSPConvergedReason reason;
523     PetscReal rnorm;
524     PetscInt its;
525     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
526     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
527     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
528     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
529     ierr = PCMGGetType(pc, &pcmgtype); CHKERRQ(ierr);
530     if (!test_mode || reason < 0 || rnorm > 1e-8) {
531       ierr = PetscPrintf(comm,
532                          "  KSP:\n"
533                          "    KSP Type                           : %s\n"
534                          "    KSP Convergence                    : %s\n"
535                          "    Total KSP Iterations               : %D\n"
536                          "    Final rnorm                        : %e\n",
537                          ksptype, KSPConvergedReasons[reason], its,
538                          (double)rnorm); CHKERRQ(ierr);
539       ierr = PetscPrintf(comm,
540                          "  PCMG:\n"
541                          "    PCMG Type                          : %s\n"
542                          "    PCMG Cycle Type                    : %s\n",
543                          PCMGTypes[pcmgtype],
544                          PCMGCycleTypes[pcgmcycletype]); CHKERRQ(ierr);
545     }
546     if (!test_mode) {
547       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
548     }
549     {
550       PetscReal maxerror;
551       ierr = ComputeErrorMax(userO[numlevels-1], op_error, X[numlevels-1], target,
552                              &maxerror); CHKERRQ(ierr);
553       PetscReal tol = 5e-2;
554       if (!test_mode || maxerror > tol) {
555         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
556         CHKERRQ(ierr);
557         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
558         CHKERRQ(ierr);
559         ierr = PetscPrintf(comm,
560                            "    Pointwise Error (max)              : %e\n"
561                            "    CG Solve Time                      : %g (%g) sec\n",
562                            (double)maxerror, rt_max, rt_min); CHKERRQ(ierr);
563       }
564     }
565     if (benchmark_mode && (!test_mode)) {
566       ierr = PetscPrintf(comm,
567                          "    DoFs/Sec in CG                     : %g (%g) million\n",
568                          1e-6*gsize[numlevels-1]*its/rt_max,
569                          1e-6*gsize[numlevels-1]*its/rt_min);
570       CHKERRQ(ierr);
571     }
572   }
573 
574   if (write_solution) {
575     PetscViewer vtkviewersoln;
576 
577     ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr);
578     ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr);
579     ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr);
580     ierr = VecView(X[numlevels-1], vtkviewersoln); CHKERRQ(ierr);
581     ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr);
582   }
583 
584   // Cleanup
585   for (int i=0; i<numlevels; i++) {
586     ierr = VecDestroy(&X[i]); CHKERRQ(ierr);
587     ierr = VecDestroy(&Xloc[i]); CHKERRQ(ierr);
588     ierr = VecDestroy(&mult[i]); CHKERRQ(ierr);
589     ierr = VecDestroy(&userO[i]->Yloc); CHKERRQ(ierr);
590     ierr = VecDestroy(&userO[i]->diag); CHKERRQ(ierr);
591     ierr = MatDestroy(&matO[i]); CHKERRQ(ierr);
592     ierr = PetscFree(userO[i]); CHKERRQ(ierr);
593     if (i > 0) {
594       ierr = MatDestroy(&matI[i]); CHKERRQ(ierr);
595       ierr = PetscFree(userI[i]); CHKERRQ(ierr);
596       ierr = MatDestroy(&matR[i]); CHKERRQ(ierr);
597       ierr = PetscFree(userR[i]); CHKERRQ(ierr);
598     }
599     ierr = CeedDataDestroy(i, ceeddata[i]); CHKERRQ(ierr);
600     ierr = DMDestroy(&dm[i]); CHKERRQ(ierr);
601   }
602   ierr = PetscFree(leveldegrees); CHKERRQ(ierr);
603   ierr = PetscFree(dm); CHKERRQ(ierr);
604   ierr = PetscFree(X); CHKERRQ(ierr);
605   ierr = PetscFree(Xloc); CHKERRQ(ierr);
606   ierr = PetscFree(mult); CHKERRQ(ierr);
607   ierr = PetscFree(matO); CHKERRQ(ierr);
608   ierr = PetscFree(matI); CHKERRQ(ierr);
609   ierr = PetscFree(matR); CHKERRQ(ierr);
610   ierr = PetscFree(ceeddata); CHKERRQ(ierr);
611   ierr = PetscFree(userO); CHKERRQ(ierr);
612   ierr = PetscFree(userI); CHKERRQ(ierr);
613   ierr = PetscFree(userR); CHKERRQ(ierr);
614   ierr = PetscFree(lsize); CHKERRQ(ierr);
615   ierr = PetscFree(xlsize); CHKERRQ(ierr);
616   ierr = PetscFree(gsize); CHKERRQ(ierr);
617   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
618   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
619   ierr = MatDestroy(&matCoarse); CHKERRQ(ierr);
620   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
621   ierr = SNESDestroy(&snes_dummy); CHKERRQ(ierr);
622   ierr = DMDestroy(&dmOrig); CHKERRQ(ierr);
623   CeedVectorDestroy(&target);
624   CeedQFunctionDestroy(&qf_error);
625   CeedQFunctionDestroy(&qf_restrict);
626   CeedQFunctionDestroy(&qf_prolong);
627   CeedOperatorDestroy(&op_error);
628   CeedDestroy(&ceed);
629   return PetscFinalize();
630 }
631