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