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 9 // 10 // This example demonstrates a simple usage of libCEED with PETSc to solve the CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps. 11 // 12 // The code uses higher level communication protocols in DMPlex. 13 // 14 // Build with: 15 // 16 // make bps [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] 17 // 18 // Sample runs: 19 // 20 // ./bps -problem bp1 -degree 3 21 // ./bps -problem bp2 -degree 3 22 // ./bps -problem bp3 -degree 3 23 // ./bps -problem bp4 -degree 3 24 // ./bps -problem bp5 -degree 3 -ceed /cpu/self 25 // ./bps -problem bp6 -degree 3 -ceed /gpu/cuda 26 // 27 //TESTARGS -ceed {ceed_resource} -test -problem bp5 -degree 3 -ksp_max_it_clip 15,15 28 //TESTARGS -ceed {ceed_resource} -test -problem bp3 -degree 3 -ksp_max_it_clip 50,50 -simplex 29 30 /// @file 31 /// CEED BPs example using PETSc with DMPlex 32 /// See bpsraw.c for a "raw" implementation using a structured grid. 33 const char help[] = "Solve CEED BPs using PETSc with DMPlex\n"; 34 35 #include "bps.h" 36 37 #include <ceed.h> 38 #include <petscdmplex.h> 39 #include <petscksp.h> 40 #include <stdbool.h> 41 #include <string.h> 42 43 #include "include/bpsproblemdata.h" 44 #include "include/libceedsetup.h" 45 #include "include/matops.h" 46 #include "include/petscutils.h" 47 #include "include/petscversion.h" 48 #include "include/structs.h" 49 50 // ----------------------------------------------------------------------------- 51 // Main body of program, called in a loop for performance benchmarking purposes 52 // ----------------------------------------------------------------------------- 53 static PetscErrorCode RunWithDM(RunParams rp, DM dm, const char *ceed_resource) { 54 double my_rt_start, my_rt, rt_min, rt_max; 55 PetscInt xl_size, l_size, g_size; 56 PetscScalar *r; 57 Vec X, X_loc, rhs, rhs_loc; 58 Mat mat_O; 59 KSP ksp; 60 OperatorApplyContext op_apply_ctx, op_error_ctx; 61 Ceed ceed; 62 CeedData ceed_data; 63 CeedQFunction qf_error; 64 CeedOperator op_error; 65 CeedVector rhs_ceed, target; 66 VecType vec_type; 67 PetscMemType mem_type; 68 69 PetscFunctionBeginUser; 70 // Set up libCEED 71 CeedInit(ceed_resource, &ceed); 72 CeedMemType mem_type_backend; 73 CeedGetPreferredMemType(ceed, &mem_type_backend); 74 75 PetscCall(DMGetVecType(dm, &vec_type)); 76 if (!vec_type) { // Not yet set by user -dm_vec_type 77 switch (mem_type_backend) { 78 case CEED_MEM_HOST: 79 vec_type = VECSTANDARD; 80 break; 81 case CEED_MEM_DEVICE: { 82 const char *resolved; 83 CeedGetResource(ceed, &resolved); 84 if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA; 85 else if (strstr(resolved, "/gpu/hip/occa")) vec_type = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678 86 else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP; 87 else vec_type = VECSTANDARD; 88 } 89 } 90 PetscCall(DMSetVecType(dm, vec_type)); 91 } 92 93 // Create global and local solution vectors 94 PetscCall(DMCreateGlobalVector(dm, &X)); 95 PetscCall(VecGetLocalSize(X, &l_size)); 96 PetscCall(VecGetSize(X, &g_size)); 97 PetscCall(DMCreateLocalVector(dm, &X_loc)); 98 PetscCall(VecGetSize(X_loc, &xl_size)); 99 PetscCall(VecDuplicate(X, &rhs)); 100 101 // Operator 102 PetscCall(PetscMalloc1(1, &op_apply_ctx)); 103 PetscCall(PetscMalloc1(1, &op_error_ctx)); 104 PetscCall(MatCreateShell(rp->comm, l_size, l_size, g_size, g_size, op_apply_ctx, &mat_O)); 105 PetscCall(MatShellSetOperation(mat_O, MATOP_MULT, (void (*)(void))MatMult_Ceed)); 106 PetscCall(MatShellSetOperation(mat_O, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag)); 107 PetscCall(MatShellSetVecType(mat_O, vec_type)); 108 109 // Print summary 110 if (!rp->test_mode) { 111 PetscInt P = rp->degree + 1, Q = P + rp->q_extra; 112 113 const char *used_resource; 114 CeedGetResource(ceed, &used_resource); 115 116 VecType vec_type; 117 PetscCall(VecGetType(X, &vec_type)); 118 119 PetscInt c_start, c_end; 120 PetscCall(DMPlexGetHeightStratum(dm, 0, &c_start, &c_end)); 121 DMPolytopeType cell_type; 122 PetscCall(DMPlexGetCellType(dm, c_start, &cell_type)); 123 CeedElemTopology elem_topo = ElemTopologyP2C(cell_type); 124 PetscMPIInt comm_size; 125 PetscCall(MPI_Comm_size(rp->comm, &comm_size)); 126 PetscCall(PetscPrintf(rp->comm, 127 "\n-- CEED Benchmark Problem %" CeedInt_FMT " -- libCEED + PETSc --\n" 128 " MPI:\n" 129 " Hostname : %s\n" 130 " Total ranks : %d\n" 131 " Ranks per compute node : %d\n" 132 " PETSc:\n" 133 " PETSc Vec Type : %s\n" 134 " libCEED:\n" 135 " libCEED Backend : %s\n" 136 " libCEED Backend MemType : %s\n" 137 " Mesh:\n" 138 " Solution Order (P) : %" CeedInt_FMT "\n" 139 " Quadrature Order (Q) : %" CeedInt_FMT "\n" 140 " Additional quadrature points (q_extra) : %" CeedInt_FMT "\n" 141 " Global nodes : %" PetscInt_FMT "\n" 142 " Local Elements : %" PetscInt_FMT "\n" 143 " Element topology : %s\n" 144 " Owned nodes : %" PetscInt_FMT "\n" 145 " DoF per node : %" PetscInt_FMT "\n", 146 rp->bp_choice + 1, rp->hostname, comm_size, rp->ranks_per_node, vec_type, used_resource, CeedMemTypes[mem_type_backend], P, 147 Q, rp->q_extra, g_size / rp->num_comp_u, c_end - c_start, CeedElemTopologies[elem_topo], l_size / rp->num_comp_u, 148 rp->num_comp_u)); 149 } 150 151 // Create RHS vector 152 PetscCall(VecDuplicate(X_loc, &rhs_loc)); 153 PetscCall(VecZeroEntries(rhs_loc)); 154 PetscCall(VecGetArrayAndMemType(rhs_loc, &r, &mem_type)); 155 CeedVectorCreate(ceed, xl_size, &rhs_ceed); 156 CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r); 157 158 PetscCall(PetscMalloc1(1, &ceed_data)); 159 PetscCall(SetupLibceedByDegree(dm, ceed, rp->degree, rp->dim, rp->q_extra, rp->dim, rp->num_comp_u, g_size, xl_size, bp_options[rp->bp_choice], 160 ceed_data, true, rhs_ceed, &target)); 161 162 // Gather RHS 163 CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL); 164 PetscCall(VecRestoreArrayAndMemType(rhs_loc, &r)); 165 PetscCall(VecZeroEntries(rhs)); 166 PetscCall(DMLocalToGlobal(dm, rhs_loc, ADD_VALUES, rhs)); 167 CeedVectorDestroy(&rhs_ceed); 168 169 // Create the error QFunction 170 CeedQFunctionCreateInterior(ceed, 1, bp_options[rp->bp_choice].error, bp_options[rp->bp_choice].error_loc, &qf_error); 171 CeedQFunctionAddInput(qf_error, "u", rp->num_comp_u, CEED_EVAL_INTERP); 172 CeedQFunctionAddInput(qf_error, "true_soln", rp->num_comp_u, CEED_EVAL_NONE); 173 CeedQFunctionAddInput(qf_error, "qdata", ceed_data->q_data_size, CEED_EVAL_NONE); 174 CeedQFunctionAddOutput(qf_error, "error", rp->num_comp_u, CEED_EVAL_INTERP); 175 176 // Create the error operator 177 CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_error); 178 CeedOperatorSetField(op_error, "u", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE); 179 CeedOperatorSetField(op_error, "true_soln", ceed_data->elem_restr_u_i, CEED_BASIS_COLLOCATED, target); 180 CeedOperatorSetField(op_error, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 181 CeedOperatorSetField(op_error, "error", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE); 182 183 // Set up apply operator context 184 PetscCall(SetupApplyOperatorCtx(rp->comm, dm, ceed, ceed_data, X_loc, op_apply_ctx)); 185 PetscCall(KSPCreate(rp->comm, &ksp)); 186 { 187 PC pc; 188 PetscCall(KSPGetPC(ksp, &pc)); 189 if (rp->bp_choice == CEED_BP1 || rp->bp_choice == CEED_BP2) { 190 PetscCall(PCSetType(pc, PCJACOBI)); 191 if (rp->simplex) { 192 PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL)); 193 } else { 194 PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 195 } 196 } else { 197 PetscCall(PCSetType(pc, PCNONE)); 198 } 199 PetscCall(KSPSetType(ksp, KSPCG)); 200 PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL)); 201 PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 202 } 203 PetscCall(KSPSetOperators(ksp, mat_O, mat_O)); 204 205 // First run's performance log is not considered for benchmarking purposes 206 PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1)); 207 my_rt_start = MPI_Wtime(); 208 PetscCall(KSPSolve(ksp, rhs, X)); 209 my_rt = MPI_Wtime() - my_rt_start; 210 PetscCall(MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, rp->comm)); 211 // Set maxits based on first iteration timing 212 if (my_rt > 0.02) { 213 PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, rp->ksp_max_it_clip[0])); 214 } else { 215 PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, rp->ksp_max_it_clip[1])); 216 } 217 PetscCall(KSPSetFromOptions(ksp)); 218 219 // Timed solve 220 PetscCall(VecZeroEntries(X)); 221 PetscCall(PetscBarrier((PetscObject)ksp)); 222 223 // -- Performance logging 224 PetscCall(PetscLogStagePush(rp->solve_stage)); 225 226 // -- Solve 227 my_rt_start = MPI_Wtime(); 228 PetscCall(KSPSolve(ksp, rhs, X)); 229 my_rt = MPI_Wtime() - my_rt_start; 230 231 // -- Performance logging 232 PetscCall(PetscLogStagePop()); 233 234 // Output results 235 { 236 KSPType ksp_type; 237 KSPConvergedReason reason; 238 PetscReal rnorm; 239 PetscInt its; 240 PetscCall(KSPGetType(ksp, &ksp_type)); 241 PetscCall(KSPGetConvergedReason(ksp, &reason)); 242 PetscCall(KSPGetIterationNumber(ksp, &its)); 243 PetscCall(KSPGetResidualNorm(ksp, &rnorm)); 244 if (!rp->test_mode || reason < 0 || rnorm > 1e-8) { 245 PetscCall(PetscPrintf(rp->comm, 246 " KSP:\n" 247 " KSP Type : %s\n" 248 " KSP Convergence : %s\n" 249 " Total KSP Iterations : %" PetscInt_FMT "\n" 250 " Final rnorm : %e\n", 251 ksp_type, KSPConvergedReasons[reason], its, (double)rnorm)); 252 } 253 if (!rp->test_mode) { 254 PetscCall(PetscPrintf(rp->comm, " Performance:\n")); 255 } 256 { 257 // Set up error operator context 258 PetscCall(SetupErrorOperatorCtx(rp->comm, dm, ceed, ceed_data, X_loc, op_error, op_error_ctx)); 259 PetscScalar l2_error; 260 PetscCall(ComputeL2Error(X, &l2_error, op_error_ctx)); 261 PetscReal tol = 5e-2; 262 if (!rp->test_mode || l2_error > tol) { 263 PetscCall(MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, rp->comm)); 264 PetscCall(MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, rp->comm)); 265 PetscCall(PetscPrintf(rp->comm, 266 " L2 Error : %e\n" 267 " CG Solve Time : %g (%g) sec\n", 268 (double)l2_error, rt_max, rt_min)); 269 } 270 } 271 if (!rp->test_mode) { 272 PetscCall(PetscPrintf(rp->comm, " DoFs/Sec in CG : %g (%g) million\n", 1e-6 * g_size * its / rt_max, 273 1e-6 * g_size * its / rt_min)); 274 } 275 } 276 277 if (rp->write_solution) { 278 PetscViewer vtk_viewer_soln; 279 280 PetscCall(PetscViewerCreate(rp->comm, &vtk_viewer_soln)); 281 PetscCall(PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK)); 282 PetscCall(PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu")); 283 PetscCall(VecView(X, vtk_viewer_soln)); 284 PetscCall(PetscViewerDestroy(&vtk_viewer_soln)); 285 } 286 287 // Cleanup 288 PetscCall(VecDestroy(&X)); 289 PetscCall(VecDestroy(&X_loc)); 290 PetscCall(VecDestroy(&op_apply_ctx->Y_loc)); 291 PetscCall(VecDestroy(&op_error_ctx->Y_loc)); 292 PetscCall(MatDestroy(&mat_O)); 293 PetscCall(PetscFree(op_apply_ctx)); 294 PetscCall(PetscFree(op_error_ctx)); 295 PetscCall(CeedDataDestroy(0, ceed_data)); 296 297 PetscCall(VecDestroy(&rhs)); 298 PetscCall(VecDestroy(&rhs_loc)); 299 PetscCall(KSPDestroy(&ksp)); 300 CeedVectorDestroy(&target); 301 CeedQFunctionDestroy(&qf_error); 302 CeedOperatorDestroy(&op_error); 303 CeedDestroy(&ceed); 304 PetscFunctionReturn(0); 305 } 306 307 static PetscErrorCode Run(RunParams rp, PetscInt num_resources, char *const *ceed_resources, PetscInt num_bp_choices, const BPType *bp_choices) { 308 DM dm; 309 310 PetscFunctionBeginUser; 311 // Setup DM 312 PetscCall(CreateDistributedDM(rp, &dm)); 313 314 for (PetscInt b = 0; b < num_bp_choices; b++) { 315 DM dm_deg; 316 VecType vec_type; 317 PetscInt q_extra = rp->q_extra; 318 rp->bp_choice = bp_choices[b]; 319 rp->num_comp_u = bp_options[rp->bp_choice].num_comp_u; 320 rp->q_extra = q_extra < 0 ? bp_options[rp->bp_choice].q_extra : q_extra; 321 PetscCall(DMClone(dm, &dm_deg)); 322 PetscCall(DMGetVecType(dm, &vec_type)); 323 PetscCall(DMSetVecType(dm_deg, vec_type)); 324 // Create DM 325 PetscInt dim; 326 PetscCall(DMGetDimension(dm_deg, &dim)); 327 PetscCall(SetupDMByDegree(dm_deg, rp->degree, rp->q_extra, rp->num_comp_u, dim, bp_options[rp->bp_choice].enforce_bc)); 328 for (PetscInt r = 0; r < num_resources; r++) { 329 PetscCall(RunWithDM(rp, dm_deg, ceed_resources[r])); 330 } 331 PetscCall(DMDestroy(&dm_deg)); 332 rp->q_extra = q_extra; 333 } 334 335 PetscCall(DMDestroy(&dm)); 336 PetscFunctionReturn(0); 337 } 338 339 int main(int argc, char **argv) { 340 PetscInt comm_size; 341 RunParams rp; 342 MPI_Comm comm; 343 char filename[PETSC_MAX_PATH_LEN]; 344 char *ceed_resources[30]; 345 PetscInt num_ceed_resources = 30; 346 char hostname[PETSC_MAX_PATH_LEN]; 347 348 PetscInt dim = 3, mesh_elem[3] = {3, 3, 3}; 349 PetscInt num_degrees = 30, degree[30] = {}, num_local_nodes = 2, local_nodes[2] = {}; 350 PetscMPIInt ranks_per_node; 351 PetscBool degree_set; 352 BPType bp_choices[10]; 353 PetscInt num_bp_choices = 10; 354 355 // Initialize PETSc 356 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 357 comm = PETSC_COMM_WORLD; 358 PetscCall(MPI_Comm_size(comm, &comm_size)); 359 #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) 360 { 361 MPI_Comm splitcomm; 362 PetscCall(MPI_Comm_split_type(comm, MPI_COMM_TYPE_SHARED, 0, MPI_INFO_NULL, &splitcomm)); 363 PetscCall(MPI_Comm_size(splitcomm, &ranks_per_node)); 364 PetscCall(MPI_Comm_free(&splitcomm)); 365 } 366 #else 367 ranks_per_node = -1; // Unknown 368 #endif 369 370 // Setup all parameters needed in Run() 371 PetscCall(PetscMalloc1(1, &rp)); 372 rp->comm = comm; 373 374 // Read command line options 375 PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); 376 { 377 PetscBool set; 378 PetscCall(PetscOptionsEnumArray("-problem", "CEED benchmark problem to solve", NULL, bp_types, (PetscEnum *)bp_choices, &num_bp_choices, &set)); 379 if (!set) { 380 bp_choices[0] = CEED_BP1; 381 num_bp_choices = 1; 382 } 383 } 384 rp->test_mode = PETSC_FALSE; 385 PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, rp->test_mode, &rp->test_mode, NULL)); 386 rp->write_solution = PETSC_FALSE; 387 PetscCall(PetscOptionsBool("-write_solution", "Write solution for visualization", NULL, rp->write_solution, &rp->write_solution, NULL)); 388 rp->simplex = PETSC_FALSE; 389 PetscCall(PetscOptionsBool("-simplex", "Element topology (default:hex)", NULL, rp->simplex, &rp->simplex, NULL)); 390 if ((bp_choices[0] == CEED_BP5 || bp_choices[0] == CEED_BP6) && (rp->simplex)) { 391 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "BP5/6 is not supported with simplex"); 392 } 393 degree[0] = rp->test_mode ? 3 : 2; 394 PetscCall(PetscOptionsIntArray("-degree", "Polynomial degree of tensor product basis", NULL, degree, &num_degrees, °ree_set)); 395 if (!degree_set) num_degrees = 1; 396 rp->q_extra = PETSC_DECIDE; 397 PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points (-1 for auto)", NULL, rp->q_extra, &rp->q_extra, NULL)); 398 { 399 PetscBool set; 400 PetscCall(PetscOptionsStringArray("-ceed", "CEED resource specifier (comma-separated list)", NULL, ceed_resources, &num_ceed_resources, &set)); 401 if (!set) { 402 PetscCall(PetscStrallocpy("/cpu/self", &ceed_resources[0])); 403 num_ceed_resources = 1; 404 } 405 } 406 PetscCall(PetscGetHostName(hostname, sizeof hostname)); 407 PetscCall(PetscOptionsString("-hostname", "Hostname for output", NULL, hostname, hostname, sizeof(hostname), NULL)); 408 rp->read_mesh = PETSC_FALSE; 409 PetscCall(PetscOptionsString("-mesh", "Read mesh from file", NULL, filename, filename, sizeof(filename), &rp->read_mesh)); 410 rp->filename = filename; 411 if (!rp->read_mesh) { 412 PetscInt tmp = dim; 413 PetscCall(PetscOptionsIntArray("-cells", "Number of cells per dimension", NULL, mesh_elem, &tmp, NULL)); 414 } 415 local_nodes[0] = 1000; 416 PetscCall(PetscOptionsIntArray("-local_nodes", 417 "Target number of locally owned nodes per " 418 "process (single value or min,max)", 419 NULL, local_nodes, &num_local_nodes, &rp->user_l_nodes)); 420 if (num_local_nodes < 2) local_nodes[1] = 2 * local_nodes[0]; 421 { 422 PetscInt two = 2; 423 rp->ksp_max_it_clip[0] = 5; 424 rp->ksp_max_it_clip[1] = 20; 425 PetscCall(PetscOptionsIntArray("-ksp_max_it_clip", "Min and max number of iterations to use during benchmarking", NULL, rp->ksp_max_it_clip, &two, 426 NULL)); 427 } 428 if (!degree_set) { 429 PetscInt max_degree = 8; 430 PetscCall(PetscOptionsInt("-max_degree", "Range of degrees [1, max_degree] to run with", NULL, max_degree, &max_degree, NULL)); 431 for (PetscInt i = 0; i < max_degree; i++) degree[i] = i + 1; 432 num_degrees = max_degree; 433 } 434 { 435 PetscBool flg; 436 PetscInt p = ranks_per_node; 437 PetscCall(PetscOptionsInt("-p", "Number of MPI ranks per node", NULL, p, &p, &flg)); 438 if (flg) ranks_per_node = p; 439 } 440 441 PetscOptionsEnd(); 442 443 // Register PETSc logging stage 444 PetscCall(PetscLogStageRegister("Solve Stage", &rp->solve_stage)); 445 446 rp->hostname = hostname; 447 rp->dim = dim; 448 rp->mesh_elem = mesh_elem; 449 rp->ranks_per_node = ranks_per_node; 450 451 for (PetscInt d = 0; d < num_degrees; d++) { 452 PetscInt deg = degree[d]; 453 for (PetscInt n = local_nodes[0]; n < local_nodes[1]; n *= 2) { 454 rp->degree = deg; 455 rp->local_nodes = n; 456 PetscCall(Run(rp, num_ceed_resources, ceed_resources, num_bp_choices, bp_choices)); 457 } 458 } 459 // Clear memory 460 PetscCall(PetscFree(rp)); 461 for (PetscInt i = 0; i < num_ceed_resources; i++) { 462 PetscCall(PetscFree(ceed_resources[i])); 463 } 464 return PetscFinalize(); 465 } 466