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