xref: /libCEED/examples/petsc/multigrid.c (revision 16911fdad6ca4b1dd37f9d3206958ee664667dbe)
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   KSP ksp;
56   PC pc;
57   Mat *matO, *matI, *matR;
58   Vec *X, *Xloc, *mult, rhs, rhsloc, diagloc;
59   UserO *userO;
60   UserIR *userI, *userR;
61   Ceed ceed;
62   CeedData *ceeddata;
63   CeedVector rhsceed, diagceed, target;
64   CeedQFunction qf_error, qf_restrict, qf_prolong;
65   CeedOperator op_error;
66   bpType bpChoice;
67   coarsenType coarsen;
68 
69   ierr = PetscInitialize(&argc, &argv, NULL, help);
70   if (ierr) return ierr;
71   comm = PETSC_COMM_WORLD;
72 
73   // Parse command line options
74   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr);
75   bpChoice = CEED_BP3;
76   ierr = PetscOptionsEnum("-problem",
77                           "CEED benchmark problem to solve", NULL,
78                           bpTypes, (PetscEnum)bpChoice, (PetscEnum *)&bpChoice,
79                           NULL); CHKERRQ(ierr);
80   ncompu = bpOptions[bpChoice].ncompu;
81   test_mode = PETSC_FALSE;
82   ierr = PetscOptionsBool("-test",
83                           "Testing mode (do not print unless error is large)",
84                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
85   benchmark_mode = PETSC_FALSE;
86   ierr = PetscOptionsBool("-benchmark",
87                           "Benchmarking mode (prints benchmark statistics)",
88                           NULL, benchmark_mode, &benchmark_mode, NULL);
89   CHKERRQ(ierr);
90   write_solution = PETSC_FALSE;
91   ierr = PetscOptionsBool("-write_solution",
92                           "Write solution for visualization",
93                           NULL, write_solution, &write_solution, NULL);
94   CHKERRQ(ierr);
95   degree = test_mode ? 3 : 2;
96   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
97                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
98   if (degree < 1) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
99                              "-degree %D must be at least 1", degree);
100   qextra = bpOptions[bpChoice].qextra;
101   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
102                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
103   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
104                             NULL, ceedresource, ceedresource,
105                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
106   coarsen = COARSEN_UNIFORM;
107   ierr = PetscOptionsEnum("-coarsen",
108                           "Coarsening strategy to use", NULL,
109                           coarsenTypes, (PetscEnum)coarsen,
110                           (PetscEnum *)&coarsen, NULL); CHKERRQ(ierr);
111   read_mesh = PETSC_FALSE;
112   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL,
113                             filename, filename, sizeof(filename), &read_mesh);
114   CHKERRQ(ierr);
115   if (!read_mesh) {
116     PetscInt tmp = dim;
117     ierr = PetscOptionsIntArray("-cells","Number of cells per dimension", NULL,
118                                 melem, &tmp, NULL); CHKERRQ(ierr);
119   }
120   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
121 
122   // Setup DM
123   if (read_mesh) {
124     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dmOrig);
125     CHKERRQ(ierr);
126   } else {
127     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_FALSE, melem, NULL,
128                                NULL, NULL, PETSC_TRUE,&dmOrig); CHKERRQ(ierr);
129   }
130 
131   {
132     DM dmDist = NULL;
133     PetscPartitioner part;
134 
135     ierr = DMPlexGetPartitioner(dmOrig, &part); CHKERRQ(ierr);
136     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
137     ierr = DMPlexDistribute(dmOrig, 0, NULL, &dmDist); CHKERRQ(ierr);
138     if (dmDist) {
139       ierr = DMDestroy(&dmOrig); CHKERRQ(ierr);
140       dmOrig = dmDist;
141     }
142   }
143 
144   // Allocate arrays for PETSc objects for each level
145   switch (coarsen) {
146   case COARSEN_UNIFORM:
147     numlevels = degree;
148     break;
149   case COARSEN_LOGARITHMIC:
150     numlevels = ceil(log(degree)/log(2)) + 1;
151     break;
152   }
153   ierr = PetscMalloc1(numlevels, &leveldegrees); CHKERRQ(ierr);
154   switch (coarsen) {
155   case COARSEN_UNIFORM:
156     for (int i=0; i<numlevels; i++) leveldegrees[i] = i + 1;
157     break;
158   case COARSEN_LOGARITHMIC:
159     for (int i=0; i<numlevels-1; i++) leveldegrees[i] = pow(2,i);
160     leveldegrees[numlevels-1] = degree;
161     break;
162   }
163   ierr = PetscMalloc1(numlevels, &dm); CHKERRQ(ierr);
164   ierr = PetscMalloc1(numlevels, &X); CHKERRQ(ierr);
165   ierr = PetscMalloc1(numlevels, &Xloc); CHKERRQ(ierr);
166   ierr = PetscMalloc1(numlevels, &mult); CHKERRQ(ierr);
167   ierr = PetscMalloc1(numlevels, &userO); CHKERRQ(ierr);
168   ierr = PetscMalloc1(numlevels, &userI); CHKERRQ(ierr);
169   ierr = PetscMalloc1(numlevels, &userR); CHKERRQ(ierr);
170   ierr = PetscMalloc1(numlevels, &matO); CHKERRQ(ierr);
171   ierr = PetscMalloc1(numlevels, &matI); CHKERRQ(ierr);
172   ierr = PetscMalloc1(numlevels, &matR); CHKERRQ(ierr);
173   ierr = PetscMalloc1(numlevels, &lsize); CHKERRQ(ierr);
174   ierr = PetscMalloc1(numlevels, &xlsize); CHKERRQ(ierr);
175   ierr = PetscMalloc1(numlevels, &gsize); CHKERRQ(ierr);
176 
177   // Setup DM and Operator Mat Shells for each level
178   for (CeedInt i=0; i<numlevels; i++) {
179     // Create DM
180     ierr = DMClone(dmOrig, &dm[i]); CHKERRQ(ierr);
181     ierr = SetupDMByDegree(dm[i], leveldegrees[i], ncompu, bpChoice);
182     CHKERRQ(ierr);
183 
184     // Create vectors
185     ierr = DMCreateGlobalVector(dm[i], &X[i]); CHKERRQ(ierr);
186     ierr = VecGetLocalSize(X[i], &lsize[i]); CHKERRQ(ierr);
187     ierr = VecGetSize(X[i], &gsize[i]); CHKERRQ(ierr);
188     ierr = DMCreateLocalVector(dm[i], &Xloc[i]); CHKERRQ(ierr);
189     ierr = VecGetSize(Xloc[i], &xlsize[i]); CHKERRQ(ierr);
190 
191     // Operator
192     ierr = PetscMalloc1(1, &userO[i]); CHKERRQ(ierr);
193     ierr = MatCreateShell(comm, lsize[i], lsize[i], gsize[i], gsize[i],
194                           userO[i], &matO[i]); CHKERRQ(ierr);
195     ierr = MatShellSetOperation(matO[i], MATOP_MULT,
196                                 (void(*)(void))MatMult_Ceed);
197     ierr = MatShellSetOperation(matO[i], MATOP_GET_DIAGONAL,
198                                 (void(*)(void))MatGetDiag);
199     CHKERRQ(ierr);
200 
201     // Level transfers
202     if (i > 0) {
203       // Interp
204       ierr = PetscMalloc1(1, &userI[i]); CHKERRQ(ierr);
205       ierr = MatCreateShell(comm, lsize[i], lsize[i-1], gsize[i], gsize[i-1],
206                             userI[i], &matI[i]); CHKERRQ(ierr);
207       ierr = MatShellSetOperation(matI[i], MATOP_MULT,
208                                   (void(*)(void))MatMult_Interp);
209       CHKERRQ(ierr);
210 
211       // Restrict
212       ierr = PetscMalloc1(1, &userR[i]); CHKERRQ(ierr);
213       ierr = MatCreateShell(comm, lsize[i-1], lsize[i], gsize[i-1], gsize[i],
214                             userR[i], &matR[i]); CHKERRQ(ierr);
215       ierr = MatShellSetOperation(matR[i], MATOP_MULT,
216                                   (void(*)(void))MatMult_Restrict);
217       CHKERRQ(ierr);
218     }
219   }
220   ierr = VecDuplicate(X[numlevels-1], &rhs); CHKERRQ(ierr);
221 
222   // Set up libCEED
223   CeedInit(ceedresource, &ceed);
224 
225   // Print global grid information
226   if (!test_mode) {
227     PetscInt P = degree + 1, Q = P + qextra;
228     const char *usedresource;
229     CeedGetResource(ceed, &usedresource);
230     ierr = PetscPrintf(comm,
231                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc + PCMG --\n"
232                        "  libCEED:\n"
233                        "    libCEED Backend                    : %s\n"
234                        "  Mesh:\n"
235                        "    Number of 1D Basis Nodes (p)       : %d\n"
236                        "    Number of 1D Quadrature Points (q) : %d\n"
237                        "    Global Nodes                       : %D\n"
238                        "    Owned Nodes                        : %D\n"
239                        "  Multigrid:\n"
240                        "    Number of Levels                   : %d\n",
241                        bpChoice+1, usedresource, P, Q,
242                        gsize[numlevels-1]/ncompu, lsize[numlevels-1]/ncompu,
243                        numlevels); CHKERRQ(ierr);
244   }
245 
246   // Create RHS vector
247   ierr = VecDuplicate(Xloc[numlevels-1], &rhsloc); CHKERRQ(ierr);
248   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
249   ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr);
250   CeedVectorCreate(ceed, xlsize[numlevels-1], &rhsceed);
251   CeedVectorSetArray(rhsceed, CEED_MEM_HOST, CEED_USE_POINTER, r);
252 
253   // Set up libCEED operators on each level
254   ierr = PetscMalloc1(numlevels, &ceeddata); CHKERRQ(ierr);
255   for (int i=0; i<numlevels; i++) {
256     // Print level information
257     if (!test_mode && (i == 0 || i == numlevels-1)) {
258       ierr = PetscPrintf(comm,"    Level %D (%s):\n"
259                          "      Number of 1D Basis Nodes (p)     : %d\n"
260                          "      Global Nodes                     : %D\n"
261                          "      Owned Nodes                      : %D\n",
262                          i, (i? "fine" : "coarse"), leveldegrees[i] + 1,
263                          gsize[i]/ncompu, lsize[i]/ncompu); CHKERRQ(ierr);
264     }
265     ierr = PetscMalloc1(1, &ceeddata[i]); CHKERRQ(ierr);
266     ierr = SetupLibceedByDegree(dm[i], ceed, leveldegrees[i], dim, qextra,
267                                 ncompu, gsize[i], xlsize[i], bpChoice,
268                                 ceeddata[i], i==(numlevels-1), rhsceed,
269                                 &target); CHKERRQ(ierr);
270   }
271 
272   // Gather RHS
273   ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr);
274   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
275   ierr = DMLocalToGlobalBegin(dm[numlevels-1], rhsloc, ADD_VALUES, rhs);
276   CHKERRQ(ierr);
277   ierr = DMLocalToGlobalEnd(dm[numlevels-1], rhsloc, ADD_VALUES, rhs);
278   CHKERRQ(ierr);
279   CeedVectorDestroy(&rhsceed);
280 
281   // Create the restriction/interpolation Q-function
282   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].ident,
283                               bpOptions[bpChoice].identfname, &qf_restrict);
284   CeedQFunctionAddInput(qf_restrict, "uin", ncompu, CEED_EVAL_NONE);
285   CeedQFunctionAddOutput(qf_restrict, "uout", ncompu, CEED_EVAL_INTERP);
286   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].ident,
287                               bpOptions[bpChoice].identfname, &qf_prolong);
288   CeedQFunctionAddInput(qf_prolong, "uin", ncompu, CEED_EVAL_INTERP);
289   CeedQFunctionAddOutput(qf_prolong, "uout", ncompu, CEED_EVAL_NONE);
290 
291   // Set up libCEED level transfer operators
292   ierr = CeedLevelTransferSetup(ceed, numlevels, ncompu, bpChoice, ceeddata,
293                                 leveldegrees, qf_restrict, qf_prolong);
294   CHKERRQ(ierr);
295 
296   // Create the error Q-function
297   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].error,
298                               bpOptions[bpChoice].errorfname, &qf_error);
299   CeedQFunctionAddInput(qf_error, "u", ncompu, CEED_EVAL_INTERP);
300   CeedQFunctionAddInput(qf_error, "true_soln", ncompu, CEED_EVAL_NONE);
301   CeedQFunctionAddOutput(qf_error, "error", ncompu, CEED_EVAL_NONE);
302 
303   // Create the error operator
304   CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error);
305   CeedOperatorSetField(op_error, "u", ceeddata[numlevels-1]->Erestrictu,
306                        CEED_TRANSPOSE, ceeddata[numlevels-1]->basisu,
307                        CEED_VECTOR_ACTIVE);
308   CeedOperatorSetField(op_error, "true_soln", ceeddata[numlevels-1]->Erestrictui,
309                        CEED_NOTRANSPOSE, CEED_BASIS_COLLOCATED, target);
310   CeedOperatorSetField(op_error, "error", ceeddata[numlevels-1]->Erestrictui,
311                        CEED_NOTRANSPOSE, CEED_BASIS_COLLOCATED,
312                        CEED_VECTOR_ACTIVE);
313 
314   // Calculate multiplicity
315   for (int i=0; i<numlevels; i++) {
316     PetscScalar *x;
317 
318     // CEED vector
319     ierr = VecGetArray(Xloc[i], &x); CHKERRQ(ierr);
320     CeedVectorSetArray(ceeddata[i]->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
321 
322     // Multiplicity
323     CeedElemRestrictionGetMultiplicity(ceeddata[i]->Erestrictu,
324                                        ceeddata[i]->xceed);
325 
326     // Restore vector
327     ierr = VecRestoreArray(Xloc[i], &x); CHKERRQ(ierr);
328 
329     // Creat mult vector
330     ierr = VecDuplicate(Xloc[i], &mult[i]); CHKERRQ(ierr);
331 
332     // Local-to-global
333     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
334     ierr = DMLocalToGlobalBegin(dm[i], Xloc[i], ADD_VALUES, X[i]);
335     CHKERRQ(ierr);
336     ierr = DMLocalToGlobalEnd(dm[i], Xloc[i], ADD_VALUES, X[i]);
337     CHKERRQ(ierr);
338     ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr);
339 
340     // Global-to-local
341     ierr = DMGlobalToLocalBegin(dm[i], X[i], INSERT_VALUES, mult[i]);
342     CHKERRQ(ierr);
343     ierr = DMGlobalToLocalEnd(dm[i], X[i], INSERT_VALUES, mult[i]);
344     CHKERRQ(ierr);
345     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
346 
347     // Multiplicity scaling
348     ierr = VecReciprocal(mult[i]);
349   }
350 
351   // Set up Mat
352   for (int i=0; i<numlevels; i++) {
353     // User Operator
354     userO[i]->comm = comm;
355     userO[i]->dm = dm[i];
356     userO[i]->Xloc = Xloc[i];
357     ierr = VecDuplicate(Xloc[i], &userO[i]->Yloc); CHKERRQ(ierr);
358     userO[i]->xceed = ceeddata[i]->xceed;
359     userO[i]->yceed = ceeddata[i]->yceed;
360     userO[i]->op = ceeddata[i]->op_apply;
361     userO[i]->ceed = ceed;
362 
363     // Set up diagonal
364     const CeedScalar *ceedarray;
365     PetscScalar *petscarray;
366     CeedInt length;
367 
368     ierr = VecDuplicate(X[i], &userO[i]->diag); CHKERRQ(ierr);
369     ierr = VecDuplicate(Xloc[i], &diagloc); CHKERRQ(ierr);
370 
371     // -- Local diagonal
372     CeedOperatorAssembleLinearDiagonal(userO[i]->op, &diagceed,
373                                        CEED_REQUEST_IMMEDIATE);
374 
375     // -- Copy values
376     CeedVectorGetArrayRead(diagceed, CEED_MEM_HOST, &ceedarray);
377     ierr = VecGetArray(diagloc, &petscarray); CHKERRQ(ierr);
378     CeedVectorGetLength(diagceed, &length);
379     for (CeedInt i=0; i<length; i++)
380       petscarray[i] = ceedarray[i];
381     CeedVectorRestoreArrayRead(diagceed, &ceedarray);
382     ierr = VecRestoreArray(diagloc, &petscarray); CHKERRQ(ierr);
383 
384     // -- Global diagonal
385     ierr = VecZeroEntries(userO[i]->diag); CHKERRQ(ierr);
386     ierr = DMLocalToGlobalBegin(userO[i]->dm, diagloc, ADD_VALUES,
387                                 userO[i]->diag); CHKERRQ(ierr);
388     ierr = DMLocalToGlobalEnd(userO[i]->dm, diagloc, ADD_VALUES,
389                               userO[i]->diag); CHKERRQ(ierr);
390 
391     // -- Cleanup
392     ierr = VecDestroy(&diagloc); CHKERRQ(ierr);
393     CeedVectorDestroy(&diagceed);
394 
395     if (i > 0) {
396       // Interp Operator
397       userI[i]->comm = comm;
398       userI[i]->dmc = dm[i-1];
399       userI[i]->dmf = dm[i];
400       userI[i]->Xloc = Xloc[i-1];
401       userI[i]->Yloc = userO[i]->Yloc;
402       userI[i]->mult = mult[i];
403       userI[i]->ceedvecc = userO[i-1]->xceed;
404       userI[i]->ceedvecf = userO[i]->yceed;
405       userI[i]->op = ceeddata[i]->op_interp;
406       userI[i]->ceed = ceed;
407 
408       // Restrict Operator
409       userR[i]->comm = comm;
410       userR[i]->dmc = dm[i-1];
411       userR[i]->dmf = dm[i];
412       userR[i]->Xloc = Xloc[i];
413       userR[i]->Yloc = userO[i-1]->Yloc;
414       userR[i]->mult = mult[i];
415       userR[i]->ceedvecf = userO[i]->xceed;
416       userR[i]->ceedvecc = userO[i-1]->yceed;
417       userR[i]->op = ceeddata[i]->op_restrict;
418       userR[i]->ceed = ceed;
419     }
420   }
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, KSPCG); CHKERRQ(ierr);
474       ierr = KSPSetOperators(coarse, matO[0], matO[0]); CHKERRQ(ierr);
475       ierr = KSPSetTolerances(coarse, 1e-10, 1e-10, PETSC_DEFAULT,
476                               PETSC_DEFAULT); CHKERRQ(ierr);
477       ierr = KSPGetPC(coarse, &coarse_pc); CHKERRQ(ierr);
478       ierr = PCSetType(coarse_pc, PCJACOBI); CHKERRQ(ierr);
479       ierr = PCJacobiSetType(coarse_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
480     }
481 
482     // PCMG options
483     ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr);
484     ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr);
485     ierr = PCMGSetCycleType(pc, pcgmcycletype); CHKERRQ(ierr);
486   }
487 
488   // First run, if benchmarking
489   if (benchmark_mode) {
490     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
491     CHKERRQ(ierr);
492     ierr = VecZeroEntries(X[numlevels-1]); CHKERRQ(ierr);
493     my_rt_start = MPI_Wtime();
494     ierr = KSPSolve(ksp, rhs, X[numlevels-1]); CHKERRQ(ierr);
495     my_rt = MPI_Wtime() - my_rt_start;
496     ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
497     CHKERRQ(ierr);
498     // Set maxits based on first iteration timing
499     if (my_rt > 0.02) {
500       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
501       CHKERRQ(ierr);
502     } else {
503       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
504       CHKERRQ(ierr);
505     }
506   }
507 
508   // Timed solve
509   ierr = VecZeroEntries(X[numlevels-1]); CHKERRQ(ierr);
510   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
511   my_rt_start = MPI_Wtime();
512   ierr = KSPSolve(ksp, rhs, X[numlevels-1]); CHKERRQ(ierr);
513   my_rt = MPI_Wtime() - my_rt_start;
514 
515   // Output results
516   {
517     KSPType ksptype;
518     PCMGType pcmgtype;
519     KSPConvergedReason reason;
520     PetscReal rnorm;
521     PetscInt its;
522     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
523     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
524     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
525     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
526     ierr = PCMGGetType(pc, &pcmgtype); CHKERRQ(ierr);
527     if (!test_mode || reason < 0 || rnorm > 1e-8) {
528       ierr = PetscPrintf(comm,
529                          "  KSP:\n"
530                          "    KSP Type                           : %s\n"
531                          "    KSP Convergence                    : %s\n"
532                          "    Total KSP Iterations               : %D\n"
533                          "    Final rnorm                        : %e\n",
534                          ksptype, KSPConvergedReasons[reason], its,
535                          (double)rnorm); CHKERRQ(ierr);
536       ierr = PetscPrintf(comm,
537                          "  PCMG:\n"
538                          "    PCMG Type                          : %s\n"
539                          "    PCMG Cycle Type                    : %s\n",
540                          PCMGTypes[pcmgtype],
541                          PCMGCycleTypes[pcgmcycletype]); CHKERRQ(ierr);
542     }
543     if (!test_mode) {
544       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
545     }
546     {
547       PetscReal maxerror;
548       ierr = ComputeErrorMax(userO[numlevels-1], op_error, X[numlevels-1], target,
549                              &maxerror); CHKERRQ(ierr);
550       PetscReal tol = 5e-2;
551       if (!test_mode || maxerror > tol) {
552         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
553         CHKERRQ(ierr);
554         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
555         CHKERRQ(ierr);
556         ierr = PetscPrintf(comm,
557                            "    Pointwise Error (max)              : %e\n"
558                            "    CG Solve Time                      : %g (%g) sec\n",
559                            (double)maxerror, rt_max, rt_min); CHKERRQ(ierr);
560       }
561     }
562     if (benchmark_mode && (!test_mode)) {
563       ierr = PetscPrintf(comm,
564                          "    DoFs/Sec in CG                     : %g (%g) million\n",
565                          1e-6*gsize[numlevels-1]*its/rt_max,
566                          1e-6*gsize[numlevels-1]*its/rt_min);
567       CHKERRQ(ierr);
568     }
569   }
570 
571   if (write_solution) {
572     PetscViewer vtkviewersoln;
573 
574     ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr);
575     ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr);
576     ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr);
577     ierr = VecView(X[numlevels-1], vtkviewersoln); CHKERRQ(ierr);
578     ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr);
579   }
580 
581   // Cleanup
582   for (int i=0; i<numlevels; i++) {
583     ierr = VecDestroy(&X[i]); CHKERRQ(ierr);
584     ierr = VecDestroy(&Xloc[i]); CHKERRQ(ierr);
585     ierr = VecDestroy(&mult[i]); CHKERRQ(ierr);
586     ierr = VecDestroy(&userO[i]->Yloc); CHKERRQ(ierr);
587     ierr = VecDestroy(&userO[i]->diag); CHKERRQ(ierr);
588     ierr = MatDestroy(&matO[i]); CHKERRQ(ierr);
589     ierr = PetscFree(userO[i]); CHKERRQ(ierr);
590     if (i > 0) {
591       ierr = MatDestroy(&matI[i]); CHKERRQ(ierr);
592       ierr = PetscFree(userI[i]); CHKERRQ(ierr);
593       ierr = MatDestroy(&matR[i]); CHKERRQ(ierr);
594       ierr = PetscFree(userR[i]); CHKERRQ(ierr);
595     }
596     ierr = CeedDataDestroy(i, ceeddata[i]); CHKERRQ(ierr);
597     ierr = DMDestroy(&dm[i]); CHKERRQ(ierr);
598   }
599   ierr = PetscFree(leveldegrees); CHKERRQ(ierr);
600   ierr = PetscFree(dm); CHKERRQ(ierr);
601   ierr = PetscFree(X); CHKERRQ(ierr);
602   ierr = PetscFree(Xloc); CHKERRQ(ierr);
603   ierr = PetscFree(mult); CHKERRQ(ierr);
604   ierr = PetscFree(matO); CHKERRQ(ierr);
605   ierr = PetscFree(matI); CHKERRQ(ierr);
606   ierr = PetscFree(matR); CHKERRQ(ierr);
607   ierr = PetscFree(ceeddata); CHKERRQ(ierr);
608   ierr = PetscFree(userO); CHKERRQ(ierr);
609   ierr = PetscFree(userI); CHKERRQ(ierr);
610   ierr = PetscFree(userR); CHKERRQ(ierr);
611   ierr = PetscFree(lsize); CHKERRQ(ierr);
612   ierr = PetscFree(xlsize); CHKERRQ(ierr);
613   ierr = PetscFree(gsize); CHKERRQ(ierr);
614   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
615   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
616   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
617   ierr = DMDestroy(&dmOrig); CHKERRQ(ierr);
618   CeedVectorDestroy(&target);
619   CeedQFunctionDestroy(&qf_error);
620   CeedQFunctionDestroy(&qf_restrict);
621   CeedQFunctionDestroy(&qf_prolong);
622   CeedOperatorDestroy(&op_error);
623   CeedDestroy(&ceed);
624   return PetscFinalize();
625 }
626