xref: /libCEED/examples/petsc/multigrid.c (revision de1229c50b0fa287ef0d97965becc5ec55ca5a76)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 //                        libCEED + PETSc Example: CEED BPs 3-6 with Multigrid
9 //
10 // This example demonstrates a simple usage of libCEED with PETSc to solve the
11 // CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps.
12 //
13 // The code uses higher level communication protocols in DMPlex.
14 //
15 // Build with:
16 //
17 //     make multigrid [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
18 //
19 // Sample runs:
20 //
21 //     multigrid -problem bp3
22 //     multigrid -problem bp4
23 //     multigrid -problem bp5 -ceed /cpu/self
24 //     multigrid -problem bp6 -ceed /gpu/cuda
25 //
26 //TESTARGS -ceed {ceed_resource} -test -problem bp3 -degree 3
27 
28 /// @file
29 /// CEED BPs 1-6 multigrid example using PETSc
30 const char help[] = "Solve CEED BPs using p-multigrid with PETSc and DMPlex\n";
31 
32 #include <stdbool.h>
33 #include <string.h>
34 #include <ceed.h>
35 #include <petsc.h>
36 #include <petscdmplex.h>
37 #include <petscksp.h>
38 #include <petscsys.h>
39 
40 #include "bps.h"
41 #include "include/bpsproblemdata.h"
42 #include "include/petscutils.h"
43 #include "include/petscversion.h"
44 #include "include/matops.h"
45 #include "include/structs.h"
46 #include "include/libceedsetup.h"
47 
48 #if PETSC_VERSION_LT(3,12,0)
49 #ifdef PETSC_HAVE_CUDA
50 #include <petsccuda.h>
51 // Note: With PETSc prior to version 3.12.0, providing the source path to
52 //       include 'cublas_v2.h' will be needed to use 'petsccuda.h'.
53 #endif
54 #endif
55 
56 int main(int argc, char **argv) {
57   PetscInt ierr;
58   MPI_Comm comm;
59   char filename[PETSC_MAX_PATH_LEN],
60        ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self";
61   double my_rt_start, my_rt, rt_min, rt_max;
62   PetscInt degree = 3, q_extra, *l_size, *xl_size, *g_size, dim = 3, fine_level,
63            mesh_elem[3] = {3, 3, 3}, num_comp_u = 1, num_levels = degree, *level_degrees;
64   PetscScalar *r;
65   PetscScalar eps = 1.0;
66   PetscBool test_mode, benchmark_mode, read_mesh, write_solution, simplex;
67   PetscLogStage solve_stage;
68   PetscLogEvent assemble_event;
69   DM  *dm, dm_orig;
70   KSP ksp;
71   PC pc;
72   Mat *mat_O, *mat_pr, mat_coarse;
73   Vec *X, *X_loc, *mult, rhs, rhs_loc;
74   PetscMemType mem_type;
75   UserO *user_O;
76   UserProlongRestr *user_pr;
77   Ceed ceed;
78   CeedData *ceed_data;
79   CeedVector rhs_ceed, target;
80   CeedQFunction qf_error, qf_restrict, qf_prolong;
81   CeedOperator op_error;
82   BPType bp_choice;
83   CoarsenType coarsen;
84 
85   ierr = PetscInitialize(&argc, &argv, NULL, help);
86   if (ierr) return ierr;
87   comm = PETSC_COMM_WORLD;
88 
89   // Parse command line options
90   PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
91   bp_choice = CEED_BP3;
92   ierr = PetscOptionsEnum("-problem",
93                           "CEED benchmark problem to solve", NULL,
94                           bp_types, (PetscEnum)bp_choice, (PetscEnum *)&bp_choice,
95                           NULL); CHKERRQ(ierr);
96   num_comp_u = bp_options[bp_choice].num_comp_u;
97   test_mode = PETSC_FALSE;
98   ierr = PetscOptionsBool("-test",
99                           "Testing mode (do not print unless error is large)",
100                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
101   benchmark_mode = PETSC_FALSE;
102   ierr = PetscOptionsBool("-benchmark",
103                           "Benchmarking mode (prints benchmark statistics)",
104                           NULL, benchmark_mode, &benchmark_mode, NULL);
105   CHKERRQ(ierr);
106   write_solution = PETSC_FALSE;
107   ierr = PetscOptionsBool("-write_solution",
108                           "Write solution for visualization",
109                           NULL, write_solution, &write_solution, NULL);
110   CHKERRQ(ierr);
111   simplex = PETSC_FALSE;
112   ierr = PetscOptionsBool("-simplex",
113                           "Element topology (default:hex)",
114                           NULL, simplex, &simplex, NULL);
115   CHKERRQ(ierr);
116   ierr = PetscOptionsScalar("-eps",
117                             "Epsilon parameter for Kershaw mesh transformation",
118                             NULL, eps, &eps, NULL);
119   if (eps > 1 || eps <= 0) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
120                                      "-eps %g must be (0,1]", (double)PetscRealPart(eps));
121   degree = test_mode ? 3 : 2;
122   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
123                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
124   if (degree < 1) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
125                             "-degree %" PetscInt_FMT " must be at least 1", degree);
126   q_extra = bp_options[bp_choice].q_extra;
127   ierr = PetscOptionsInt("-q_extra", "Number of extra quadrature points",
128                          NULL, q_extra, &q_extra, NULL); CHKERRQ(ierr);
129   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
130                             NULL, ceed_resource, ceed_resource,
131                             sizeof(ceed_resource), NULL); CHKERRQ(ierr);
132   coarsen = COARSEN_UNIFORM;
133   ierr = PetscOptionsEnum("-coarsen",
134                           "Coarsening strategy to use", NULL,
135                           coarsen_types, (PetscEnum)coarsen,
136                           (PetscEnum *)&coarsen, NULL); CHKERRQ(ierr);
137   read_mesh = PETSC_FALSE;
138   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL,
139                             filename, filename, sizeof(filename), &read_mesh);
140   CHKERRQ(ierr);
141   if (!read_mesh) {
142     PetscInt tmp = dim;
143     ierr = PetscOptionsIntArray("-cells","Number of cells per dimension", NULL,
144                                 mesh_elem, &tmp, NULL); CHKERRQ(ierr);
145   }
146   PetscOptionsEnd();
147 
148   // Set up libCEED
149   CeedInit(ceed_resource, &ceed);
150   CeedMemType mem_type_backend;
151   CeedGetPreferredMemType(ceed, &mem_type_backend);
152 
153   // Setup DM
154   if (read_mesh) {
155     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, NULL, PETSC_TRUE,
156                                 &dm_orig);
157     CHKERRQ(ierr);
158   } else {
159     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, simplex, mesh_elem, NULL,
160                                NULL, NULL, PETSC_TRUE, &dm_orig); CHKERRQ(ierr);
161   }
162 
163   VecType vec_type;
164   switch (mem_type_backend) {
165   case CEED_MEM_HOST: vec_type = VECSTANDARD; break;
166   case CEED_MEM_DEVICE: {
167     const char *resolved;
168     CeedGetResource(ceed, &resolved);
169     if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA;
170     else if (strstr(resolved, "/gpu/hip/occa"))
171       vec_type = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678
172     else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP;
173     else vec_type = VECSTANDARD;
174   }
175   }
176   ierr = DMSetVecType(dm_orig, vec_type); CHKERRQ(ierr);
177   ierr = DMSetFromOptions(dm_orig); CHKERRQ(ierr);
178   ierr = DMViewFromOptions(dm_orig, NULL, "-dm_view"); CHKERRQ(ierr);
179 
180   // Apply Kershaw mesh transformation
181   ierr = Kershaw(dm_orig, eps); CHKERRQ(ierr);
182 
183   // Allocate arrays for PETSc objects for each level
184   switch (coarsen) {
185   case COARSEN_UNIFORM:
186     num_levels = degree;
187     break;
188   case COARSEN_LOGARITHMIC:
189     num_levels = ceil(log(degree)/log(2)) + 1;
190     break;
191   }
192   ierr = PetscMalloc1(num_levels, &level_degrees); CHKERRQ(ierr);
193   fine_level = num_levels - 1;
194 
195   switch (coarsen) {
196   case COARSEN_UNIFORM:
197     for (int i=0; i<num_levels; i++) level_degrees[i] = i + 1;
198     break;
199   case COARSEN_LOGARITHMIC:
200     for (int i=0; i<num_levels - 1; i++) level_degrees[i] = pow(2,i);
201     level_degrees[fine_level] = degree;
202     break;
203   }
204   ierr = PetscMalloc1(num_levels, &dm); CHKERRQ(ierr);
205   ierr = PetscMalloc1(num_levels, &X); CHKERRQ(ierr);
206   ierr = PetscMalloc1(num_levels, &X_loc); CHKERRQ(ierr);
207   ierr = PetscMalloc1(num_levels, &mult); CHKERRQ(ierr);
208   ierr = PetscMalloc1(num_levels, &user_O); CHKERRQ(ierr);
209   ierr = PetscMalloc1(num_levels, &user_pr); CHKERRQ(ierr);
210   ierr = PetscMalloc1(num_levels, &mat_O); CHKERRQ(ierr);
211   ierr = PetscMalloc1(num_levels, &mat_pr); CHKERRQ(ierr);
212   ierr = PetscMalloc1(num_levels, &l_size); CHKERRQ(ierr);
213   ierr = PetscMalloc1(num_levels, &xl_size); CHKERRQ(ierr);
214   ierr = PetscMalloc1(num_levels, &g_size); CHKERRQ(ierr);
215 
216   // Setup DM and Operator Mat Shells for each level
217   for (CeedInt i=0; i<num_levels; i++) {
218     // Create DM
219     ierr = DMClone(dm_orig, &dm[i]); CHKERRQ(ierr);
220     ierr = DMGetVecType(dm_orig, &vec_type); CHKERRQ(ierr);
221     ierr = DMSetVecType(dm[i], vec_type); CHKERRQ(ierr);
222     PetscInt dim;
223     ierr = DMGetDimension(dm[i], &dim); CHKERRQ(ierr);
224     ierr = SetupDMByDegree(dm[i], level_degrees[i], q_extra, num_comp_u, dim,
225                            bp_options[bp_choice].enforce_bc, bp_options[bp_choice].bc_func);
226     CHKERRQ(ierr);
227 
228     // Create vectors
229     ierr = DMCreateGlobalVector(dm[i], &X[i]); CHKERRQ(ierr);
230     ierr = VecGetLocalSize(X[i], &l_size[i]); CHKERRQ(ierr);
231     ierr = VecGetSize(X[i], &g_size[i]); CHKERRQ(ierr);
232     ierr = DMCreateLocalVector(dm[i], &X_loc[i]); CHKERRQ(ierr);
233     ierr = VecGetSize(X_loc[i], &xl_size[i]); CHKERRQ(ierr);
234 
235     // Operator
236     ierr = PetscMalloc1(1, &user_O[i]); CHKERRQ(ierr);
237     ierr = MatCreateShell(comm, l_size[i], l_size[i], g_size[i], g_size[i],
238                           user_O[i], &mat_O[i]); CHKERRQ(ierr);
239     ierr = MatShellSetOperation(mat_O[i], MATOP_MULT,
240                                 (void(*)(void))MatMult_Ceed); CHKERRQ(ierr);
241     ierr = MatShellSetOperation(mat_O[i], MATOP_GET_DIAGONAL,
242                                 (void(*)(void))MatGetDiag); CHKERRQ(ierr);
243     ierr = MatShellSetVecType(mat_O[i], vec_type); CHKERRQ(ierr);
244 
245     // Level transfers
246     if (i > 0) {
247       // Interp
248       ierr = PetscMalloc1(1, &user_pr[i]); CHKERRQ(ierr);
249       ierr = MatCreateShell(comm, l_size[i], l_size[i-1], g_size[i], g_size[i-1],
250                             user_pr[i], &mat_pr[i]); CHKERRQ(ierr);
251       ierr = MatShellSetOperation(mat_pr[i], MATOP_MULT,
252                                   (void(*)(void))MatMult_Prolong);
253       CHKERRQ(ierr);
254       ierr = MatShellSetOperation(mat_pr[i], MATOP_MULT_TRANSPOSE,
255                                   (void(*)(void))MatMult_Restrict);
256       CHKERRQ(ierr);
257       ierr = MatShellSetVecType(mat_pr[i], vec_type); CHKERRQ(ierr);
258     }
259   }
260   ierr = VecDuplicate(X[fine_level], &rhs); CHKERRQ(ierr);
261 
262   // Print global grid information
263   if (!test_mode) {
264     PetscInt P = degree + 1, Q = P + q_extra;
265 
266     const char *used_resource;
267     CeedGetResource(ceed, &used_resource);
268 
269     ierr = VecGetType(X[0], &vec_type); CHKERRQ(ierr);
270 
271     ierr = PetscPrintf(comm,
272                        "\n-- CEED Benchmark Problem %" CeedInt_FMT " -- libCEED + PETSc + PCMG --\n"
273                        "  PETSc:\n"
274                        "    PETSc Vec Type                          : %s\n"
275                        "  libCEED:\n"
276                        "    libCEED Backend                         : %s\n"
277                        "    libCEED Backend MemType                 : %s\n"
278                        "  Mesh:\n"
279                        "    Number of 1D Basis Nodes (p)            : %" CeedInt_FMT "\n"
280                        "    Number of 1D Quadrature Points (q)      : %" CeedInt_FMT "\n"
281                        "    Additional quadrature points (q_extra)  : %" CeedInt_FMT "\n"
282                        "    Global Nodes                            : %" PetscInt_FMT "\n"
283                        "    Owned Nodes                             : %" PetscInt_FMT "\n"
284                        "    DoF per node                            : %" PetscInt_FMT "\n"
285                        "  Multigrid:\n"
286                        "    Number of Levels                        : %" CeedInt_FMT "\n",
287                        bp_choice+1, vec_type, used_resource,
288                        CeedMemTypes[mem_type_backend], P, Q, q_extra,
289                        g_size[fine_level]/num_comp_u, l_size[fine_level]/num_comp_u,
290                        num_comp_u, num_levels); CHKERRQ(ierr);
291   }
292 
293   // Create RHS vector
294   ierr = VecDuplicate(X_loc[fine_level], &rhs_loc); CHKERRQ(ierr);
295   ierr = VecZeroEntries(rhs_loc); CHKERRQ(ierr);
296   ierr = VecGetArrayAndMemType(rhs_loc, &r, &mem_type); CHKERRQ(ierr);
297   CeedVectorCreate(ceed, xl_size[fine_level], &rhs_ceed);
298   CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r);
299 
300   // Set up libCEED operators on each level
301   ierr = PetscMalloc1(num_levels, &ceed_data); CHKERRQ(ierr);
302   for (PetscInt i=0; i<num_levels; i++) {
303     // Print level information
304     if (!test_mode && (i == 0 || i == fine_level)) {
305       ierr = PetscPrintf(comm,"    Level %" PetscInt_FMT " (%s):\n"
306                          "      Number of 1D Basis Nodes (p)     : %" CeedInt_FMT "\n"
307                          "      Global Nodes                     : %" PetscInt_FMT "\n"
308                          "      Owned Nodes                      : %" PetscInt_FMT "\n",
309                          i, (i? "fine" : "coarse"), level_degrees[i] + 1,
310                          g_size[i]/num_comp_u, l_size[i]/num_comp_u); CHKERRQ(ierr);
311     }
312     ierr = PetscMalloc1(1, &ceed_data[i]); CHKERRQ(ierr);
313     ierr = SetupLibceedByDegree(dm[i], ceed, level_degrees[i], dim, q_extra,
314                                 dim, num_comp_u, g_size[i], xl_size[i], bp_options[bp_choice],
315                                 ceed_data[i], i==(fine_level), rhs_ceed, &target);
316     CHKERRQ(ierr);
317   }
318 
319   // Gather RHS
320   CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL);
321   ierr = VecRestoreArrayAndMemType(rhs_loc, &r); CHKERRQ(ierr);
322   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
323   ierr = DMLocalToGlobal(dm[fine_level], rhs_loc, ADD_VALUES, rhs); CHKERRQ(ierr);
324   CeedVectorDestroy(&rhs_ceed);
325 
326   // Create the restriction/interpolation QFunction
327   CeedQFunctionCreateIdentity(ceed, num_comp_u, CEED_EVAL_NONE, CEED_EVAL_INTERP,
328                               &qf_restrict);
329   CeedQFunctionCreateIdentity(ceed, num_comp_u, CEED_EVAL_INTERP, CEED_EVAL_NONE,
330                               &qf_prolong);
331 
332   // Set up libCEED level transfer operators
333   ierr = CeedLevelTransferSetup(ceed, num_levels, num_comp_u, ceed_data,
334                                 level_degrees,
335                                 qf_restrict, qf_prolong); CHKERRQ(ierr);
336 
337   // Create the error QFunction
338   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].error,
339                               bp_options[bp_choice].error_loc, &qf_error);
340   CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP);
341   CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE);
342   CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_NONE);
343 
344   // Create the error operator
345   CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
346                      &op_error);
347   CeedOperatorSetField(op_error, "u", ceed_data[fine_level]->elem_restr_u,
348                        ceed_data[fine_level]->basis_u, CEED_VECTOR_ACTIVE);
349   CeedOperatorSetField(op_error, "true_soln",
350                        ceed_data[fine_level]->elem_restr_u_i,
351                        CEED_BASIS_COLLOCATED, target);
352   CeedOperatorSetField(op_error, "error", ceed_data[fine_level]->elem_restr_u_i,
353                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
354 
355   // Calculate multiplicity
356   for (int i=0; i<num_levels; i++) {
357     PetscScalar *x;
358 
359     // CEED vector
360     ierr = VecZeroEntries(X_loc[i]); CHKERRQ(ierr);
361     ierr = VecGetArray(X_loc[i], &x); CHKERRQ(ierr);
362     CeedVectorSetArray(ceed_data[i]->x_ceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
363 
364     // Multiplicity
365     CeedElemRestrictionGetMultiplicity(ceed_data[i]->elem_restr_u,
366                                        ceed_data[i]->x_ceed);
367     CeedVectorSyncArray(ceed_data[i]->x_ceed, CEED_MEM_HOST);
368 
369     // Restore vector
370     ierr = VecRestoreArray(X_loc[i], &x); CHKERRQ(ierr);
371 
372     // Creat mult vector
373     ierr = VecDuplicate(X_loc[i], &mult[i]); CHKERRQ(ierr);
374 
375     // Local-to-global
376     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
377     ierr = DMLocalToGlobal(dm[i], X_loc[i], ADD_VALUES, X[i]);
378     CHKERRQ(ierr);
379     ierr = VecZeroEntries(X_loc[i]); CHKERRQ(ierr);
380 
381     // Global-to-local
382     ierr = DMGlobalToLocal(dm[i], X[i], INSERT_VALUES, mult[i]);
383     CHKERRQ(ierr);
384     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
385 
386     // Multiplicity scaling
387     ierr = VecReciprocal(mult[i]);
388   }
389 
390   // Set up Mat
391   for (int i=0; i<num_levels; i++) {
392     // User Operator
393     user_O[i]->comm = comm;
394     user_O[i]->dm = dm[i];
395     user_O[i]->X_loc = X_loc[i];
396     ierr = VecDuplicate(X_loc[i], &user_O[i]->Y_loc); CHKERRQ(ierr);
397     user_O[i]->x_ceed = ceed_data[i]->x_ceed;
398     user_O[i]->y_ceed = ceed_data[i]->y_ceed;
399     user_O[i]->op = ceed_data[i]->op_apply;
400     user_O[i]->ceed = ceed;
401 
402     if (i > 0) {
403       // Prolongation/Restriction Operator
404       user_pr[i]->comm = comm;
405       user_pr[i]->dmf = dm[i];
406       user_pr[i]->dmc = dm[i-1];
407       user_pr[i]->loc_vec_c = X_loc[i-1];
408       user_pr[i]->loc_vec_f = user_O[i]->Y_loc;
409       user_pr[i]->mult_vec = mult[i];
410       user_pr[i]->ceed_vec_c = user_O[i-1]->x_ceed;
411       user_pr[i]->ceed_vec_f = user_O[i]->y_ceed;
412       user_pr[i]->op_prolong = ceed_data[i]->op_prolong;
413       user_pr[i]->op_restrict = ceed_data[i]->op_restrict;
414       user_pr[i]->ceed = ceed;
415     }
416   }
417 
418   // Assemble coarse grid Jacobian for AMG (or other sparse matrix) solve
419   ierr = DMCreateMatrix(dm[0], &mat_coarse); CHKERRQ(ierr);
420 
421   ierr = PetscLogEventRegister("AssembleMatrix", MAT_CLASSID, &assemble_event);
422   CHKERRQ(ierr);
423   {
424     // Assemble matrix analytically
425     PetscCount num_entries;
426     CeedInt *rows, *cols;
427     CeedVector coo_values;
428     CeedOperatorLinearAssembleSymbolic(user_O[0]->op, &num_entries, &rows, &cols);
429     ISLocalToGlobalMapping ltog_row, ltog_col;
430     ierr = MatGetLocalToGlobalMapping(mat_coarse, &ltog_row, &ltog_col);
431     CHKERRQ(ierr);
432     ierr = ISLocalToGlobalMappingApply(ltog_row, num_entries, rows, rows);
433     CHKERRQ(ierr);
434     ierr = ISLocalToGlobalMappingApply(ltog_col, num_entries, cols, cols);
435     CHKERRQ(ierr);
436     ierr = MatSetPreallocationCOO(mat_coarse, num_entries, rows, cols);
437     CHKERRQ(ierr);
438     free(rows);
439     free(cols);
440     CeedVectorCreate(ceed, num_entries, &coo_values);
441     ierr = PetscLogEventBegin(assemble_event, mat_coarse, 0, 0, 0); CHKERRQ(ierr);
442     CeedOperatorLinearAssemble(user_O[0]->op, coo_values);
443     const CeedScalar *values;
444     CeedVectorGetArrayRead(coo_values, CEED_MEM_HOST, &values);
445     ierr = MatSetValuesCOO(mat_coarse, values, ADD_VALUES); CHKERRQ(ierr);
446     CeedVectorRestoreArrayRead(coo_values, &values);
447     ierr = PetscLogEventEnd(assemble_event, mat_coarse, 0, 0, 0); CHKERRQ(ierr);
448     CeedVectorDestroy(&coo_values);
449   }
450 
451   // Set up KSP
452   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
453   {
454     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
455     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
456     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
457                             PETSC_DEFAULT); CHKERRQ(ierr);
458   }
459   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
460   ierr = KSPSetOperators(ksp, mat_O[fine_level], mat_O[fine_level]);
461   CHKERRQ(ierr);
462 
463   // Set up PCMG
464   ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
465   PCMGCycleType pcmg_cycle_type = PC_MG_CYCLE_V;
466   {
467     ierr = PCSetType(pc, PCMG); CHKERRQ(ierr);
468 
469     // PCMG levels
470     ierr = PCMGSetLevels(pc, num_levels, NULL); CHKERRQ(ierr);
471     for (int i=0; i<num_levels; i++) {
472       // Smoother
473       KSP smoother;
474       PC smoother_pc;
475       ierr = PCMGGetSmoother(pc, i, &smoother); CHKERRQ(ierr);
476       ierr = KSPSetType(smoother, KSPCHEBYSHEV); CHKERRQ(ierr);
477       ierr = KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1); CHKERRQ(ierr);
478       ierr = KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE); CHKERRQ(ierr);
479       ierr = KSPSetOperators(smoother, mat_O[i], mat_O[i]); CHKERRQ(ierr);
480       ierr = KSPGetPC(smoother, &smoother_pc); CHKERRQ(ierr);
481       ierr = PCSetType(smoother_pc, PCJACOBI); CHKERRQ(ierr);
482       ierr = PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
483 
484       // Work vector
485       if (i < num_levels - 1) {
486         ierr = PCMGSetX(pc, i, X[i]); CHKERRQ(ierr);
487       }
488 
489       // Level transfers
490       if (i > 0) {
491         // Interpolation
492         ierr = PCMGSetInterpolation(pc, i, mat_pr[i]); CHKERRQ(ierr);
493       }
494 
495       // Coarse solve
496       KSP coarse;
497       PC coarse_pc;
498       ierr = PCMGGetCoarseSolve(pc, &coarse); CHKERRQ(ierr);
499       ierr = KSPSetType(coarse, KSPPREONLY); CHKERRQ(ierr);
500       ierr = KSPSetOperators(coarse, mat_coarse, mat_coarse); CHKERRQ(ierr);
501 
502       ierr = KSPGetPC(coarse, &coarse_pc); CHKERRQ(ierr);
503       ierr = PCSetType(coarse_pc, PCGAMG); CHKERRQ(ierr);
504 
505       ierr = KSPSetOptionsPrefix(coarse, "coarse_"); CHKERRQ(ierr);
506       ierr = PCSetOptionsPrefix(coarse_pc, "coarse_"); CHKERRQ(ierr);
507       ierr = KSPSetFromOptions(coarse); CHKERRQ(ierr);
508       ierr = PCSetFromOptions(coarse_pc); CHKERRQ(ierr);
509     }
510 
511     // PCMG options
512     ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr);
513     ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr);
514     ierr = PCMGSetCycleType(pc, pcmg_cycle_type); CHKERRQ(ierr);
515   }
516 
517   // First run, if benchmarking
518   if (benchmark_mode) {
519     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
520     CHKERRQ(ierr);
521     ierr = VecZeroEntries(X[fine_level]); CHKERRQ(ierr);
522     my_rt_start = MPI_Wtime();
523     ierr = KSPSolve(ksp, rhs, X[fine_level]); CHKERRQ(ierr);
524     my_rt = MPI_Wtime() - my_rt_start;
525     ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
526     CHKERRQ(ierr);
527     // Set maxits based on first iteration timing
528     if (my_rt > 0.02) {
529       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
530       CHKERRQ(ierr);
531     } else {
532       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
533       CHKERRQ(ierr);
534     }
535   }
536 
537   // Timed solve
538   ierr = VecZeroEntries(X[fine_level]); CHKERRQ(ierr);
539   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
540 
541   // -- Performance logging
542   ierr = PetscLogStageRegister("Solve Stage", &solve_stage); CHKERRQ(ierr);
543   ierr = PetscLogStagePush(solve_stage); CHKERRQ(ierr);
544 
545   // -- Solve
546   my_rt_start = MPI_Wtime();
547   ierr = KSPSolve(ksp, rhs, X[fine_level]); CHKERRQ(ierr);
548   my_rt = MPI_Wtime() - my_rt_start;
549 
550 
551   // -- Performance logging
552   ierr = PetscLogStagePop();
553 
554   // Output results
555   {
556     KSPType ksp_type;
557     PCMGType pcmg_type;
558     KSPConvergedReason reason;
559     PetscReal rnorm;
560     PetscInt its;
561     ierr = KSPGetType(ksp, &ksp_type); CHKERRQ(ierr);
562     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
563     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
564     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
565     ierr = PCMGGetType(pc, &pcmg_type); CHKERRQ(ierr);
566     if (!test_mode || reason < 0 || rnorm > 1e-8) {
567       ierr = PetscPrintf(comm,
568                          "  KSP:\n"
569                          "    KSP Type                                : %s\n"
570                          "    KSP Convergence                         : %s\n"
571                          "    Total KSP Iterations                    : %D\n"
572                          "    Final rnorm                             : %e\n",
573                          ksp_type, KSPConvergedReasons[reason], its,
574                          (double)rnorm); CHKERRQ(ierr);
575       ierr = PetscPrintf(comm,
576                          "  PCMG:\n"
577                          "    PCMG Type                               : %s\n"
578                          "    PCMG Cycle Type                         : %s\n",
579                          PCMGTypes[pcmg_type],
580                          PCMGCycleTypes[pcmg_cycle_type]); CHKERRQ(ierr);
581     }
582     if (!test_mode) {
583       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
584     }
585     {
586       PetscReal max_error;
587       ierr = ComputeErrorMax(user_O[fine_level], op_error, X[fine_level], target,
588                              &max_error); CHKERRQ(ierr);
589       PetscReal tol = 5e-2;
590       if (!test_mode || max_error > tol) {
591         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
592         CHKERRQ(ierr);
593         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
594         CHKERRQ(ierr);
595         ierr = PetscPrintf(comm,
596                            "    Pointwise Error (max)                   : %e\n"
597                            "    CG Solve Time                           : %g (%g) sec\n",
598                            (double)max_error, rt_max, rt_min); CHKERRQ(ierr);
599       }
600     }
601     if (benchmark_mode && (!test_mode)) {
602       ierr = PetscPrintf(comm,
603                          "    DoFs/Sec in CG                     : %g (%g) million\n",
604                          1e-6*g_size[fine_level]*its/rt_max,
605                          1e-6*g_size[fine_level]*its/rt_min);
606       CHKERRQ(ierr);
607     }
608   }
609 
610   if (write_solution) {
611     PetscViewer vtk_viewer_soln;
612 
613     ierr = PetscViewerCreate(comm, &vtk_viewer_soln); CHKERRQ(ierr);
614     ierr = PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK); CHKERRQ(ierr);
615     ierr = PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"); CHKERRQ(ierr);
616     ierr = VecView(X[fine_level], vtk_viewer_soln); CHKERRQ(ierr);
617     ierr = PetscViewerDestroy(&vtk_viewer_soln); CHKERRQ(ierr);
618   }
619 
620   // Cleanup
621   for (int i=0; i<num_levels; i++) {
622     ierr = VecDestroy(&X[i]); CHKERRQ(ierr);
623     ierr = VecDestroy(&X_loc[i]); CHKERRQ(ierr);
624     ierr = VecDestroy(&mult[i]); CHKERRQ(ierr);
625     ierr = VecDestroy(&user_O[i]->Y_loc); CHKERRQ(ierr);
626     ierr = MatDestroy(&mat_O[i]); CHKERRQ(ierr);
627     ierr = PetscFree(user_O[i]); CHKERRQ(ierr);
628     if (i > 0) {
629       ierr = MatDestroy(&mat_pr[i]); CHKERRQ(ierr);
630       ierr = PetscFree(user_pr[i]); CHKERRQ(ierr);
631     }
632     ierr = CeedDataDestroy(i, ceed_data[i]); CHKERRQ(ierr);
633     ierr = DMDestroy(&dm[i]); CHKERRQ(ierr);
634   }
635   ierr = PetscFree(level_degrees); CHKERRQ(ierr);
636   ierr = PetscFree(dm); CHKERRQ(ierr);
637   ierr = PetscFree(X); CHKERRQ(ierr);
638   ierr = PetscFree(X_loc); CHKERRQ(ierr);
639   ierr = PetscFree(mult); CHKERRQ(ierr);
640   ierr = PetscFree(mat_O); CHKERRQ(ierr);
641   ierr = PetscFree(mat_pr); CHKERRQ(ierr);
642   ierr = PetscFree(ceed_data); CHKERRQ(ierr);
643   ierr = PetscFree(user_O); CHKERRQ(ierr);
644   ierr = PetscFree(user_pr); CHKERRQ(ierr);
645   ierr = PetscFree(l_size); CHKERRQ(ierr);
646   ierr = PetscFree(xl_size); CHKERRQ(ierr);
647   ierr = PetscFree(g_size); CHKERRQ(ierr);
648   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
649   ierr = VecDestroy(&rhs_loc); CHKERRQ(ierr);
650   ierr = MatDestroy(&mat_coarse); CHKERRQ(ierr);
651   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
652   ierr = DMDestroy(&dm_orig); CHKERRQ(ierr);
653   CeedVectorDestroy(&target);
654   CeedQFunctionDestroy(&qf_error);
655   CeedQFunctionDestroy(&qf_restrict);
656   CeedQFunctionDestroy(&qf_prolong);
657   CeedOperatorDestroy(&op_error);
658   CeedDestroy(&ceed);
659   return PetscFinalize();
660 }
661