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