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