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