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