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