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