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