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