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 " Solution Order (P) : %" CeedInt_FMT "\n" 152 " Quadrature Order (Q) : %" CeedInt_FMT "\n" 153 " Additional quadrature points (q_extra) : %" CeedInt_FMT "\n" 154 " Global nodes : %" PetscInt_FMT "\n" 155 " Local Elements : %" PetscInt_FMT "\n" 156 " Element topology : %s\n" 157 " Owned nodes : %" PetscInt_FMT "\n" 158 " DoF per node : %" PetscInt_FMT "\n", 159 rp->bp_choice+1, rp->hostname, comm_size, 160 rp->ranks_per_node, vec_type, used_resource, 161 CeedMemTypes[mem_type_backend], P, Q, rp->q_extra, 162 g_size/rp->num_comp_u, c_end - c_start, 163 CeedElemTopologies[elem_topo], 164 l_size/rp->num_comp_u, rp->num_comp_u); 165 CHKERRQ(ierr); 166 } 167 168 // Create RHS vector 169 ierr = VecDuplicate(X_loc, &rhs_loc); CHKERRQ(ierr); 170 ierr = VecZeroEntries(rhs_loc); CHKERRQ(ierr); 171 ierr = VecGetArrayAndMemType(rhs_loc, &r, &mem_type); CHKERRQ(ierr); 172 CeedVectorCreate(ceed, xl_size, &rhs_ceed); 173 CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r); 174 175 ierr = PetscMalloc1(1, &ceed_data); CHKERRQ(ierr); 176 ierr = SetupLibceedByDegree(dm, ceed, rp->degree, rp->dim, rp->q_extra, 177 rp->dim, rp->num_comp_u, g_size, xl_size, bp_options[rp->bp_choice], 178 ceed_data, true, rhs_ceed, &target); CHKERRQ(ierr); 179 180 // Gather RHS 181 CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL); 182 ierr = VecRestoreArrayAndMemType(rhs_loc, &r); CHKERRQ(ierr); 183 ierr = VecZeroEntries(rhs); CHKERRQ(ierr); 184 ierr = DMLocalToGlobal(dm, rhs_loc, ADD_VALUES, rhs); CHKERRQ(ierr); 185 CeedVectorDestroy(&rhs_ceed); 186 187 // Create the error QFunction 188 CeedQFunctionCreateInterior(ceed, 1, bp_options[rp->bp_choice].error, 189 bp_options[rp->bp_choice].error_loc, &qf_error); 190 CeedQFunctionAddInput(qf_error, "u", rp->num_comp_u, CEED_EVAL_INTERP); 191 CeedQFunctionAddInput(qf_error, "true_soln", rp->num_comp_u, CEED_EVAL_NONE); 192 CeedQFunctionAddOutput(qf_error, "error", rp->num_comp_u, CEED_EVAL_NONE); 193 194 // Create the error operator 195 CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, 196 &op_error); 197 CeedOperatorSetField(op_error, "u", ceed_data->elem_restr_u, 198 ceed_data->basis_u, CEED_VECTOR_ACTIVE); 199 CeedOperatorSetField(op_error, "true_soln", ceed_data->elem_restr_u_i, 200 CEED_BASIS_COLLOCATED, target); 201 CeedOperatorSetField(op_error, "error", ceed_data->elem_restr_u_i, 202 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 203 204 // Set up Mat 205 user_O->comm = rp->comm; 206 user_O->dm = dm; 207 user_O->X_loc = X_loc; 208 ierr = VecDuplicate(X_loc, &user_O->Y_loc); CHKERRQ(ierr); 209 user_O->x_ceed = ceed_data->x_ceed; 210 user_O->y_ceed = ceed_data->y_ceed; 211 user_O->op = ceed_data->op_apply; 212 user_O->ceed = ceed; 213 214 ierr = KSPCreate(rp->comm, &ksp); CHKERRQ(ierr); 215 { 216 PC pc; 217 ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 218 if (rp->bp_choice == CEED_BP1 || rp->bp_choice == CEED_BP2) { 219 ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr); 220 if (rp->simplex) { 221 ierr = PCJacobiSetType(pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 222 } else { 223 ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr); 224 } 225 } else { 226 ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr); 227 } 228 ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 229 ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 230 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 231 PETSC_DEFAULT); CHKERRQ(ierr); 232 } 233 ierr = KSPSetOperators(ksp, mat_O, mat_O); CHKERRQ(ierr); 234 235 // First run's performance log is not considered for benchmarking purposes 236 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1); 237 CHKERRQ(ierr); 238 my_rt_start = MPI_Wtime(); 239 ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 240 my_rt = MPI_Wtime() - my_rt_start; 241 ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, rp->comm); 242 CHKERRQ(ierr); 243 // Set maxits based on first iteration timing 244 if (my_rt > 0.02) { 245 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 246 rp->ksp_max_it_clip[0]); 247 CHKERRQ(ierr); 248 } else { 249 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 250 rp->ksp_max_it_clip[1]); 251 CHKERRQ(ierr); 252 } 253 ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr); 254 255 // Timed solve 256 ierr = VecZeroEntries(X); CHKERRQ(ierr); 257 ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr); 258 259 // -- Performance logging 260 ierr = PetscLogStagePush(rp->solve_stage); CHKERRQ(ierr); 261 262 // -- Solve 263 my_rt_start = MPI_Wtime(); 264 ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 265 my_rt = MPI_Wtime() - my_rt_start; 266 267 // -- Performance logging 268 ierr = PetscLogStagePop(); 269 270 // Output results 271 { 272 KSPType ksp_type; 273 KSPConvergedReason reason; 274 PetscReal rnorm; 275 PetscInt its; 276 ierr = KSPGetType(ksp, &ksp_type); CHKERRQ(ierr); 277 ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr); 278 ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr); 279 ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr); 280 if (!rp->test_mode || reason < 0 || rnorm > 1e-8) { 281 ierr = PetscPrintf(rp->comm, 282 " KSP:\n" 283 " KSP Type : %s\n" 284 " KSP Convergence : %s\n" 285 " Total KSP Iterations : %" PetscInt_FMT "\n" 286 " Final rnorm : %e\n", 287 ksp_type, KSPConvergedReasons[reason], its, 288 (double)rnorm); CHKERRQ(ierr); 289 } 290 if (!rp->test_mode) { 291 ierr = PetscPrintf(rp->comm," Performance:\n"); CHKERRQ(ierr); 292 } 293 { 294 PetscReal max_error; 295 ierr = ComputeErrorMax(user_O, op_error, X, target, &max_error); 296 CHKERRQ(ierr); 297 PetscReal tol = 5e-2; 298 if (!rp->test_mode || max_error > tol) { 299 ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, rp->comm); 300 CHKERRQ(ierr); 301 ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, rp->comm); 302 CHKERRQ(ierr); 303 ierr = PetscPrintf(rp->comm, 304 " Pointwise Error (max) : %e\n" 305 " CG Solve Time : %g (%g) sec\n", 306 (double)max_error, rt_max, rt_min); CHKERRQ(ierr); 307 } 308 } 309 if (!rp->test_mode) { 310 ierr = PetscPrintf(rp->comm, 311 " DoFs/Sec in CG : %g (%g) million\n", 312 1e-6*g_size*its/rt_max, 313 1e-6*g_size*its/rt_min); CHKERRQ(ierr); 314 } 315 } 316 317 if (rp->write_solution) { 318 PetscViewer vtk_viewer_soln; 319 320 ierr = PetscViewerCreate(rp->comm, &vtk_viewer_soln); CHKERRQ(ierr); 321 ierr = PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK); CHKERRQ(ierr); 322 ierr = PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"); CHKERRQ(ierr); 323 ierr = VecView(X, vtk_viewer_soln); CHKERRQ(ierr); 324 ierr = PetscViewerDestroy(&vtk_viewer_soln); CHKERRQ(ierr); 325 } 326 327 // Cleanup 328 ierr = VecDestroy(&X); CHKERRQ(ierr); 329 ierr = VecDestroy(&X_loc); CHKERRQ(ierr); 330 ierr = VecDestroy(&user_O->Y_loc); CHKERRQ(ierr); 331 ierr = MatDestroy(&mat_O); CHKERRQ(ierr); 332 ierr = PetscFree(user_O); CHKERRQ(ierr); 333 ierr = CeedDataDestroy(0, ceed_data); CHKERRQ(ierr); 334 335 ierr = VecDestroy(&rhs); CHKERRQ(ierr); 336 ierr = VecDestroy(&rhs_loc); CHKERRQ(ierr); 337 ierr = KSPDestroy(&ksp); CHKERRQ(ierr); 338 CeedVectorDestroy(&target); 339 CeedQFunctionDestroy(&qf_error); 340 CeedOperatorDestroy(&op_error); 341 CeedDestroy(&ceed); 342 PetscFunctionReturn(0); 343 } 344 345 static PetscErrorCode Run(RunParams rp, PetscInt num_resources, 346 char *const *ceed_resources, PetscInt num_bp_choices, 347 const BPType *bp_choices) { 348 PetscInt ierr; 349 DM dm; 350 351 PetscFunctionBeginUser; 352 // Setup DM 353 ierr = CreateDistributedDM(rp, &dm); CHKERRQ(ierr); 354 355 for (PetscInt b = 0; b < num_bp_choices; b++) { 356 DM dm_deg; 357 VecType vec_type; 358 PetscInt q_extra = rp->q_extra; 359 rp->bp_choice = bp_choices[b]; 360 rp->num_comp_u = bp_options[rp->bp_choice].num_comp_u; 361 rp->q_extra = q_extra < 0 ? bp_options[rp->bp_choice].q_extra : q_extra; 362 ierr = DMClone(dm, &dm_deg); CHKERRQ(ierr); 363 ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr); 364 ierr = DMSetVecType(dm_deg, vec_type); CHKERRQ(ierr); 365 // Create DM 366 PetscInt dim; 367 ierr = DMGetDimension(dm_deg, &dim); CHKERRQ(ierr); 368 ierr = SetupDMByDegree(dm_deg, rp->degree, q_extra, rp->num_comp_u, dim, 369 bp_options[rp->bp_choice].enforce_bc, 370 bp_options[rp->bp_choice].bc_func); CHKERRQ(ierr); 371 for (PetscInt r = 0; r < num_resources; r++) { 372 ierr = RunWithDM(rp, dm_deg, ceed_resources[r]); CHKERRQ(ierr); 373 } 374 ierr = DMDestroy(&dm_deg); CHKERRQ(ierr); 375 rp->q_extra = q_extra; 376 } 377 378 ierr = DMDestroy(&dm); CHKERRQ(ierr); 379 PetscFunctionReturn(0); 380 } 381 382 int main(int argc, char **argv) { 383 PetscInt ierr, comm_size; 384 RunParams rp; 385 MPI_Comm comm; 386 char filename[PETSC_MAX_PATH_LEN]; 387 char *ceed_resources[30]; 388 PetscInt num_ceed_resources = 30; 389 char hostname[PETSC_MAX_PATH_LEN]; 390 391 PetscInt dim = 3, mesh_elem[3] = {3, 3, 3}; 392 PetscInt num_degrees = 30, degree[30] = {}, num_local_nodes = 2, 393 local_nodes[2] = {}; 394 PetscMPIInt ranks_per_node; 395 PetscBool degree_set; 396 BPType bp_choices[10]; 397 PetscInt num_bp_choices = 10; 398 399 // Initialize PETSc 400 ierr = PetscInitialize(&argc, &argv, NULL, help); 401 if (ierr) return ierr; 402 comm = PETSC_COMM_WORLD; 403 ierr = MPI_Comm_size(comm, &comm_size); 404 if (ierr != MPI_SUCCESS) return ierr; 405 #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) 406 { 407 MPI_Comm splitcomm; 408 ierr = MPI_Comm_split_type(comm, MPI_COMM_TYPE_SHARED, 0, MPI_INFO_NULL, 409 &splitcomm); 410 CHKERRQ(ierr); 411 ierr = MPI_Comm_size(splitcomm, &ranks_per_node); CHKERRQ(ierr); 412 ierr = MPI_Comm_free(&splitcomm); CHKERRQ(ierr); 413 } 414 #else 415 ranks_per_node = -1; // Unknown 416 #endif 417 418 // Setup all parameters needed in Run() 419 ierr = PetscMalloc1(1, &rp); CHKERRQ(ierr); 420 rp->comm = comm; 421 422 // Read command line options 423 PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); 424 { 425 PetscBool set; 426 ierr = PetscOptionsEnumArray("-problem", "CEED benchmark problem to solve", 427 NULL, 428 bp_types, (PetscEnum *)bp_choices, &num_bp_choices, &set); 429 CHKERRQ(ierr); 430 if (!set) { 431 bp_choices[0] = CEED_BP1; 432 num_bp_choices = 1; 433 } 434 } 435 rp->test_mode = PETSC_FALSE; 436 ierr = PetscOptionsBool("-test", 437 "Testing mode (do not print unless error is large)", 438 NULL, rp->test_mode, &rp->test_mode, NULL); CHKERRQ(ierr); 439 rp->write_solution = PETSC_FALSE; 440 ierr = PetscOptionsBool("-write_solution", "Write solution for visualization", 441 NULL, rp->write_solution, &rp->write_solution, NULL); 442 CHKERRQ(ierr); 443 rp->simplex = PETSC_FALSE; 444 ierr = PetscOptionsBool("-simplex", "Element topology (default:hex)", 445 NULL, rp->simplex, &rp->simplex, NULL); 446 CHKERRQ(ierr); 447 degree[0] = rp->test_mode ? 3 : 2; 448 ierr = PetscOptionsIntArray("-degree", 449 "Polynomial degree of tensor product basis", NULL, 450 degree, &num_degrees, °ree_set); CHKERRQ(ierr); 451 if (!degree_set) 452 num_degrees = 1; 453 rp->q_extra = PETSC_DECIDE; 454 ierr = PetscOptionsInt("-q_extra", 455 "Number of extra quadrature points (-1 for auto)", NULL, 456 rp->q_extra, &rp->q_extra, NULL); CHKERRQ(ierr); 457 { 458 PetscBool set; 459 ierr = PetscOptionsStringArray("-ceed", 460 "CEED resource specifier (comma-separated list)", NULL, 461 ceed_resources, &num_ceed_resources, &set); CHKERRQ(ierr); 462 if (!set) { 463 ierr = PetscStrallocpy( "/cpu/self", &ceed_resources[0]); CHKERRQ(ierr); 464 num_ceed_resources = 1; 465 } 466 } 467 ierr = PetscGetHostName(hostname, sizeof hostname); CHKERRQ(ierr); 468 ierr = PetscOptionsString("-hostname", "Hostname for output", NULL, hostname, 469 hostname, sizeof(hostname), NULL); CHKERRQ(ierr); 470 rp->read_mesh = PETSC_FALSE; 471 ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL, filename, 472 filename, sizeof(filename), &rp->read_mesh); 473 CHKERRQ(ierr); 474 rp->filename = filename; 475 if (!rp->read_mesh) { 476 PetscInt tmp = dim; 477 ierr = PetscOptionsIntArray("-cells", "Number of cells per dimension", NULL, 478 mesh_elem, &tmp, NULL); CHKERRQ(ierr); 479 } 480 local_nodes[0] = 1000; 481 ierr = PetscOptionsIntArray("-local_nodes", 482 "Target number of locally owned nodes per " 483 "process (single value or min,max)", 484 NULL, local_nodes, &num_local_nodes, &rp->user_l_nodes); 485 CHKERRQ(ierr); 486 if (num_local_nodes < 2) 487 local_nodes[1] = 2 * local_nodes[0]; 488 { 489 PetscInt two = 2; 490 rp->ksp_max_it_clip[0] = 5; 491 rp->ksp_max_it_clip[1] = 20; 492 ierr = PetscOptionsIntArray("-ksp_max_it_clip", 493 "Min and max number of iterations to use during benchmarking", 494 NULL, rp->ksp_max_it_clip, &two, NULL); CHKERRQ(ierr); 495 } 496 if (!degree_set) { 497 PetscInt max_degree = 8; 498 ierr = PetscOptionsInt("-max_degree", 499 "Range of degrees [1, max_degree] to run with", 500 NULL, max_degree, &max_degree, NULL); 501 CHKERRQ(ierr); 502 for (PetscInt i = 0; i < max_degree; i++) 503 degree[i] = i + 1; 504 num_degrees = max_degree; 505 } 506 { 507 PetscBool flg; 508 PetscInt p = ranks_per_node; 509 ierr = PetscOptionsInt("-p", "Number of MPI ranks per node", NULL, 510 p, &p, &flg); 511 CHKERRQ(ierr); 512 if (flg) ranks_per_node = p; 513 } 514 515 PetscOptionsEnd(); 516 517 // Register PETSc logging stage 518 ierr = PetscLogStageRegister("Solve Stage", &rp->solve_stage); 519 CHKERRQ(ierr); 520 521 rp->hostname = hostname; 522 rp->dim = dim; 523 rp->mesh_elem = mesh_elem; 524 rp->ranks_per_node = ranks_per_node; 525 526 for (PetscInt d = 0; d < num_degrees; d++) { 527 PetscInt deg = degree[d]; 528 for (PetscInt n = local_nodes[0]; n < local_nodes[1]; n *= 2) { 529 rp->degree = deg; 530 rp->local_nodes = n; 531 ierr = Run(rp, num_ceed_resources, ceed_resources, 532 num_bp_choices, bp_choices); CHKERRQ(ierr); 533 } 534 } 535 // Clear memory 536 ierr = PetscFree(rp); CHKERRQ(ierr); 537 for (PetscInt i=0; i<num_ceed_resources; i++) { 538 ierr = PetscFree(ceed_resources[i]); CHKERRQ(ierr); 539 } 540 return PetscFinalize(); 541 } 542