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