xref: /libCEED/examples/petsc/multigrid.c (revision 2d03409ca85173e0565dadbac5c9383b99bac8d1)
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 ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self",
48        filename[PETSC_MAX_PATH_LEN];
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   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL,
112                             filename, filename, sizeof(filename), &read_mesh);
113   CHKERRQ(ierr);
114   if (!read_mesh) {
115     PetscInt tmp = dim;
116     ierr = PetscOptionsIntArray("-cells","Number of cells per dimension", NULL,
117                                 melem, &tmp, NULL); CHKERRQ(ierr);
118   }
119   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
120 
121   // Setup DM
122   if (read_mesh) {
123     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dmOrig);
124     CHKERRQ(ierr);
125   } else {
126     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_FALSE, melem, NULL,
127                                NULL, NULL, PETSC_TRUE,&dmOrig); CHKERRQ(ierr);
128   }
129 
130   {
131     DM dmDist = NULL;
132     PetscPartitioner part;
133 
134     ierr = DMPlexGetPartitioner(dmOrig, &part); CHKERRQ(ierr);
135     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
136     ierr = DMPlexDistribute(dmOrig, 0, NULL, &dmDist); CHKERRQ(ierr);
137     if (dmDist) {
138       ierr = DMDestroy(&dmOrig); CHKERRQ(ierr);
139       dmOrig = dmDist;
140     }
141   }
142 
143   // Allocate arrays for PETSc objects for each level
144   switch (coarsen) {
145   case COARSEN_UNIFORM:
146     numlevels = degree;
147     break;
148   case COARSEN_LOGRITHMIC:
149     numlevels = ceil(log(degree)/log(2)) + 1;
150     break;
151   }
152   ierr = PetscMalloc1(numlevels, &leveldegrees); CHKERRQ(ierr);
153   switch (coarsen) {
154   case COARSEN_UNIFORM:
155     for (int i=0; i<numlevels; i++) leveldegrees[i] = i + 1;
156     break;
157   case COARSEN_LOGRITHMIC:
158     for (int i=0; i<numlevels-1; i++) leveldegrees[i] = pow(2,i);
159     leveldegrees[numlevels-1] = degree;
160     break;
161   }
162   ierr = PetscMalloc1(numlevels, &dm); CHKERRQ(ierr);
163   ierr = PetscMalloc1(numlevels, &X); CHKERRQ(ierr);
164   ierr = PetscMalloc1(numlevels, &Xloc); CHKERRQ(ierr);
165   ierr = PetscMalloc1(numlevels, &mult); CHKERRQ(ierr);
166   ierr = PetscMalloc1(numlevels, &userO); CHKERRQ(ierr);
167   ierr = PetscMalloc1(numlevels, &userI); CHKERRQ(ierr);
168   ierr = PetscMalloc1(numlevels, &userR); CHKERRQ(ierr);
169   ierr = PetscMalloc1(numlevels, &matO); CHKERRQ(ierr);
170   ierr = PetscMalloc1(numlevels, &matI); CHKERRQ(ierr);
171   ierr = PetscMalloc1(numlevels, &matR); 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);
196     ierr = MatShellSetOperation(matO[i], MATOP_GET_DIAGONAL,
197                                 (void(*)(void))MatGetDiag);
198     CHKERRQ(ierr);
199 
200     // Level transfers
201     if (i > 0) {
202       // Interp
203       ierr = PetscMalloc1(1, &userI[i]); CHKERRQ(ierr);
204       ierr = MatCreateShell(comm, lsize[i], lsize[i-1], gsize[i], gsize[i-1],
205                             userI[i], &matI[i]); CHKERRQ(ierr);
206       ierr = MatShellSetOperation(matI[i], MATOP_MULT,
207                                   (void(*)(void))MatMult_Interp);
208       CHKERRQ(ierr);
209 
210       // Restrict
211       ierr = PetscMalloc1(1, &userR[i]); CHKERRQ(ierr);
212       ierr = MatCreateShell(comm, lsize[i-1], lsize[i], gsize[i-1], gsize[i],
213                             userR[i], &matR[i]); CHKERRQ(ierr);
214       ierr = MatShellSetOperation(matR[i], MATOP_MULT,
215                                   (void(*)(void))MatMult_Restrict);
216       CHKERRQ(ierr);
217     }
218   }
219   ierr = VecDuplicate(X[numlevels-1], &rhs); CHKERRQ(ierr);
220 
221   // Set up libCEED
222   CeedInit(ceedresource, &ceed);
223 
224   // Print global grid information
225   if (!test_mode) {
226     PetscInt P = degree + 1, Q = P + qextra;
227     const char *usedresource;
228     CeedGetResource(ceed, &usedresource);
229     ierr = PetscPrintf(comm,
230                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc + PCMG --\n"
231                        "  libCEED:\n"
232                        "    libCEED Backend                    : %s\n"
233                        "  Mesh:\n"
234                        "    Number of 1D Basis Nodes (p)       : %d\n"
235                        "    Number of 1D Quadrature Points (q) : %d\n"
236                        "    Global Nodes                       : %D\n"
237                        "    Owned Nodes                        : %D\n"
238                        "  Multigrid:\n"
239                        "    Number of Levels                   : %d\n",
240                        bpChoice+1, usedresource, P, Q,
241                        gsize[numlevels-1]/ncompu, lsize[numlevels-1]/ncompu,
242                        numlevels); CHKERRQ(ierr);
243   }
244 
245   // Create RHS vector
246   ierr = VecDuplicate(Xloc[numlevels-1], &rhsloc); CHKERRQ(ierr);
247   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
248   ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr);
249   CeedVectorCreate(ceed, xlsize[numlevels-1], &rhsceed);
250   CeedVectorSetArray(rhsceed, CEED_MEM_HOST, CEED_USE_POINTER, r);
251 
252   // Set up libCEED operators on each level
253   ierr = PetscMalloc1(numlevels, &ceeddata); CHKERRQ(ierr);
254   for (int i=0; i<numlevels; i++) {
255     // Print level information
256     if (!test_mode && (i == 0 || i == numlevels-1)) {
257       ierr = PetscPrintf(comm,"    Level %D (%s):\n"
258                          "      Number of 1D Basis Nodes (p)     : %d\n"
259                          "      Global Nodes                     : %D\n"
260                          "      Owned Nodes                      : %D\n",
261                          i, (i? "fine" : "coarse"), leveldegrees[i] + 1,
262                          gsize[i]/ncompu, lsize[i]/ncompu); CHKERRQ(ierr);
263     }
264     ierr = PetscMalloc1(1, &ceeddata[i]); CHKERRQ(ierr);
265     ierr = SetupLibceedByDegree(dm[i], ceed, leveldegrees[i], dim, qextra,
266                                 ncompu, gsize[i], xlsize[i], bpChoice,
267                                 ceeddata[i], i==(numlevels-1), rhsceed,
268                                 &target); CHKERRQ(ierr);
269   }
270 
271   // Gather RHS
272   ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr);
273   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
274   ierr = DMLocalToGlobalBegin(dm[numlevels-1], rhsloc, ADD_VALUES, rhs);
275   CHKERRQ(ierr);
276   ierr = DMLocalToGlobalEnd(dm[numlevels-1], rhsloc, ADD_VALUES, rhs);
277   CHKERRQ(ierr);
278   CeedVectorDestroy(&rhsceed);
279 
280   // Create the restriction/interpolation Q-function
281   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].ident,
282                               bpOptions[bpChoice].identfname, &qf_restrict);
283   CeedQFunctionAddInput(qf_restrict, "uin", ncompu, CEED_EVAL_NONE);
284   CeedQFunctionAddOutput(qf_restrict, "uout", ncompu, CEED_EVAL_INTERP);
285   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].ident,
286                               bpOptions[bpChoice].identfname, &qf_prolong);
287   CeedQFunctionAddInput(qf_prolong, "uin", ncompu, CEED_EVAL_INTERP);
288   CeedQFunctionAddOutput(qf_prolong, "uout", ncompu, CEED_EVAL_NONE);
289 
290   // Set up libCEED level transfer operators
291   ierr = CeedLevelTransferSetup(ceed, numlevels, ncompu, bpChoice, ceeddata,
292                                 leveldegrees, qf_restrict, qf_prolong);
293   CHKERRQ(ierr);
294 
295   // Create the error Q-function
296   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].error,
297                               bpOptions[bpChoice].errorfname, &qf_error);
298   CeedQFunctionAddInput(qf_error, "u", ncompu, CEED_EVAL_INTERP);
299   CeedQFunctionAddInput(qf_error, "true_soln", ncompu, CEED_EVAL_NONE);
300   CeedQFunctionAddOutput(qf_error, "error", ncompu, CEED_EVAL_NONE);
301 
302   // Create the error operator
303   CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error);
304   CeedOperatorSetField(op_error, "u", ceeddata[numlevels-1]->Erestrictu,
305                        CEED_TRANSPOSE, ceeddata[numlevels-1]->basisu,
306                        CEED_VECTOR_ACTIVE);
307   CeedOperatorSetField(op_error, "true_soln", ceeddata[numlevels-1]->Erestrictui,
308                        CEED_NOTRANSPOSE, CEED_BASIS_COLLOCATED, target);
309   CeedOperatorSetField(op_error, "error", ceeddata[numlevels-1]->Erestrictui,
310                        CEED_NOTRANSPOSE, CEED_BASIS_COLLOCATED,
311                        CEED_VECTOR_ACTIVE);
312 
313   // Calculate multiplicity
314   for (int i=0; i<numlevels; i++) {
315     PetscScalar *x;
316 
317     // CEED vector
318     ierr = VecGetArray(Xloc[i], &x); CHKERRQ(ierr);
319     CeedVectorSetArray(ceeddata[i]->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
320 
321     // Multiplicity
322     CeedElemRestrictionGetMultiplicity(ceeddata[i]->Erestrictu,
323                                        ceeddata[i]->xceed);
324 
325     // Restore vector
326     ierr = VecRestoreArray(Xloc[i], &x); CHKERRQ(ierr);
327 
328     // Creat mult vector
329     ierr = VecDuplicate(Xloc[i], &mult[i]); CHKERRQ(ierr);
330 
331     // Local-to-global
332     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
333     ierr = DMLocalToGlobalBegin(dm[i], Xloc[i], ADD_VALUES, X[i]);
334     CHKERRQ(ierr);
335     ierr = DMLocalToGlobalEnd(dm[i], Xloc[i], ADD_VALUES, X[i]);
336     CHKERRQ(ierr);
337     ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr);
338 
339     // Global-to-local
340     ierr = DMGlobalToLocalBegin(dm[i], X[i], INSERT_VALUES, mult[i]);
341     CHKERRQ(ierr);
342     ierr = DMGlobalToLocalEnd(dm[i], X[i], INSERT_VALUES, mult[i]);
343     CHKERRQ(ierr);
344     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
345 
346     // Multiplicity scaling
347     ierr = VecReciprocal(mult[i]);
348   }
349 
350   // Set up Mat
351   for (int i=0; i<numlevels; i++) {
352     // User Operator
353     userO[i]->comm = comm;
354     userO[i]->dm = dm[i];
355     userO[i]->Xloc = Xloc[i];
356     ierr = VecDuplicate(Xloc[i], &userO[i]->Yloc); CHKERRQ(ierr);
357     userO[i]->xceed = ceeddata[i]->xceed;
358     userO[i]->yceed = ceeddata[i]->yceed;
359     userO[i]->op = ceeddata[i]->op_apply;
360     userO[i]->ceed = ceed;
361 
362     // Set up diagonal
363     const CeedScalar *ceedarray;
364     PetscScalar *petscarray;
365     CeedInt length;
366 
367     ierr = VecDuplicate(X[i], &userO[i]->diag); CHKERRQ(ierr);
368     ierr = VecDuplicate(Xloc[i], &diagloc); CHKERRQ(ierr);
369 
370     // -- Local diagonal
371     CeedOperatorAssembleLinearDiagonal(userO[i]->op, &diagceed,
372                                        CEED_REQUEST_IMMEDIATE);
373 
374     // -- Copy values
375     CeedVectorGetArrayRead(diagceed, CEED_MEM_HOST, &ceedarray);
376     ierr = VecGetArray(diagloc, &petscarray); CHKERRQ(ierr);
377     CeedVectorGetLength(diagceed, &length);
378     for (CeedInt i=0; i<length; i++)
379       petscarray[i] = ceedarray[i];
380     CeedVectorRestoreArrayRead(diagceed, &ceedarray);
381     ierr = VecRestoreArray(diagloc, &petscarray); CHKERRQ(ierr);
382 
383     // -- Global diagonal
384     ierr = VecZeroEntries(userO[i]->diag); CHKERRQ(ierr);
385     ierr = DMLocalToGlobalBegin(userO[i]->dm, diagloc, ADD_VALUES,
386                                 userO[i]->diag); CHKERRQ(ierr);
387     ierr = DMLocalToGlobalEnd(userO[i]->dm, diagloc, ADD_VALUES,
388                               userO[i]->diag); CHKERRQ(ierr);
389 
390     // -- Cleanup
391     ierr = VecDestroy(&diagloc); CHKERRQ(ierr);
392     CeedVectorDestroy(&diagceed);
393 
394     if (i > 0) {
395       // Interp Operator
396       userI[i]->comm = comm;
397       userI[i]->dmc = dm[i-1];
398       userI[i]->dmf = dm[i];
399       userI[i]->Xloc = Xloc[i-1];
400       userI[i]->Yloc = userO[i]->Yloc;
401       userI[i]->mult = mult[i];
402       userI[i]->ceedvecc = userO[i-1]->xceed;
403       userI[i]->ceedvecf = userO[i]->yceed;
404       userI[i]->op = ceeddata[i]->op_interp;
405       userI[i]->ceed = ceed;
406 
407       // Restrict Operator
408       userR[i]->comm = comm;
409       userR[i]->dmc = dm[i-1];
410       userR[i]->dmf = dm[i];
411       userR[i]->Xloc = Xloc[i];
412       userR[i]->Yloc = userO[i-1]->Yloc;
413       userR[i]->mult = mult[i];
414       userR[i]->ceedvecf = userO[i]->xceed;
415       userR[i]->ceedvecc = userO[i-1]->yceed;
416       userR[i]->op = ceeddata[i]->op_restrict;
417       userR[i]->ceed = ceed;
418     }
419   }
420 
421   // Set up KSP
422   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
423   {
424     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
425     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
426     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
427                             PETSC_DEFAULT); CHKERRQ(ierr);
428   }
429   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
430   ierr = KSPSetOperators(ksp, matO[numlevels-1], matO[numlevels-1]);
431   CHKERRQ(ierr);
432 
433   // Set up PCMG
434   ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
435   PCMGCycleType pcgmcycletype = PC_MG_CYCLE_V;
436   {
437     ierr = PCSetType(pc, PCMG); CHKERRQ(ierr);
438 
439     // PCMG levels
440     ierr = PCMGSetLevels(pc, numlevels, NULL); CHKERRQ(ierr);
441     for (int i=0; i<numlevels; i++) {
442       // Smoother
443       KSP smoother;
444       PC smoother_pc;
445       ierr = PCMGGetSmoother(pc, i, &smoother); CHKERRQ(ierr);
446       ierr = KSPSetType(smoother, KSPCHEBYSHEV); CHKERRQ(ierr);
447       ierr = KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1); CHKERRQ(ierr);
448       ierr = KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE); CHKERRQ(ierr);
449       ierr = KSPSetOperators(smoother, matO[i], matO[i]); CHKERRQ(ierr);
450       ierr = KSPGetPC(smoother, &smoother_pc); CHKERRQ(ierr);
451       ierr = PCSetType(smoother_pc, PCJACOBI); CHKERRQ(ierr);
452       ierr = PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
453 
454       // Work vector
455       if (i < numlevels-1) {
456         ierr = PCMGSetX(pc, i, X[i]); CHKERRQ(ierr);
457       }
458 
459       // Level transfers
460       if (i > 0) {
461         // Interpolation
462         ierr = PCMGSetInterpolation(pc, i, matI[i]); CHKERRQ(ierr);
463 
464         // Restriction
465         ierr = PCMGSetRestriction(pc, i, matR[i]); CHKERRQ(ierr);
466       }
467 
468       // Coarse solve
469       KSP coarse;
470       PC coarse_pc;
471       ierr = PCMGGetCoarseSolve(pc, &coarse); CHKERRQ(ierr);
472       ierr = KSPSetType(coarse, KSPCG); CHKERRQ(ierr);
473       ierr = KSPSetOperators(coarse, matO[0], matO[0]); CHKERRQ(ierr);
474       ierr = KSPSetTolerances(coarse, 1e-10, 1e-10, PETSC_DEFAULT,
475                               PETSC_DEFAULT); CHKERRQ(ierr);
476       ierr = KSPGetPC(coarse, &coarse_pc); CHKERRQ(ierr);
477       ierr = PCSetType(coarse_pc, PCJACOBI); CHKERRQ(ierr);
478       ierr = PCJacobiSetType(coarse_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
479     }
480 
481     // PCMG options
482     ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr);
483     ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr);
484     ierr = PCMGSetCycleType(pc, pcgmcycletype); CHKERRQ(ierr);
485   }
486 
487   // First run, if benchmarking
488   if (benchmark_mode) {
489     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
490     CHKERRQ(ierr);
491     ierr = VecZeroEntries(X[numlevels-1]); CHKERRQ(ierr);
492     my_rt_start = MPI_Wtime();
493     ierr = KSPSolve(ksp, rhs, X[numlevels-1]); CHKERRQ(ierr);
494     my_rt = MPI_Wtime() - my_rt_start;
495     ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
496     CHKERRQ(ierr);
497     // Set maxits based on first iteration timing
498     if (my_rt > 0.02) {
499       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
500       CHKERRQ(ierr);
501     } else {
502       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
503       CHKERRQ(ierr);
504     }
505   }
506 
507   // Timed solve
508   ierr = VecZeroEntries(X[numlevels-1]); CHKERRQ(ierr);
509   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
510   my_rt_start = MPI_Wtime();
511   ierr = KSPSolve(ksp, rhs, X[numlevels-1]); CHKERRQ(ierr);
512   my_rt = MPI_Wtime() - my_rt_start;
513 
514   // Output results
515   {
516     KSPType ksptype;
517     PCMGType pcmgtype;
518     KSPConvergedReason reason;
519     PetscReal rnorm;
520     PetscInt its;
521     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
522     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
523     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
524     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
525     ierr = PCMGGetType(pc, &pcmgtype); CHKERRQ(ierr);
526     if (!test_mode || reason < 0 || rnorm > 1e-8) {
527       ierr = PetscPrintf(comm,
528                          "  KSP:\n"
529                          "    KSP Type                           : %s\n"
530                          "    KSP Convergence                    : %s\n"
531                          "    Total KSP Iterations               : %D\n"
532                          "    Final rnorm                        : %e\n",
533                          ksptype, KSPConvergedReasons[reason], its,
534                          (double)rnorm); CHKERRQ(ierr);
535       ierr = PetscPrintf(comm,
536                          "  PCMG:\n"
537                          "    PCMG Type                          : %s\n"
538                          "    PCMG Cycle Type                    : %s\n",
539                          PCMGTypes[pcmgtype],
540                          PCMGCycleTypes[pcgmcycletype]); CHKERRQ(ierr);
541     }
542     if (!test_mode) {
543       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
544     }
545     {
546       PetscReal maxerror;
547       ierr = ComputeErrorMax(userO[numlevels-1], op_error, X[numlevels-1], target,
548                              &maxerror); CHKERRQ(ierr);
549       PetscReal tol = 5e-2;
550       if (!test_mode || maxerror > tol) {
551         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
552         CHKERRQ(ierr);
553         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
554         CHKERRQ(ierr);
555         ierr = PetscPrintf(comm,
556                            "    Pointwise Error (max)              : %e\n"
557                            "    CG Solve Time                      : %g (%g) sec\n",
558                            (double)maxerror, rt_max, rt_min); CHKERRQ(ierr);
559       }
560     }
561     if (benchmark_mode && (!test_mode)) {
562       ierr = PetscPrintf(comm,
563                          "    DoFs/Sec in CG                     : %g (%g) million\n",
564                          1e-6*gsize[numlevels-1]*its/rt_max,
565                          1e-6*gsize[numlevels-1]*its/rt_min);
566       CHKERRQ(ierr);
567     }
568   }
569 
570   if (write_solution) {
571     PetscViewer vtkviewersoln;
572 
573     ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr);
574     ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr);
575     ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr);
576     ierr = VecView(X[numlevels-1], vtkviewersoln); CHKERRQ(ierr);
577     ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr);
578   }
579 
580   // Cleanup
581   for (int i=0; i<numlevels; i++) {
582     ierr = VecDestroy(&X[i]); CHKERRQ(ierr);
583     ierr = VecDestroy(&Xloc[i]); CHKERRQ(ierr);
584     ierr = VecDestroy(&mult[i]); CHKERRQ(ierr);
585     ierr = VecDestroy(&userO[i]->Yloc); CHKERRQ(ierr);
586     ierr = VecDestroy(&userO[i]->diag); CHKERRQ(ierr);
587     ierr = MatDestroy(&matO[i]); CHKERRQ(ierr);
588     ierr = PetscFree(userO[i]); CHKERRQ(ierr);
589     if (i > 0) {
590       ierr = MatDestroy(&matI[i]); CHKERRQ(ierr);
591       ierr = PetscFree(userI[i]); CHKERRQ(ierr);
592       ierr = MatDestroy(&matR[i]); CHKERRQ(ierr);
593       ierr = PetscFree(userR[i]); CHKERRQ(ierr);
594     }
595     ierr = CeedDataDestroy(i, ceeddata[i]); CHKERRQ(ierr);
596     ierr = DMDestroy(&dm[i]); CHKERRQ(ierr);
597   }
598   ierr = PetscFree(leveldegrees); CHKERRQ(ierr);
599   ierr = PetscFree(dm); CHKERRQ(ierr);
600   ierr = PetscFree(X); CHKERRQ(ierr);
601   ierr = PetscFree(Xloc); CHKERRQ(ierr);
602   ierr = PetscFree(mult); CHKERRQ(ierr);
603   ierr = PetscFree(matO); CHKERRQ(ierr);
604   ierr = PetscFree(matI); CHKERRQ(ierr);
605   ierr = PetscFree(matR); CHKERRQ(ierr);
606   ierr = PetscFree(ceeddata); CHKERRQ(ierr);
607   ierr = PetscFree(userO); CHKERRQ(ierr);
608   ierr = PetscFree(userI); CHKERRQ(ierr);
609   ierr = PetscFree(userR); CHKERRQ(ierr);
610   ierr = PetscFree(lsize); CHKERRQ(ierr);
611   ierr = PetscFree(xlsize); CHKERRQ(ierr);
612   ierr = PetscFree(gsize); CHKERRQ(ierr);
613   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
614   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
615   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
616   ierr = DMDestroy(&dmOrig); CHKERRQ(ierr);
617   CeedVectorDestroy(&target);
618   CeedQFunctionDestroy(&qf_error);
619   CeedQFunctionDestroy(&qf_restrict);
620   CeedQFunctionDestroy(&qf_prolong);
621   CeedOperatorDestroy(&op_error);
622   CeedDestroy(&ceed);
623   return PetscFinalize();
624 }
625