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