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