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