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