xref: /libCEED/examples/petsc/multigrid.c (revision 9396343da0997500e192eb5de23b22627660fffe)
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   CeedMemType memtyperequested;
65   CeedVector rhsceed, target;
66   CeedQFunction qferror, qfrestrict, qfprolong;
67   CeedOperator operror;
68   bpType bpchoice;
69   coarsenType coarsen;
70 
71   ierr = PetscInitialize(&argc, &argv, NULL, help);
72   if (ierr) return ierr;
73   comm = PETSC_COMM_WORLD;
74 
75   // Check PETSc CUDA avaliability
76   PetscBool petschavecuda, setmemtyperequest = PETSC_FALSE;
77   // *INDENT-OFF*
78   #ifdef PETSC_HAVE_CUDA
79   petschavecuda = PETSC_TRUE;
80   #else
81   petschavecuda = PETSC_FALSE;
82   #endif
83   // *INDENT-ON*
84 
85   // Parse command line options
86   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr);
87   bpchoice = CEED_BP3;
88   ierr = PetscOptionsEnum("-problem",
89                           "CEED benchmark problem to solve", NULL,
90                           bpTypes, (PetscEnum)bpchoice, (PetscEnum *)&bpchoice,
91                           NULL); CHKERRQ(ierr);
92   ncompu = bpOptions[bpchoice].ncompu;
93   test_mode = PETSC_FALSE;
94   ierr = PetscOptionsBool("-test",
95                           "Testing mode (do not print unless error is large)",
96                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
97   benchmark_mode = PETSC_FALSE;
98   ierr = PetscOptionsBool("-benchmark",
99                           "Benchmarking mode (prints benchmark statistics)",
100                           NULL, benchmark_mode, &benchmark_mode, NULL);
101   CHKERRQ(ierr);
102   write_solution = PETSC_FALSE;
103   ierr = PetscOptionsBool("-write_solution",
104                           "Write solution for visualization",
105                           NULL, write_solution, &write_solution, NULL);
106   CHKERRQ(ierr);
107   degree = test_mode ? 3 : 2;
108   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
109                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
110   if (degree < 1) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
111                              "-degree %D must be at least 1", degree);
112   qextra = bpOptions[bpchoice].qextra;
113   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
114                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
115   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
116                             NULL, ceedresource, ceedresource,
117                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
118   coarsen = COARSEN_UNIFORM;
119   ierr = PetscOptionsEnum("-coarsen",
120                           "Coarsening strategy to use", NULL,
121                           coarsenTypes, (PetscEnum)coarsen,
122                           (PetscEnum *)&coarsen, NULL); CHKERRQ(ierr);
123   read_mesh = PETSC_FALSE;
124   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL,
125                             filename, filename, sizeof(filename), &read_mesh);
126   CHKERRQ(ierr);
127   if (!read_mesh) {
128     PetscInt tmp = dim;
129     ierr = PetscOptionsIntArray("-cells","Number of cells per dimension", NULL,
130                                 melem, &tmp, NULL); CHKERRQ(ierr);
131   }
132   memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST;
133   ierr = PetscOptionsEnum("-memtype",
134                           "CEED MemType requested", NULL,
135                           memTypes, (PetscEnum)memtyperequested,
136                           (PetscEnum *)&memtyperequested, &setmemtyperequest);
137   CHKERRQ(ierr);
138   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
139 
140   // Set up libCEED
141   CeedInit(ceedresource, &ceed);
142   CeedMemType memtypebackend;
143   CeedGetPreferredMemType(ceed, &memtypebackend);
144 
145   // Check memtype compatibility
146   if (!setmemtyperequest)
147     memtyperequested = memtypebackend;
148   else if (!petschavecuda && memtyperequested == CEED_MEM_DEVICE)
149     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS,
150              "PETSc was not built with CUDA. "
151              "Requested MemType CEED_MEM_DEVICE is not supported.", NULL);
152 
153   // Setup DM
154   if (read_mesh) {
155     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dmorig);
156     CHKERRQ(ierr);
157   } else {
158     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_FALSE, melem, NULL,
159                                NULL, NULL, PETSC_TRUE,&dmorig); CHKERRQ(ierr);
160   }
161 
162   {
163     DM dmDist = NULL;
164     PetscPartitioner part;
165 
166     ierr = DMPlexGetPartitioner(dmorig, &part); CHKERRQ(ierr);
167     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
168     ierr = DMPlexDistribute(dmorig, 0, NULL, &dmDist); CHKERRQ(ierr);
169     if (dmDist) {
170       ierr = DMDestroy(&dmorig); CHKERRQ(ierr);
171       dmorig = dmDist;
172     }
173   }
174 
175   // Allocate arrays for PETSc objects for each level
176   switch (coarsen) {
177   case COARSEN_UNIFORM:
178     numlevels = degree;
179     break;
180   case COARSEN_LOGARITHMIC:
181     numlevels = ceil(log(degree)/log(2)) + 1;
182     break;
183   }
184   ierr = PetscMalloc1(numlevels, &leveldegrees); CHKERRQ(ierr);
185   fineLevel = numlevels - 1;
186 
187   switch (coarsen) {
188   case COARSEN_UNIFORM:
189     for (int i=0; i<numlevels; i++) leveldegrees[i] = i + 1;
190     break;
191   case COARSEN_LOGARITHMIC:
192     for (int i=0; i<numlevels - 1; i++) leveldegrees[i] = pow(2,i);
193     leveldegrees[fineLevel] = degree;
194     break;
195   }
196   ierr = PetscMalloc1(numlevels, &dm); CHKERRQ(ierr);
197   ierr = PetscMalloc1(numlevels, &X); CHKERRQ(ierr);
198   ierr = PetscMalloc1(numlevels, &Xloc); CHKERRQ(ierr);
199   ierr = PetscMalloc1(numlevels, &mult); CHKERRQ(ierr);
200   ierr = PetscMalloc1(numlevels, &userO); CHKERRQ(ierr);
201   ierr = PetscMalloc1(numlevels, &userPR); CHKERRQ(ierr);
202   ierr = PetscMalloc1(numlevels, &matO); CHKERRQ(ierr);
203   ierr = PetscMalloc1(numlevels, &matPR); CHKERRQ(ierr);
204   ierr = PetscMalloc1(numlevels, &lsize); CHKERRQ(ierr);
205   ierr = PetscMalloc1(numlevels, &xlsize); CHKERRQ(ierr);
206   ierr = PetscMalloc1(numlevels, &gsize); CHKERRQ(ierr);
207 
208   // Setup DM and Operator Mat Shells for each level
209   for (CeedInt i=0; i<numlevels; i++) {
210     // Create DM
211     ierr = DMClone(dmorig, &dm[i]); CHKERRQ(ierr);
212     ierr = SetupDMByDegree(dm[i], leveldegrees[i], ncompu, bpchoice);
213     CHKERRQ(ierr);
214 
215     // Create vectors
216     if (memtyperequested == CEED_MEM_DEVICE) {
217       ierr = DMSetVecType(dm[i], VECCUDA); CHKERRQ(ierr);
218     }
219     ierr = DMCreateGlobalVector(dm[i], &X[i]); CHKERRQ(ierr);
220     ierr = VecGetLocalSize(X[i], &lsize[i]); CHKERRQ(ierr);
221     ierr = VecGetSize(X[i], &gsize[i]); CHKERRQ(ierr);
222     ierr = DMCreateLocalVector(dm[i], &Xloc[i]); CHKERRQ(ierr);
223     ierr = VecGetSize(Xloc[i], &xlsize[i]); CHKERRQ(ierr);
224 
225     // Operator
226     ierr = PetscMalloc1(1, &userO[i]); CHKERRQ(ierr);
227     ierr = MatCreateShell(comm, lsize[i], lsize[i], gsize[i], gsize[i],
228                           userO[i], &matO[i]); CHKERRQ(ierr);
229     ierr = MatShellSetOperation(matO[i], MATOP_MULT,
230                                 (void(*)(void))MatMult_Ceed); CHKERRQ(ierr);
231     ierr = MatShellSetOperation(matO[i], MATOP_GET_DIAGONAL,
232                                 (void(*)(void))MatGetDiag); CHKERRQ(ierr);
233     if (memtyperequested == CEED_MEM_DEVICE) {
234       ierr = MatShellSetVecType(matO[i], VECCUDA); CHKERRQ(ierr);
235     }
236 
237     // Level transfers
238     if (i > 0) {
239       // Interp
240       ierr = PetscMalloc1(1, &userPR[i]); CHKERRQ(ierr);
241       ierr = MatCreateShell(comm, lsize[i], lsize[i-1], gsize[i], gsize[i-1],
242                             userPR[i], &matPR[i]); CHKERRQ(ierr);
243       ierr = MatShellSetOperation(matPR[i], MATOP_MULT,
244                                   (void(*)(void))MatMult_Prolong);
245       CHKERRQ(ierr);
246       ierr = MatShellSetOperation(matPR[i], MATOP_MULT_TRANSPOSE,
247                                   (void(*)(void))MatMult_Restrict);
248       CHKERRQ(ierr);
249       if (memtyperequested == CEED_MEM_DEVICE) {
250         ierr = MatShellSetVecType(matPR[i], VECCUDA); CHKERRQ(ierr);
251       }
252     }
253   }
254   ierr = VecDuplicate(X[fineLevel], &rhs); CHKERRQ(ierr);
255 
256   // Set up libCEED
257   CeedInit(ceedresource, &ceed);
258 
259   // Print global grid information
260   if (!test_mode) {
261     PetscInt P = degree + 1, Q = P + qextra;
262 
263     const char *usedresource;
264     CeedGetResource(ceed, &usedresource);
265 
266     VecType vectype;
267     ierr = VecGetType(X[0], &vectype); CHKERRQ(ierr);
268 
269     ierr = PetscPrintf(comm,
270                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc + PCMG --\n"
271                        "  PETSc:\n"
272                        "    PETSc Vec Type                     : %s\n"
273                        "  libCEED:\n"
274                        "    libCEED Backend                    : %s\n"
275                        "    libCEED Backend MemType            : %s\n"
276                        "    libCEED User Requested MemType     : %s\n"
277                        "  Mesh:\n"
278                        "    Number of 1D Basis Nodes (p)       : %d\n"
279                        "    Number of 1D Quadrature Points (q) : %d\n"
280                        "    Global Nodes                       : %D\n"
281                        "    Owned Nodes                        : %D\n"
282                        "    DoF per node                       : %D\n"
283                        "  Multigrid:\n"
284                        "    Number of Levels                   : %d\n",
285                        bpchoice+1, vectype, usedresource,
286                        CeedMemTypes[memtypebackend],
287                        (setmemtyperequest) ?
288                        CeedMemTypes[memtyperequested] : "none",
289                        P, Q, gsize[fineLevel]/ncompu, lsize[fineLevel]/ncompu,
290                        ncompu, numlevels); CHKERRQ(ierr);
291   }
292 
293   // Create RHS vector
294   ierr = VecDuplicate(Xloc[fineLevel], &rhsloc); CHKERRQ(ierr);
295   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
296   if (memtyperequested == CEED_MEM_HOST) {
297     ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr);
298   } else {
299     ierr = VecCUDAGetArray(rhsloc, &r); CHKERRQ(ierr);
300   }
301   CeedVectorCreate(ceed, xlsize[fineLevel], &rhsceed);
302   CeedVectorSetArray(rhsceed, CEED_MEM_HOST, CEED_USE_POINTER, r);
303 
304   // Set up libCEED operators on each level
305   ierr = PetscMalloc1(numlevels, &ceeddata); CHKERRQ(ierr);
306   for (int i=0; i<numlevels; i++) {
307     // Print level information
308     if (!test_mode && (i == 0 || i == fineLevel)) {
309       ierr = PetscPrintf(comm,"    Level %D (%s):\n"
310                          "      Number of 1D Basis Nodes (p)     : %d\n"
311                          "      Global Nodes                     : %D\n"
312                          "      Owned Nodes                      : %D\n",
313                          i, (i? "fine" : "coarse"), leveldegrees[i] + 1,
314                          gsize[i]/ncompu, lsize[i]/ncompu); CHKERRQ(ierr);
315     }
316     ierr = PetscMalloc1(1, &ceeddata[i]); CHKERRQ(ierr);
317     ierr = SetupLibceedByDegree(dm[i], ceed, leveldegrees[i], dim, qextra,
318                                 ncompu, gsize[i], xlsize[i], bpchoice,
319                                 ceeddata[i], i==(fineLevel), rhsceed,
320                                 &target); CHKERRQ(ierr);
321   }
322 
323   // Gather RHS
324   CeedVectorSyncArray(rhsceed, memtyperequested);
325   if (memtyperequested == CEED_MEM_HOST) {
326     ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr);
327   } else {
328     ierr = VecCUDARestoreArray(rhsloc, &r); CHKERRQ(ierr);
329   }
330   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
331   ierr = DMLocalToGlobal(dm[fineLevel], rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr);
332   CeedVectorDestroy(&rhsceed);
333 
334   // Create the restriction/interpolation Q-function
335   CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_NONE, CEED_EVAL_INTERP,
336                               &qfrestrict);
337   CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_INTERP, CEED_EVAL_NONE,
338                               &qfprolong);
339 
340   // Set up libCEED level transfer operators
341   ierr = CeedLevelTransferSetup(ceed, numlevels, ncompu, bpchoice, ceeddata,
342                                 leveldegrees, qfrestrict, qfprolong);
343   CHKERRQ(ierr);
344 
345   // Create the error Q-function
346   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].error,
347                               bpOptions[bpchoice].errorfname, &qferror);
348   CeedQFunctionAddInput(qferror, "u", ncompu, CEED_EVAL_INTERP);
349   CeedQFunctionAddInput(qferror, "true_soln", ncompu, CEED_EVAL_NONE);
350   CeedQFunctionAddOutput(qferror, "error", ncompu, CEED_EVAL_NONE);
351 
352   // Create the error operator
353   CeedOperatorCreate(ceed, qferror, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
354                      &operror);
355   CeedOperatorSetField(operror, "u", ceeddata[fineLevel]->Erestrictu,
356                        ceeddata[fineLevel]->basisu, CEED_VECTOR_ACTIVE);
357   CeedOperatorSetField(operror, "true_soln", ceeddata[fineLevel]->Erestrictui,
358                        CEED_BASIS_COLLOCATED, target);
359   CeedOperatorSetField(operror, "error", ceeddata[fineLevel]->Erestrictui,
360                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
361 
362   // Calculate multiplicity
363   for (int i=0; i<numlevels; i++) {
364     PetscScalar *x;
365 
366     // CEED vector
367     ierr = VecGetArray(Xloc[i], &x); CHKERRQ(ierr);
368     CeedVectorSetArray(ceeddata[i]->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
369 
370     // Multiplicity
371     CeedElemRestrictionGetMultiplicity(ceeddata[i]->Erestrictu,
372                                        ceeddata[i]->xceed);
373 
374     // Restore vector
375     ierr = VecRestoreArray(Xloc[i], &x); CHKERRQ(ierr);
376 
377     // Creat mult vector
378     ierr = VecDuplicate(Xloc[i], &mult[i]); CHKERRQ(ierr);
379 
380     // Local-to-global
381     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
382     ierr = DMLocalToGlobal(dm[i], Xloc[i], ADD_VALUES, X[i]);
383     CHKERRQ(ierr);
384     ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr);
385 
386     // Global-to-local
387     ierr = DMGlobalToLocal(dm[i], X[i], INSERT_VALUES, mult[i]);
388     CHKERRQ(ierr);
389     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
390 
391     // Multiplicity scaling
392     ierr = VecReciprocal(mult[i]);
393   }
394 
395   // Set up Mat
396   for (int i=0; i<numlevels; i++) {
397     // User Operator
398     userO[i]->comm = comm;
399     userO[i]->dm = dm[i];
400     userO[i]->Xloc = Xloc[i];
401     ierr = VecDuplicate(Xloc[i], &userO[i]->Yloc); CHKERRQ(ierr);
402     userO[i]->xceed = ceeddata[i]->xceed;
403     userO[i]->yceed = ceeddata[i]->yceed;
404     userO[i]->op = ceeddata[i]->opapply;
405     userO[i]->ceed = ceed;
406     userO[i]->memtype = memtyperequested;
407     if (memtyperequested == CEED_MEM_HOST) {
408       userO[i]->VecGetArray = VecGetArray;
409       userO[i]->VecGetArrayRead = VecGetArrayRead;
410       userO[i]->VecRestoreArray = VecRestoreArray;
411       userO[i]->VecRestoreArrayRead = VecRestoreArrayRead;
412     } else {
413       userO[i]->VecGetArray = VecCUDAGetArray;
414       userO[i]->VecGetArrayRead = VecCUDAGetArrayRead;
415       userO[i]->VecRestoreArray = VecCUDARestoreArray;
416       userO[i]->VecRestoreArrayRead = VecCUDARestoreArrayRead;
417     }
418 
419     if (i > 0) {
420       // Prolongation/Restriction Operator
421       userPR[i]->comm = comm;
422       userPR[i]->dmf = dm[i];
423       userPR[i]->dmc = dm[i-1];
424       userPR[i]->locvecc = Xloc[i-1];
425       userPR[i]->locvecf = userO[i]->Yloc;
426       userPR[i]->multvec = mult[i];
427       userPR[i]->ceedvecc = userO[i-1]->xceed;
428       userPR[i]->ceedvecf = userO[i]->yceed;
429       userPR[i]->opprolong = ceeddata[i]->opprolong;
430       userPR[i]->oprestrict = ceeddata[i]->oprestrict;
431       userPR[i]->ceed = ceed;
432       userPR[i]->memtype = userO[i]->memtype;
433       userPR[i]->VecGetArray = userO[i]->VecGetArray;
434       userPR[i]->VecGetArrayRead = userO[i]->VecGetArrayRead;
435       userPR[i]->VecRestoreArray = userO[i]->VecRestoreArray;
436       userPR[i]->VecRestoreArrayRead = userO[i]->VecRestoreArrayRead;
437     }
438   }
439 
440   // Setup dummy SNES for AMG coarse solve
441   ierr = SNESCreate(comm, &snesdummy); CHKERRQ(ierr);
442   ierr = SNESSetDM(snesdummy, dm[0]); CHKERRQ(ierr);
443   ierr = SNESSetSolution(snesdummy, X[0]); CHKERRQ(ierr);
444 
445   // -- Jacobian matrix
446   ierr = DMSetMatType(dm[0], MATAIJ); CHKERRQ(ierr);
447   ierr = DMCreateMatrix(dm[0], &matcoarse); CHKERRQ(ierr);
448   ierr = SNESSetJacobian(snesdummy, matcoarse, matcoarse, NULL,
449                          NULL); CHKERRQ(ierr);
450 
451   // -- Residual evaluation function
452   ierr = SNESSetFunction(snesdummy, X[0], FormResidual_Ceed,
453                          userO[0]); CHKERRQ(ierr);
454 
455   // -- Form Jacobian
456   ierr = SNESComputeJacobianDefaultColor(snesdummy, X[0], matO[0],
457                                          matcoarse, NULL); CHKERRQ(ierr);
458 
459   // Set up KSP
460   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
461   {
462     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
463     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
464     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
465                             PETSC_DEFAULT); CHKERRQ(ierr);
466   }
467   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
468   ierr = KSPSetOperators(ksp, matO[fineLevel], matO[fineLevel]);
469   CHKERRQ(ierr);
470 
471   // Set up PCMG
472   ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
473   PCMGCycleType pcgmcycletype = PC_MG_CYCLE_V;
474   {
475     ierr = PCSetType(pc, PCMG); CHKERRQ(ierr);
476 
477     // PCMG levels
478     ierr = PCMGSetLevels(pc, numlevels, NULL); CHKERRQ(ierr);
479     for (int i=0; i<numlevels; i++) {
480       // Smoother
481       KSP smoother;
482       PC smoother_pc;
483       ierr = PCMGGetSmoother(pc, i, &smoother); CHKERRQ(ierr);
484       ierr = KSPSetType(smoother, KSPCHEBYSHEV); CHKERRQ(ierr);
485       ierr = KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1); CHKERRQ(ierr);
486       ierr = KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE); CHKERRQ(ierr);
487       ierr = KSPSetOperators(smoother, matO[i], matO[i]); CHKERRQ(ierr);
488       ierr = KSPGetPC(smoother, &smoother_pc); CHKERRQ(ierr);
489       ierr = PCSetType(smoother_pc, PCJACOBI); CHKERRQ(ierr);
490       ierr = PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
491 
492       // Work vector
493       if (i < numlevels - 1) {
494         ierr = PCMGSetX(pc, i, X[i]); CHKERRQ(ierr);
495       }
496 
497       // Level transfers
498       if (i > 0) {
499         // Interpolation
500         ierr = PCMGSetInterpolation(pc, i, matPR[i]); CHKERRQ(ierr);
501       }
502 
503       // Coarse solve
504       KSP coarse;
505       PC coarse_pc;
506       ierr = PCMGGetCoarseSolve(pc, &coarse); CHKERRQ(ierr);
507       ierr = KSPSetType(coarse, KSPPREONLY); CHKERRQ(ierr);
508       ierr = KSPSetOperators(coarse, matcoarse, matcoarse); CHKERRQ(ierr);
509 
510       ierr = KSPGetPC(coarse, &coarse_pc); CHKERRQ(ierr);
511       ierr = PCSetType(coarse_pc, PCGAMG); CHKERRQ(ierr);
512 
513       ierr = KSPSetOptionsPrefix(coarse, "coarse_"); CHKERRQ(ierr);
514       ierr = PCSetOptionsPrefix(coarse_pc, "coarse_"); CHKERRQ(ierr);
515       ierr = KSPSetFromOptions(coarse); CHKERRQ(ierr);
516       ierr = PCSetFromOptions(coarse_pc); CHKERRQ(ierr);
517     }
518 
519     // PCMG options
520     ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr);
521     ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr);
522     ierr = PCMGSetCycleType(pc, pcgmcycletype); CHKERRQ(ierr);
523   }
524 
525   // First run, if benchmarking
526   if (benchmark_mode) {
527     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
528     CHKERRQ(ierr);
529     ierr = VecZeroEntries(X[fineLevel]); CHKERRQ(ierr);
530     my_rt_start = MPI_Wtime();
531     ierr = KSPSolve(ksp, rhs, X[fineLevel]); CHKERRQ(ierr);
532     my_rt = MPI_Wtime() - my_rt_start;
533     ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
534     CHKERRQ(ierr);
535     // Set maxits based on first iteration timing
536     if (my_rt > 0.02) {
537       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
538       CHKERRQ(ierr);
539     } else {
540       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
541       CHKERRQ(ierr);
542     }
543   }
544 
545   // Timed solve
546   ierr = VecZeroEntries(X[fineLevel]); CHKERRQ(ierr);
547   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
548   my_rt_start = MPI_Wtime();
549   ierr = KSPSolve(ksp, rhs, X[fineLevel]); CHKERRQ(ierr);
550   my_rt = MPI_Wtime() - my_rt_start;
551 
552   // Output results
553   {
554     KSPType ksptype;
555     PCMGType pcmgtype;
556     KSPConvergedReason reason;
557     PetscReal rnorm;
558     PetscInt its;
559     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
560     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
561     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
562     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
563     ierr = PCMGGetType(pc, &pcmgtype); CHKERRQ(ierr);
564     if (!test_mode || reason < 0 || rnorm > 1e-8) {
565       ierr = PetscPrintf(comm,
566                          "  KSP:\n"
567                          "    KSP Type                           : %s\n"
568                          "    KSP Convergence                    : %s\n"
569                          "    Total KSP Iterations               : %D\n"
570                          "    Final rnorm                        : %e\n",
571                          ksptype, KSPConvergedReasons[reason], its,
572                          (double)rnorm); CHKERRQ(ierr);
573       ierr = PetscPrintf(comm,
574                          "  PCMG:\n"
575                          "    PCMG Type                          : %s\n"
576                          "    PCMG Cycle Type                    : %s\n",
577                          PCMGTypes[pcmgtype],
578                          PCMGCycleTypes[pcgmcycletype]); CHKERRQ(ierr);
579     }
580     if (!test_mode) {
581       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
582     }
583     {
584       PetscReal maxerror;
585       ierr = ComputeErrorMax(userO[fineLevel], operror, X[fineLevel], target,
586                              &maxerror); CHKERRQ(ierr);
587       PetscReal tol = 5e-2;
588       if (!test_mode || maxerror > tol) {
589         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
590         CHKERRQ(ierr);
591         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
592         CHKERRQ(ierr);
593         ierr = PetscPrintf(comm,
594                            "    Pointwise Error (max)              : %e\n"
595                            "    CG Solve Time                      : %g (%g) sec\n",
596                            (double)maxerror, rt_max, rt_min); CHKERRQ(ierr);
597       }
598     }
599     if (benchmark_mode && (!test_mode)) {
600       ierr = PetscPrintf(comm,
601                          "    DoFs/Sec in CG                     : %g (%g) million\n",
602                          1e-6*gsize[fineLevel]*its/rt_max,
603                          1e-6*gsize[fineLevel]*its/rt_min);
604       CHKERRQ(ierr);
605     }
606   }
607 
608   if (write_solution) {
609     PetscViewer vtkviewersoln;
610 
611     ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr);
612     ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr);
613     ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr);
614     ierr = VecView(X[fineLevel], vtkviewersoln); CHKERRQ(ierr);
615     ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr);
616   }
617 
618   // Cleanup
619   for (int i=0; i<numlevels; i++) {
620     ierr = VecDestroy(&X[i]); CHKERRQ(ierr);
621     ierr = VecDestroy(&Xloc[i]); CHKERRQ(ierr);
622     ierr = VecDestroy(&mult[i]); CHKERRQ(ierr);
623     ierr = VecDestroy(&userO[i]->Yloc); CHKERRQ(ierr);
624     ierr = MatDestroy(&matO[i]); CHKERRQ(ierr);
625     ierr = PetscFree(userO[i]); CHKERRQ(ierr);
626     if (i > 0) {
627       ierr = MatDestroy(&matPR[i]); CHKERRQ(ierr);
628       ierr = PetscFree(userPR[i]); CHKERRQ(ierr);
629     }
630     ierr = CeedDataDestroy(i, ceeddata[i]); CHKERRQ(ierr);
631     ierr = DMDestroy(&dm[i]); CHKERRQ(ierr);
632   }
633   ierr = PetscFree(leveldegrees); CHKERRQ(ierr);
634   ierr = PetscFree(dm); CHKERRQ(ierr);
635   ierr = PetscFree(X); CHKERRQ(ierr);
636   ierr = PetscFree(Xloc); CHKERRQ(ierr);
637   ierr = PetscFree(mult); CHKERRQ(ierr);
638   ierr = PetscFree(matO); CHKERRQ(ierr);
639   ierr = PetscFree(matPR); CHKERRQ(ierr);
640   ierr = PetscFree(ceeddata); CHKERRQ(ierr);
641   ierr = PetscFree(userO); CHKERRQ(ierr);
642   ierr = PetscFree(userPR); CHKERRQ(ierr);
643   ierr = PetscFree(lsize); CHKERRQ(ierr);
644   ierr = PetscFree(xlsize); CHKERRQ(ierr);
645   ierr = PetscFree(gsize); CHKERRQ(ierr);
646   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
647   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
648   ierr = MatDestroy(&matcoarse); CHKERRQ(ierr);
649   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
650   ierr = SNESDestroy(&snesdummy); CHKERRQ(ierr);
651   ierr = DMDestroy(&dmorig); CHKERRQ(ierr);
652   CeedVectorDestroy(&target);
653   CeedQFunctionDestroy(&qferror);
654   CeedQFunctionDestroy(&qfrestrict);
655   CeedQFunctionDestroy(&qfprolong);
656   CeedOperatorDestroy(&operror);
657   CeedDestroy(&ceed);
658   return PetscFinalize();
659 }
660