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