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