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