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