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