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