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