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