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 PetscLogEvent assemble_event; 78 DM *dm, dm_orig; 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) SETERRQ(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) SETERRQ(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 VecType vec_type; 168 switch (mem_type_backend) { 169 case CEED_MEM_HOST: vec_type = VECSTANDARD; break; 170 case CEED_MEM_DEVICE: { 171 const char *resolved; 172 CeedGetResource(ceed, &resolved); 173 if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA; 174 else if (strstr(resolved, "/gpu/hip/occa")) 175 vec_type = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678 176 else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP; 177 else vec_type = VECSTANDARD; 178 } 179 } 180 ierr = DMSetVecType(dm_orig, vec_type); CHKERRQ(ierr); 181 ierr = DMSetFromOptions(dm_orig); CHKERRQ(ierr); 182 ierr = DMViewFromOptions(dm_orig, NULL, "-dm_view"); CHKERRQ(ierr); 183 184 // Apply Kershaw mesh transformation 185 ierr = Kershaw(dm_orig, eps); CHKERRQ(ierr); 186 187 // Allocate arrays for PETSc objects for each level 188 switch (coarsen) { 189 case COARSEN_UNIFORM: 190 num_levels = degree; 191 break; 192 case COARSEN_LOGARITHMIC: 193 num_levels = ceil(log(degree)/log(2)) + 1; 194 break; 195 } 196 ierr = PetscMalloc1(num_levels, &level_degrees); CHKERRQ(ierr); 197 fine_level = num_levels - 1; 198 199 switch (coarsen) { 200 case COARSEN_UNIFORM: 201 for (int i=0; i<num_levels; i++) level_degrees[i] = i + 1; 202 break; 203 case COARSEN_LOGARITHMIC: 204 for (int i=0; i<num_levels - 1; i++) level_degrees[i] = pow(2,i); 205 level_degrees[fine_level] = degree; 206 break; 207 } 208 ierr = PetscMalloc1(num_levels, &dm); CHKERRQ(ierr); 209 ierr = PetscMalloc1(num_levels, &X); CHKERRQ(ierr); 210 ierr = PetscMalloc1(num_levels, &X_loc); CHKERRQ(ierr); 211 ierr = PetscMalloc1(num_levels, &mult); CHKERRQ(ierr); 212 ierr = PetscMalloc1(num_levels, &user_O); CHKERRQ(ierr); 213 ierr = PetscMalloc1(num_levels, &user_pr); CHKERRQ(ierr); 214 ierr = PetscMalloc1(num_levels, &mat_O); CHKERRQ(ierr); 215 ierr = PetscMalloc1(num_levels, &mat_pr); CHKERRQ(ierr); 216 ierr = PetscMalloc1(num_levels, &l_size); CHKERRQ(ierr); 217 ierr = PetscMalloc1(num_levels, &xl_size); CHKERRQ(ierr); 218 ierr = PetscMalloc1(num_levels, &g_size); CHKERRQ(ierr); 219 220 // Setup DM and Operator Mat Shells for each level 221 for (CeedInt i=0; i<num_levels; i++) { 222 // Create DM 223 ierr = DMClone(dm_orig, &dm[i]); CHKERRQ(ierr); 224 ierr = DMGetVecType(dm_orig, &vec_type); CHKERRQ(ierr); 225 ierr = DMSetVecType(dm[i], vec_type); CHKERRQ(ierr); 226 PetscInt dim; 227 ierr = DMGetDimension(dm[i], &dim); CHKERRQ(ierr); 228 ierr = SetupDMByDegree(dm[i], level_degrees[i], num_comp_u, dim, 229 bp_options[bp_choice].enforce_bc, bp_options[bp_choice].bc_func); 230 CHKERRQ(ierr); 231 232 // Create vectors 233 ierr = DMCreateGlobalVector(dm[i], &X[i]); CHKERRQ(ierr); 234 ierr = VecGetLocalSize(X[i], &l_size[i]); CHKERRQ(ierr); 235 ierr = VecGetSize(X[i], &g_size[i]); CHKERRQ(ierr); 236 ierr = DMCreateLocalVector(dm[i], &X_loc[i]); CHKERRQ(ierr); 237 ierr = VecGetSize(X_loc[i], &xl_size[i]); CHKERRQ(ierr); 238 239 // Operator 240 ierr = PetscMalloc1(1, &user_O[i]); CHKERRQ(ierr); 241 ierr = MatCreateShell(comm, l_size[i], l_size[i], g_size[i], g_size[i], 242 user_O[i], &mat_O[i]); CHKERRQ(ierr); 243 ierr = MatShellSetOperation(mat_O[i], MATOP_MULT, 244 (void(*)(void))MatMult_Ceed); CHKERRQ(ierr); 245 ierr = MatShellSetOperation(mat_O[i], MATOP_GET_DIAGONAL, 246 (void(*)(void))MatGetDiag); CHKERRQ(ierr); 247 ierr = MatShellSetVecType(mat_O[i], vec_type); CHKERRQ(ierr); 248 249 // Level transfers 250 if (i > 0) { 251 // Interp 252 ierr = PetscMalloc1(1, &user_pr[i]); CHKERRQ(ierr); 253 ierr = MatCreateShell(comm, l_size[i], l_size[i-1], g_size[i], g_size[i-1], 254 user_pr[i], &mat_pr[i]); CHKERRQ(ierr); 255 ierr = MatShellSetOperation(mat_pr[i], MATOP_MULT, 256 (void(*)(void))MatMult_Prolong); 257 CHKERRQ(ierr); 258 ierr = MatShellSetOperation(mat_pr[i], MATOP_MULT_TRANSPOSE, 259 (void(*)(void))MatMult_Restrict); 260 CHKERRQ(ierr); 261 ierr = MatShellSetVecType(mat_pr[i], vec_type); CHKERRQ(ierr); 262 } 263 } 264 ierr = VecDuplicate(X[fine_level], &rhs); CHKERRQ(ierr); 265 266 // Print global grid information 267 if (!test_mode) { 268 PetscInt P = degree + 1, Q = P + q_extra; 269 270 const char *used_resource; 271 CeedGetResource(ceed, &used_resource); 272 273 ierr = VecGetType(X[0], &vec_type); CHKERRQ(ierr); 274 275 ierr = PetscPrintf(comm, 276 "\n-- CEED Benchmark Problem %d -- libCEED + PETSc + PCMG --\n" 277 " PETSc:\n" 278 " PETSc Vec Type : %s\n" 279 " libCEED:\n" 280 " libCEED Backend : %s\n" 281 " libCEED Backend MemType : %s\n" 282 " Mesh:\n" 283 " Number of 1D Basis Nodes (p) : %d\n" 284 " Number of 1D Quadrature Points (q) : %d\n" 285 " Global Nodes : %D\n" 286 " Owned Nodes : %D\n" 287 " DoF per node : %D\n" 288 " Multigrid:\n" 289 " Number of Levels : %d\n", 290 bp_choice+1, vec_type, used_resource, 291 CeedMemTypes[mem_type_backend], 292 P, Q, g_size[fine_level]/num_comp_u, l_size[fine_level]/num_comp_u, 293 num_comp_u, num_levels); CHKERRQ(ierr); 294 } 295 296 // Create RHS vector 297 ierr = VecDuplicate(X_loc[fine_level], &rhs_loc); CHKERRQ(ierr); 298 ierr = VecZeroEntries(rhs_loc); CHKERRQ(ierr); 299 ierr = VecGetArrayAndMemType(rhs_loc, &r, &mem_type); CHKERRQ(ierr); 300 CeedVectorCreate(ceed, xl_size[fine_level], &rhs_ceed); 301 CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r); 302 303 // Set up libCEED operators on each level 304 ierr = PetscMalloc1(num_levels, &ceed_data); CHKERRQ(ierr); 305 for (int i=0; i<num_levels; i++) { 306 // Print level information 307 if (!test_mode && (i == 0 || i == fine_level)) { 308 ierr = PetscPrintf(comm," Level %D (%s):\n" 309 " Number of 1D Basis Nodes (p) : %d\n" 310 " Global Nodes : %D\n" 311 " Owned Nodes : %D\n", 312 i, (i? "fine" : "coarse"), level_degrees[i] + 1, 313 g_size[i]/num_comp_u, l_size[i]/num_comp_u); CHKERRQ(ierr); 314 } 315 ierr = PetscMalloc1(1, &ceed_data[i]); CHKERRQ(ierr); 316 ierr = SetupLibceedByDegree(dm[i], ceed, level_degrees[i], dim, q_extra, 317 dim, num_comp_u, g_size[i], xl_size[i], bp_options[bp_choice], 318 ceed_data[i], i==(fine_level), rhs_ceed, &target); 319 CHKERRQ(ierr); 320 } 321 322 // Gather RHS 323 CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL); 324 ierr = VecRestoreArrayAndMemType(rhs_loc, &r); CHKERRQ(ierr); 325 ierr = VecZeroEntries(rhs); CHKERRQ(ierr); 326 ierr = DMLocalToGlobal(dm[fine_level], rhs_loc, ADD_VALUES, rhs); CHKERRQ(ierr); 327 CeedVectorDestroy(&rhs_ceed); 328 329 // Create the restriction/interpolation QFunction 330 CeedQFunctionCreateIdentity(ceed, num_comp_u, CEED_EVAL_NONE, CEED_EVAL_INTERP, 331 &qf_restrict); 332 CeedQFunctionCreateIdentity(ceed, num_comp_u, CEED_EVAL_INTERP, CEED_EVAL_NONE, 333 &qf_prolong); 334 335 // Set up libCEED level transfer operators 336 ierr = CeedLevelTransferSetup(ceed, num_levels, num_comp_u, ceed_data, 337 level_degrees, 338 qf_restrict, qf_prolong); CHKERRQ(ierr); 339 340 // Create the error QFunction 341 CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].error, 342 bp_options[bp_choice].error_loc, &qf_error); 343 CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP); 344 CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE); 345 CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_NONE); 346 347 // Create the error operator 348 CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, 349 &op_error); 350 CeedOperatorSetField(op_error, "u", ceed_data[fine_level]->elem_restr_u, 351 ceed_data[fine_level]->basis_u, CEED_VECTOR_ACTIVE); 352 CeedOperatorSetField(op_error, "true_soln", 353 ceed_data[fine_level]->elem_restr_u_i, 354 CEED_BASIS_COLLOCATED, target); 355 CeedOperatorSetField(op_error, "error", ceed_data[fine_level]->elem_restr_u_i, 356 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 357 358 // Calculate multiplicity 359 for (int i=0; i<num_levels; i++) { 360 PetscScalar *x; 361 362 // CEED vector 363 ierr = VecZeroEntries(X_loc[i]); CHKERRQ(ierr); 364 ierr = VecGetArray(X_loc[i], &x); CHKERRQ(ierr); 365 CeedVectorSetArray(ceed_data[i]->x_ceed, CEED_MEM_HOST, CEED_USE_POINTER, x); 366 367 // Multiplicity 368 CeedElemRestrictionGetMultiplicity(ceed_data[i]->elem_restr_u, 369 ceed_data[i]->x_ceed); 370 CeedVectorSyncArray(ceed_data[i]->x_ceed, CEED_MEM_HOST); 371 372 // Restore vector 373 ierr = VecRestoreArray(X_loc[i], &x); CHKERRQ(ierr); 374 375 // Creat mult vector 376 ierr = VecDuplicate(X_loc[i], &mult[i]); CHKERRQ(ierr); 377 378 // Local-to-global 379 ierr = VecZeroEntries(X[i]); CHKERRQ(ierr); 380 ierr = DMLocalToGlobal(dm[i], X_loc[i], ADD_VALUES, X[i]); 381 CHKERRQ(ierr); 382 ierr = VecZeroEntries(X_loc[i]); CHKERRQ(ierr); 383 384 // Global-to-local 385 ierr = DMGlobalToLocal(dm[i], X[i], INSERT_VALUES, mult[i]); 386 CHKERRQ(ierr); 387 ierr = VecZeroEntries(X[i]); CHKERRQ(ierr); 388 389 // Multiplicity scaling 390 ierr = VecReciprocal(mult[i]); 391 } 392 393 // Set up Mat 394 for (int i=0; i<num_levels; i++) { 395 // User Operator 396 user_O[i]->comm = comm; 397 user_O[i]->dm = dm[i]; 398 user_O[i]->X_loc = X_loc[i]; 399 ierr = VecDuplicate(X_loc[i], &user_O[i]->Y_loc); CHKERRQ(ierr); 400 user_O[i]->x_ceed = ceed_data[i]->x_ceed; 401 user_O[i]->y_ceed = ceed_data[i]->y_ceed; 402 user_O[i]->op = ceed_data[i]->op_apply; 403 user_O[i]->ceed = ceed; 404 405 if (i > 0) { 406 // Prolongation/Restriction Operator 407 user_pr[i]->comm = comm; 408 user_pr[i]->dmf = dm[i]; 409 user_pr[i]->dmc = dm[i-1]; 410 user_pr[i]->loc_vec_c = X_loc[i-1]; 411 user_pr[i]->loc_vec_f = user_O[i]->Y_loc; 412 user_pr[i]->mult_vec = mult[i]; 413 user_pr[i]->ceed_vec_c = user_O[i-1]->x_ceed; 414 user_pr[i]->ceed_vec_f = user_O[i]->y_ceed; 415 user_pr[i]->op_prolong = ceed_data[i]->op_prolong; 416 user_pr[i]->op_restrict = ceed_data[i]->op_restrict; 417 user_pr[i]->ceed = ceed; 418 } 419 } 420 421 // Assemble coarse grid Jacobian for AMG (or other sparse matrix) solve 422 ierr = DMCreateMatrix(dm[0], &mat_coarse); CHKERRQ(ierr); 423 424 ierr = PetscLogEventRegister("AssembleMatrix", MAT_CLASSID, &assemble_event); 425 CHKERRQ(ierr); 426 { 427 // Assemble matrix analytically 428 PetscCount num_entries; 429 CeedInt *rows, *cols; 430 CeedVector coo_values; 431 CeedOperatorLinearAssembleSymbolic(user_O[0]->op, &num_entries, &rows, &cols); 432 ISLocalToGlobalMapping ltog_row, ltog_col; 433 ierr = MatGetLocalToGlobalMapping(mat_coarse, <og_row, <og_col); 434 CHKERRQ(ierr); 435 ierr = ISLocalToGlobalMappingApply(ltog_row, num_entries, rows, rows); 436 CHKERRQ(ierr); 437 ierr = ISLocalToGlobalMappingApply(ltog_col, num_entries, cols, cols); 438 CHKERRQ(ierr); 439 ierr = MatSetPreallocationCOO(mat_coarse, num_entries, rows, cols); 440 CHKERRQ(ierr); 441 free(rows); 442 free(cols); 443 CeedVectorCreate(ceed, num_entries, &coo_values); 444 ierr = PetscLogEventBegin(assemble_event, mat_coarse, 0, 0, 0); CHKERRQ(ierr); 445 CeedOperatorLinearAssemble(user_O[0]->op, coo_values); 446 const CeedScalar *values; 447 CeedVectorGetArrayRead(coo_values, CEED_MEM_HOST, &values); 448 ierr = MatSetValuesCOO(mat_coarse, values, ADD_VALUES); CHKERRQ(ierr); 449 CeedVectorRestoreArrayRead(coo_values, &values); 450 ierr = PetscLogEventEnd(assemble_event, mat_coarse, 0, 0, 0); CHKERRQ(ierr); 451 CeedVectorDestroy(&coo_values); 452 } 453 454 // Set up KSP 455 ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr); 456 { 457 ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 458 ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 459 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 460 PETSC_DEFAULT); CHKERRQ(ierr); 461 } 462 ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr); 463 ierr = KSPSetOperators(ksp, mat_O[fine_level], mat_O[fine_level]); 464 CHKERRQ(ierr); 465 466 // Set up PCMG 467 ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 468 PCMGCycleType pcmg_cycle_type = PC_MG_CYCLE_V; 469 { 470 ierr = PCSetType(pc, PCMG); CHKERRQ(ierr); 471 472 // PCMG levels 473 ierr = PCMGSetLevels(pc, num_levels, NULL); CHKERRQ(ierr); 474 for (int i=0; i<num_levels; i++) { 475 // Smoother 476 KSP smoother; 477 PC smoother_pc; 478 ierr = PCMGGetSmoother(pc, i, &smoother); CHKERRQ(ierr); 479 ierr = KSPSetType(smoother, KSPCHEBYSHEV); CHKERRQ(ierr); 480 ierr = KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1); CHKERRQ(ierr); 481 ierr = KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE); CHKERRQ(ierr); 482 ierr = KSPSetOperators(smoother, mat_O[i], mat_O[i]); CHKERRQ(ierr); 483 ierr = KSPGetPC(smoother, &smoother_pc); CHKERRQ(ierr); 484 ierr = PCSetType(smoother_pc, PCJACOBI); CHKERRQ(ierr); 485 ierr = PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 486 487 // Work vector 488 if (i < num_levels - 1) { 489 ierr = PCMGSetX(pc, i, X[i]); CHKERRQ(ierr); 490 } 491 492 // Level transfers 493 if (i > 0) { 494 // Interpolation 495 ierr = PCMGSetInterpolation(pc, i, mat_pr[i]); CHKERRQ(ierr); 496 } 497 498 // Coarse solve 499 KSP coarse; 500 PC coarse_pc; 501 ierr = PCMGGetCoarseSolve(pc, &coarse); CHKERRQ(ierr); 502 ierr = KSPSetType(coarse, KSPPREONLY); CHKERRQ(ierr); 503 ierr = KSPSetOperators(coarse, mat_coarse, mat_coarse); CHKERRQ(ierr); 504 505 ierr = KSPGetPC(coarse, &coarse_pc); CHKERRQ(ierr); 506 ierr = PCSetType(coarse_pc, PCGAMG); CHKERRQ(ierr); 507 508 ierr = KSPSetOptionsPrefix(coarse, "coarse_"); CHKERRQ(ierr); 509 ierr = PCSetOptionsPrefix(coarse_pc, "coarse_"); CHKERRQ(ierr); 510 ierr = KSPSetFromOptions(coarse); CHKERRQ(ierr); 511 ierr = PCSetFromOptions(coarse_pc); CHKERRQ(ierr); 512 } 513 514 // PCMG options 515 ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr); 516 ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr); 517 ierr = PCMGSetCycleType(pc, pcmg_cycle_type); CHKERRQ(ierr); 518 } 519 520 // First run, if benchmarking 521 if (benchmark_mode) { 522 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1); 523 CHKERRQ(ierr); 524 ierr = VecZeroEntries(X[fine_level]); CHKERRQ(ierr); 525 my_rt_start = MPI_Wtime(); 526 ierr = KSPSolve(ksp, rhs, X[fine_level]); CHKERRQ(ierr); 527 my_rt = MPI_Wtime() - my_rt_start; 528 ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm); 529 CHKERRQ(ierr); 530 // Set maxits based on first iteration timing 531 if (my_rt > 0.02) { 532 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5); 533 CHKERRQ(ierr); 534 } else { 535 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20); 536 CHKERRQ(ierr); 537 } 538 } 539 540 // Timed solve 541 ierr = VecZeroEntries(X[fine_level]); CHKERRQ(ierr); 542 ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr); 543 544 // -- Performance logging 545 ierr = PetscLogStageRegister("Solve Stage", &solve_stage); CHKERRQ(ierr); 546 ierr = PetscLogStagePush(solve_stage); CHKERRQ(ierr); 547 548 // -- Solve 549 my_rt_start = MPI_Wtime(); 550 ierr = KSPSolve(ksp, rhs, X[fine_level]); CHKERRQ(ierr); 551 my_rt = MPI_Wtime() - my_rt_start; 552 553 554 // -- Performance logging 555 ierr = PetscLogStagePop(); 556 557 // Output results 558 { 559 KSPType ksp_type; 560 PCMGType pcmg_type; 561 KSPConvergedReason reason; 562 PetscReal rnorm; 563 PetscInt its; 564 ierr = KSPGetType(ksp, &ksp_type); CHKERRQ(ierr); 565 ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr); 566 ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr); 567 ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr); 568 ierr = PCMGGetType(pc, &pcmg_type); CHKERRQ(ierr); 569 if (!test_mode || reason < 0 || rnorm > 1e-8) { 570 ierr = PetscPrintf(comm, 571 " KSP:\n" 572 " KSP Type : %s\n" 573 " KSP Convergence : %s\n" 574 " Total KSP Iterations : %D\n" 575 " Final rnorm : %e\n", 576 ksp_type, KSPConvergedReasons[reason], its, 577 (double)rnorm); CHKERRQ(ierr); 578 ierr = PetscPrintf(comm, 579 " PCMG:\n" 580 " PCMG Type : %s\n" 581 " PCMG Cycle Type : %s\n", 582 PCMGTypes[pcmg_type], 583 PCMGCycleTypes[pcmg_cycle_type]); CHKERRQ(ierr); 584 } 585 if (!test_mode) { 586 ierr = PetscPrintf(comm," Performance:\n"); CHKERRQ(ierr); 587 } 588 { 589 PetscReal max_error; 590 ierr = ComputeErrorMax(user_O[fine_level], op_error, X[fine_level], target, 591 &max_error); CHKERRQ(ierr); 592 PetscReal tol = 5e-2; 593 if (!test_mode || max_error > tol) { 594 ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm); 595 CHKERRQ(ierr); 596 ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm); 597 CHKERRQ(ierr); 598 ierr = PetscPrintf(comm, 599 " Pointwise Error (max) : %e\n" 600 " CG Solve Time : %g (%g) sec\n", 601 (double)max_error, rt_max, rt_min); CHKERRQ(ierr); 602 } 603 } 604 if (benchmark_mode && (!test_mode)) { 605 ierr = PetscPrintf(comm, 606 " DoFs/Sec in CG : %g (%g) million\n", 607 1e-6*g_size[fine_level]*its/rt_max, 608 1e-6*g_size[fine_level]*its/rt_min); 609 CHKERRQ(ierr); 610 } 611 } 612 613 if (write_solution) { 614 PetscViewer vtk_viewer_soln; 615 616 ierr = PetscViewerCreate(comm, &vtk_viewer_soln); CHKERRQ(ierr); 617 ierr = PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK); CHKERRQ(ierr); 618 ierr = PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"); CHKERRQ(ierr); 619 ierr = VecView(X[fine_level], vtk_viewer_soln); CHKERRQ(ierr); 620 ierr = PetscViewerDestroy(&vtk_viewer_soln); CHKERRQ(ierr); 621 } 622 623 // Cleanup 624 for (int i=0; i<num_levels; i++) { 625 ierr = VecDestroy(&X[i]); CHKERRQ(ierr); 626 ierr = VecDestroy(&X_loc[i]); CHKERRQ(ierr); 627 ierr = VecDestroy(&mult[i]); CHKERRQ(ierr); 628 ierr = VecDestroy(&user_O[i]->Y_loc); CHKERRQ(ierr); 629 ierr = MatDestroy(&mat_O[i]); CHKERRQ(ierr); 630 ierr = PetscFree(user_O[i]); CHKERRQ(ierr); 631 if (i > 0) { 632 ierr = MatDestroy(&mat_pr[i]); CHKERRQ(ierr); 633 ierr = PetscFree(user_pr[i]); CHKERRQ(ierr); 634 } 635 ierr = CeedDataDestroy(i, ceed_data[i]); CHKERRQ(ierr); 636 ierr = DMDestroy(&dm[i]); CHKERRQ(ierr); 637 } 638 ierr = PetscFree(level_degrees); CHKERRQ(ierr); 639 ierr = PetscFree(dm); CHKERRQ(ierr); 640 ierr = PetscFree(X); CHKERRQ(ierr); 641 ierr = PetscFree(X_loc); CHKERRQ(ierr); 642 ierr = PetscFree(mult); CHKERRQ(ierr); 643 ierr = PetscFree(mat_O); CHKERRQ(ierr); 644 ierr = PetscFree(mat_pr); CHKERRQ(ierr); 645 ierr = PetscFree(ceed_data); CHKERRQ(ierr); 646 ierr = PetscFree(user_O); CHKERRQ(ierr); 647 ierr = PetscFree(user_pr); CHKERRQ(ierr); 648 ierr = PetscFree(l_size); CHKERRQ(ierr); 649 ierr = PetscFree(xl_size); CHKERRQ(ierr); 650 ierr = PetscFree(g_size); CHKERRQ(ierr); 651 ierr = VecDestroy(&rhs); CHKERRQ(ierr); 652 ierr = VecDestroy(&rhs_loc); CHKERRQ(ierr); 653 ierr = MatDestroy(&mat_coarse); CHKERRQ(ierr); 654 ierr = KSPDestroy(&ksp); CHKERRQ(ierr); 655 ierr = DMDestroy(&dm_orig); CHKERRQ(ierr); 656 CeedVectorDestroy(&target); 657 CeedQFunctionDestroy(&qf_error); 658 CeedQFunctionDestroy(&qf_restrict); 659 CeedQFunctionDestroy(&qf_prolong); 660 CeedOperatorDestroy(&op_error); 661 CeedDestroy(&ceed); 662 return PetscFinalize(); 663 } 664