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