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