xref: /libCEED/examples/petsc/multigrid.c (revision ce74dcefbddf8db3528d9019b9e2fecff9c01924)
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, 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); CHKERRQ(ierr);
198     ierr = MatShellSetOperation(matO[i], MATOP_GET_DIAGONAL,
199                                 (void(*)(void))MatGetDiag); 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                        "    DoF per node                       : %D\n"
240                        "  Multigrid:\n"
241                        "    Number of Levels                   : %d\n",
242                        bpChoice+1, usedresource, P, Q,
243                        gsize[numlevels-1]/ncompu, lsize[numlevels-1]/ncompu,
244                        ncompu, numlevels); CHKERRQ(ierr);
245   }
246 
247   // Create RHS vector
248   ierr = VecDuplicate(Xloc[numlevels-1], &rhsloc); CHKERRQ(ierr);
249   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
250   ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr);
251   CeedVectorCreate(ceed, xlsize[numlevels-1], &rhsceed);
252   CeedVectorSetArray(rhsceed, CEED_MEM_HOST, CEED_USE_POINTER, r);
253 
254   // Set up libCEED operators on each level
255   ierr = PetscMalloc1(numlevels, &ceeddata); CHKERRQ(ierr);
256   for (int i=0; i<numlevels; i++) {
257     // Print level information
258     if (!test_mode && (i == 0 || i == numlevels-1)) {
259       ierr = PetscPrintf(comm,"    Level %D (%s):\n"
260                          "      Number of 1D Basis Nodes (p)     : %d\n"
261                          "      Global Nodes                     : %D\n"
262                          "      Owned Nodes                      : %D\n",
263                          i, (i? "fine" : "coarse"), leveldegrees[i] + 1,
264                          gsize[i]/ncompu, lsize[i]/ncompu); CHKERRQ(ierr);
265     }
266     ierr = PetscMalloc1(1, &ceeddata[i]); CHKERRQ(ierr);
267     ierr = SetupLibceedByDegree(dm[i], ceed, leveldegrees[i], dim, qextra,
268                                 ncompu, gsize[i], xlsize[i], bpChoice,
269                                 ceeddata[i], i==(numlevels-1), rhsceed,
270                                 &target); CHKERRQ(ierr);
271   }
272 
273   // Gather RHS
274   ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr);
275   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
276   ierr = DMLocalToGlobal(dm[numlevels-1], rhsloc, ADD_VALUES, rhs);
277   CHKERRQ(ierr);
278   CeedVectorDestroy(&rhsceed);
279 
280   // Create the restriction/interpolation Q-function
281   CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_NONE, CEED_EVAL_INTERP,
282                               &qf_restrict);
283   CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_INTERP, CEED_EVAL_NONE,
284                               &qf_prolong);
285 
286   // Set up libCEED level transfer operators
287   ierr = CeedLevelTransferSetup(ceed, numlevels, ncompu, bpChoice, ceeddata,
288                                 leveldegrees, qf_restrict, qf_prolong);
289   CHKERRQ(ierr);
290 
291   // Create the error Q-function
292   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].error,
293                               bpOptions[bpChoice].errorfname, &qf_error);
294   CeedQFunctionAddInput(qf_error, "u", ncompu, CEED_EVAL_INTERP);
295   CeedQFunctionAddInput(qf_error, "true_soln", ncompu, CEED_EVAL_NONE);
296   CeedQFunctionAddOutput(qf_error, "error", ncompu, CEED_EVAL_NONE);
297 
298   // Create the error operator
299   CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
300                      &op_error);
301   CeedOperatorSetField(op_error, "u", ceeddata[numlevels-1]->Erestrictu,
302                        ceeddata[numlevels-1]->basisu, CEED_VECTOR_ACTIVE);
303   CeedOperatorSetField(op_error, "true_soln", ceeddata[numlevels-1]->Erestrictui,
304                        CEED_BASIS_COLLOCATED, target);
305   CeedOperatorSetField(op_error, "error", ceeddata[numlevels-1]->Erestrictui,
306                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
307 
308   // Calculate multiplicity
309   for (int i=0; i<numlevels; i++) {
310     PetscScalar *x;
311 
312     // CEED vector
313     ierr = VecGetArray(Xloc[i], &x); CHKERRQ(ierr);
314     CeedVectorSetArray(ceeddata[i]->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
315 
316     // Multiplicity
317     CeedElemRestrictionGetMultiplicity(ceeddata[i]->Erestrictu,
318                                        ceeddata[i]->xceed);
319 
320     // Restore vector
321     ierr = VecRestoreArray(Xloc[i], &x); CHKERRQ(ierr);
322 
323     // Creat mult vector
324     ierr = VecDuplicate(Xloc[i], &mult[i]); CHKERRQ(ierr);
325 
326     // Local-to-global
327     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
328     ierr = DMLocalToGlobal(dm[i], Xloc[i], ADD_VALUES, X[i]);
329     CHKERRQ(ierr);
330     ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr);
331 
332     // Global-to-local
333     ierr = DMGlobalToLocal(dm[i], X[i], INSERT_VALUES, mult[i]);
334     CHKERRQ(ierr);
335     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
336 
337     // Multiplicity scaling
338     ierr = VecReciprocal(mult[i]);
339   }
340 
341   // Set up Mat
342   for (int i=0; i<numlevels; i++) {
343     // User Operator
344     userO[i]->comm = comm;
345     userO[i]->dm = dm[i];
346     userO[i]->Xloc = Xloc[i];
347     ierr = VecDuplicate(Xloc[i], &userO[i]->Yloc); CHKERRQ(ierr);
348     userO[i]->xceed = ceeddata[i]->xceed;
349     userO[i]->yceed = ceeddata[i]->yceed;
350     userO[i]->op = ceeddata[i]->op_apply;
351     userO[i]->ceed = ceed;
352 
353     if (i > 0) {
354       // Interp Operator
355       userI[i]->comm = comm;
356       userI[i]->dmc = dm[i-1];
357       userI[i]->dmf = dm[i];
358       userI[i]->Xloc = Xloc[i-1];
359       userI[i]->Yloc = userO[i]->Yloc;
360       userI[i]->mult = mult[i];
361       userI[i]->ceedvecc = userO[i-1]->xceed;
362       userI[i]->ceedvecf = userO[i]->yceed;
363       userI[i]->op = ceeddata[i]->op_interp;
364       userI[i]->ceed = ceed;
365 
366       // Restrict Operator
367       userR[i]->comm = comm;
368       userR[i]->dmc = dm[i-1];
369       userR[i]->dmf = dm[i];
370       userR[i]->Xloc = Xloc[i];
371       userR[i]->Yloc = userO[i-1]->Yloc;
372       userR[i]->mult = mult[i];
373       userR[i]->ceedvecf = userO[i]->xceed;
374       userR[i]->ceedvecc = userO[i-1]->yceed;
375       userR[i]->op = ceeddata[i]->op_restrict;
376       userR[i]->ceed = ceed;
377     }
378   }
379 
380   // Setup dummy SNES for AMG coarse solve
381   ierr = SNESCreate(comm, &snes_dummy); CHKERRQ(ierr);
382   ierr = SNESSetDM(snes_dummy, dm[0]); CHKERRQ(ierr);
383   ierr = SNESSetSolution(snes_dummy, X[0]); CHKERRQ(ierr);
384 
385   // -- Jacobian matrix
386   ierr = DMSetMatType(dm[0], MATAIJ); CHKERRQ(ierr);
387   ierr = DMCreateMatrix(dm[0], &matCoarse); CHKERRQ(ierr);
388   ierr = SNESSetJacobian(snes_dummy, matCoarse, matCoarse, NULL,
389                          NULL); CHKERRQ(ierr);
390 
391   // -- Residual evaluation function
392   ierr = SNESSetFunction(snes_dummy, X[0], FormResidual_Ceed,
393                          userO[0]); CHKERRQ(ierr);
394 
395   // -- Form Jacobian
396   ierr = SNESComputeJacobianDefaultColor(snes_dummy, X[0], matO[0],
397                                          matCoarse, NULL); CHKERRQ(ierr);
398 
399   // Set up KSP
400   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
401   {
402     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
403     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
404     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
405                             PETSC_DEFAULT); CHKERRQ(ierr);
406   }
407   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
408   ierr = KSPSetOperators(ksp, matO[numlevels-1], matO[numlevels-1]);
409   CHKERRQ(ierr);
410 
411   // Set up PCMG
412   ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
413   PCMGCycleType pcgmcycletype = PC_MG_CYCLE_V;
414   {
415     ierr = PCSetType(pc, PCMG); CHKERRQ(ierr);
416 
417     // PCMG levels
418     ierr = PCMGSetLevels(pc, numlevels, NULL); CHKERRQ(ierr);
419     for (int i=0; i<numlevels; i++) {
420       // Smoother
421       KSP smoother;
422       PC smoother_pc;
423       ierr = PCMGGetSmoother(pc, i, &smoother); CHKERRQ(ierr);
424       ierr = KSPSetType(smoother, KSPCHEBYSHEV); CHKERRQ(ierr);
425       ierr = KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1); CHKERRQ(ierr);
426       ierr = KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE); CHKERRQ(ierr);
427       ierr = KSPSetOperators(smoother, matO[i], matO[i]); CHKERRQ(ierr);
428       ierr = KSPGetPC(smoother, &smoother_pc); CHKERRQ(ierr);
429       ierr = PCSetType(smoother_pc, PCJACOBI); CHKERRQ(ierr);
430       ierr = PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
431 
432       // Work vector
433       if (i < numlevels-1) {
434         ierr = PCMGSetX(pc, i, X[i]); CHKERRQ(ierr);
435       }
436 
437       // Level transfers
438       if (i > 0) {
439         // Interpolation
440         ierr = PCMGSetInterpolation(pc, i, matI[i]); CHKERRQ(ierr);
441 
442         // Restriction
443         ierr = PCMGSetRestriction(pc, i, matR[i]); CHKERRQ(ierr);
444       }
445 
446       // Coarse solve
447       KSP coarse;
448       PC coarse_pc;
449       ierr = PCMGGetCoarseSolve(pc, &coarse); CHKERRQ(ierr);
450       ierr = KSPSetType(coarse, KSPPREONLY); CHKERRQ(ierr);
451       ierr = KSPSetOperators(coarse, matCoarse, matCoarse); CHKERRQ(ierr);
452 
453       ierr = KSPGetPC(coarse, &coarse_pc); CHKERRQ(ierr);
454       ierr = PCSetType(coarse_pc, PCGAMG); CHKERRQ(ierr);
455 
456       ierr = KSPSetOptionsPrefix(coarse, "coarse_"); CHKERRQ(ierr);
457       ierr = PCSetOptionsPrefix(coarse_pc, "coarse_"); CHKERRQ(ierr);
458       ierr = KSPSetFromOptions(coarse); CHKERRQ(ierr);
459       ierr = PCSetFromOptions(coarse_pc); CHKERRQ(ierr);
460     }
461 
462     // PCMG options
463     ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr);
464     ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr);
465     ierr = PCMGSetCycleType(pc, pcgmcycletype); CHKERRQ(ierr);
466   }
467 
468   // First run, if benchmarking
469   if (benchmark_mode) {
470     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
471     CHKERRQ(ierr);
472     ierr = VecZeroEntries(X[numlevels-1]); CHKERRQ(ierr);
473     my_rt_start = MPI_Wtime();
474     ierr = KSPSolve(ksp, rhs, X[numlevels-1]); CHKERRQ(ierr);
475     my_rt = MPI_Wtime() - my_rt_start;
476     ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
477     CHKERRQ(ierr);
478     // Set maxits based on first iteration timing
479     if (my_rt > 0.02) {
480       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
481       CHKERRQ(ierr);
482     } else {
483       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
484       CHKERRQ(ierr);
485     }
486   }
487 
488   // Timed solve
489   ierr = VecZeroEntries(X[numlevels-1]); CHKERRQ(ierr);
490   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
491   my_rt_start = MPI_Wtime();
492   ierr = KSPSolve(ksp, rhs, X[numlevels-1]); CHKERRQ(ierr);
493   my_rt = MPI_Wtime() - my_rt_start;
494 
495   // Output results
496   {
497     KSPType ksptype;
498     PCMGType pcmgtype;
499     KSPConvergedReason reason;
500     PetscReal rnorm;
501     PetscInt its;
502     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
503     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
504     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
505     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
506     ierr = PCMGGetType(pc, &pcmgtype); CHKERRQ(ierr);
507     if (!test_mode || reason < 0 || rnorm > 1e-8) {
508       ierr = PetscPrintf(comm,
509                          "  KSP:\n"
510                          "    KSP Type                           : %s\n"
511                          "    KSP Convergence                    : %s\n"
512                          "    Total KSP Iterations               : %D\n"
513                          "    Final rnorm                        : %e\n",
514                          ksptype, KSPConvergedReasons[reason], its,
515                          (double)rnorm); CHKERRQ(ierr);
516       ierr = PetscPrintf(comm,
517                          "  PCMG:\n"
518                          "    PCMG Type                          : %s\n"
519                          "    PCMG Cycle Type                    : %s\n",
520                          PCMGTypes[pcmgtype],
521                          PCMGCycleTypes[pcgmcycletype]); CHKERRQ(ierr);
522     }
523     if (!test_mode) {
524       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
525     }
526     {
527       PetscReal maxerror;
528       ierr = ComputeErrorMax(userO[numlevels-1], op_error, X[numlevels-1], target,
529                              &maxerror); CHKERRQ(ierr);
530       PetscReal tol = 5e-2;
531       if (!test_mode || maxerror > tol) {
532         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
533         CHKERRQ(ierr);
534         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
535         CHKERRQ(ierr);
536         ierr = PetscPrintf(comm,
537                            "    Pointwise Error (max)              : %e\n"
538                            "    CG Solve Time                      : %g (%g) sec\n",
539                            (double)maxerror, rt_max, rt_min); CHKERRQ(ierr);
540       }
541     }
542     if (benchmark_mode && (!test_mode)) {
543       ierr = PetscPrintf(comm,
544                          "    DoFs/Sec in CG                     : %g (%g) million\n",
545                          1e-6*gsize[numlevels-1]*its/rt_max,
546                          1e-6*gsize[numlevels-1]*its/rt_min);
547       CHKERRQ(ierr);
548     }
549   }
550 
551   if (write_solution) {
552     PetscViewer vtkviewersoln;
553 
554     ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr);
555     ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr);
556     ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr);
557     ierr = VecView(X[numlevels-1], vtkviewersoln); CHKERRQ(ierr);
558     ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr);
559   }
560 
561   // Cleanup
562   for (int i=0; i<numlevels; i++) {
563     ierr = VecDestroy(&X[i]); CHKERRQ(ierr);
564     ierr = VecDestroy(&Xloc[i]); CHKERRQ(ierr);
565     ierr = VecDestroy(&mult[i]); CHKERRQ(ierr);
566     ierr = VecDestroy(&userO[i]->Yloc); CHKERRQ(ierr);
567     ierr = MatDestroy(&matO[i]); CHKERRQ(ierr);
568     ierr = PetscFree(userO[i]); CHKERRQ(ierr);
569     if (i > 0) {
570       ierr = MatDestroy(&matI[i]); CHKERRQ(ierr);
571       ierr = PetscFree(userI[i]); CHKERRQ(ierr);
572       ierr = MatDestroy(&matR[i]); CHKERRQ(ierr);
573       ierr = PetscFree(userR[i]); CHKERRQ(ierr);
574     }
575     ierr = CeedDataDestroy(i, ceeddata[i]); CHKERRQ(ierr);
576     ierr = DMDestroy(&dm[i]); CHKERRQ(ierr);
577   }
578   ierr = PetscFree(leveldegrees); CHKERRQ(ierr);
579   ierr = PetscFree(dm); CHKERRQ(ierr);
580   ierr = PetscFree(X); CHKERRQ(ierr);
581   ierr = PetscFree(Xloc); CHKERRQ(ierr);
582   ierr = PetscFree(mult); CHKERRQ(ierr);
583   ierr = PetscFree(matO); CHKERRQ(ierr);
584   ierr = PetscFree(matI); CHKERRQ(ierr);
585   ierr = PetscFree(matR); CHKERRQ(ierr);
586   ierr = PetscFree(ceeddata); CHKERRQ(ierr);
587   ierr = PetscFree(userO); CHKERRQ(ierr);
588   ierr = PetscFree(userI); CHKERRQ(ierr);
589   ierr = PetscFree(userR); CHKERRQ(ierr);
590   ierr = PetscFree(lsize); CHKERRQ(ierr);
591   ierr = PetscFree(xlsize); CHKERRQ(ierr);
592   ierr = PetscFree(gsize); CHKERRQ(ierr);
593   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
594   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
595   ierr = MatDestroy(&matCoarse); CHKERRQ(ierr);
596   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
597   ierr = SNESDestroy(&snes_dummy); CHKERRQ(ierr);
598   ierr = DMDestroy(&dmOrig); CHKERRQ(ierr);
599   CeedVectorDestroy(&target);
600   CeedQFunctionDestroy(&qf_error);
601   CeedQFunctionDestroy(&qf_restrict);
602   CeedQFunctionDestroy(&qf_prolong);
603   CeedOperatorDestroy(&op_error);
604   CeedDestroy(&ceed);
605   return PetscFinalize();
606 }
607