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 // on a closed surface, such as the one of a discrete sphere. 22 // 23 // The code uses higher level communication protocols in DMPlex. 24 // 25 // Build with: 26 // 27 // make bpssphere [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] 28 // 29 // Sample runs: 30 // 31 // bpssphere -problem bp1 -degree 3 32 // bpssphere -problem bp2 -degree 3 33 // bpssphere -problem bp3 -degree 3 34 // bpssphere -problem bp4 -degree 3 35 // bpssphere -problem bp5 -degree 3 -ceed /cpu/self 36 // bpssphere -problem bp6 -degree 3 -ceed /gpu/cuda 37 // 38 //TESTARGS -ceed {ceed_resource} -test -problem bp3 -degree 3 -dm_refine 2 39 40 /// @file 41 /// CEED BPs example using PETSc with DMPlex 42 /// See bps.c for a "raw" implementation using a structured grid. 43 /// and bpsdmplex.c for an implementation using an unstructured grid. 44 static const char help[] = "Solve CEED BPs on a sphere using DMPlex in PETSc\n"; 45 46 #include <ceed.h> 47 #include <petscdmplex.h> 48 #include <petscksp.h> 49 #include <stdbool.h> 50 #include <string.h> 51 #include "bpssphere.h" 52 53 int main(int argc, char **argv) { 54 PetscInt ierr; 55 MPI_Comm comm; 56 char ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self", 57 filename[PETSC_MAX_PATH_LEN]; 58 double my_rt_start, my_rt, rt_min, rt_max; 59 PetscInt degree = 3, q_extra, l_size, g_size, topo_dim = 2, num_comp_x = 3, 60 num_comp_u = 1, xl_size; 61 PetscScalar *r; 62 PetscBool test_mode, benchmark_mode, read_mesh, write_solution, simplex; 63 PetscLogStage solve_stage; 64 Vec X, X_loc, rhs, rhs_loc; 65 Mat mat_O; 66 KSP ksp; 67 DM dm; 68 UserO user_O; 69 Ceed ceed; 70 CeedData ceed_data; 71 CeedQFunction qf_error; 72 CeedOperator op_error; 73 CeedVector rhs_ceed, target; 74 BPType bp_choice; 75 VecType vec_type; 76 PetscMemType mem_type; 77 78 ierr = PetscInitialize(&argc, &argv, NULL, help); 79 if (ierr) return ierr; 80 comm = PETSC_COMM_WORLD; 81 82 // Read command line options 83 ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr); 84 bp_choice = CEED_BP1; 85 ierr = PetscOptionsEnum("-problem", 86 "CEED benchmark problem to solve", NULL, 87 bp_types, (PetscEnum)bp_choice, (PetscEnum *)&bp_choice, 88 NULL); CHKERRQ(ierr); 89 num_comp_u = bp_options[bp_choice].num_comp_u; 90 test_mode = PETSC_FALSE; 91 ierr = PetscOptionsBool("-test", 92 "Testing mode (do not print unless error is large)", 93 NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr); 94 benchmark_mode = PETSC_FALSE; 95 ierr = PetscOptionsBool("-benchmark", 96 "Benchmarking mode (prints benchmark statistics)", 97 NULL, benchmark_mode, &benchmark_mode, NULL); 98 CHKERRQ(ierr); 99 write_solution = PETSC_FALSE; 100 ierr = PetscOptionsBool("-write_solution", 101 "Write solution for visualization", 102 NULL, write_solution, &write_solution, NULL); 103 CHKERRQ(ierr); 104 degree = test_mode ? 3 : 2; 105 ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", 106 NULL, degree, °ree, NULL); CHKERRQ(ierr); 107 q_extra = bp_options[bp_choice].q_extra; 108 ierr = PetscOptionsInt("-q_extra", "Number of extra quadrature points", 109 NULL, q_extra, &q_extra, NULL); CHKERRQ(ierr); 110 ierr = PetscOptionsString("-ceed", "CEED resource specifier", 111 NULL, ceed_resource, ceed_resource, 112 sizeof(ceed_resource), NULL); CHKERRQ(ierr); 113 read_mesh = PETSC_FALSE; 114 ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL, 115 filename, filename, sizeof(filename), &read_mesh); 116 CHKERRQ(ierr); 117 simplex = PETSC_FALSE; 118 ierr = PetscOptionsBool("-simplex", "Use simplices, or tensor product cells", 119 NULL, simplex, &simplex, NULL); CHKERRQ(ierr); 120 ierr = PetscOptionsEnd(); CHKERRQ(ierr); 121 122 // Setup DM 123 if (read_mesh) { 124 ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dm); 125 CHKERRQ(ierr); 126 } else { 127 // Create the mesh as a 0-refined sphere. This will create a cubic surface, not a box 128 ierr = DMPlexCreateSphereMesh(PETSC_COMM_WORLD, topo_dim, simplex, 1., &dm); 129 CHKERRQ(ierr); 130 // Set the object name 131 ierr = PetscObjectSetName((PetscObject)dm, "Sphere"); CHKERRQ(ierr); 132 // Distribute mesh over processes 133 { 134 DM dm_dist = NULL; 135 PetscPartitioner part; 136 137 ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr); 138 ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 139 ierr = DMPlexDistribute(dm, 0, NULL, &dm_dist); CHKERRQ(ierr); 140 if (dm_dist) { 141 ierr = DMDestroy(&dm); CHKERRQ(ierr); 142 dm = dm_dist; 143 } 144 } 145 // Refine DMPlex with uniform refinement using runtime option -dm_refine 146 ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr); 147 ierr = DMSetFromOptions(dm); CHKERRQ(ierr); 148 ierr = ProjectToUnitSphere(dm); CHKERRQ(ierr); 149 // View DMPlex via runtime option 150 ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr); 151 } 152 153 // Create DM 154 ierr = SetupDMByDegree(dm, degree, num_comp_u, topo_dim, false, 155 (BCFunction)NULL); 156 CHKERRQ(ierr); 157 158 // Create vectors 159 ierr = DMCreateGlobalVector(dm, &X); CHKERRQ(ierr); 160 ierr = VecGetLocalSize(X, &l_size); CHKERRQ(ierr); 161 ierr = VecGetSize(X, &g_size); CHKERRQ(ierr); 162 ierr = DMCreateLocalVector(dm, &X_loc); CHKERRQ(ierr); 163 ierr = VecGetSize(X_loc, &xl_size); CHKERRQ(ierr); 164 ierr = VecDuplicate(X, &rhs); CHKERRQ(ierr); 165 166 // Operator 167 ierr = PetscMalloc1(1, &user_O); CHKERRQ(ierr); 168 ierr = MatCreateShell(comm, l_size, l_size, g_size, g_size, 169 user_O, &mat_O); CHKERRQ(ierr); 170 ierr = MatShellSetOperation(mat_O, MATOP_MULT, 171 (void(*)(void))MatMult_Ceed); CHKERRQ(ierr); 172 173 // Set up libCEED 174 CeedInit(ceed_resource, &ceed); 175 CeedMemType mem_type_backend; 176 CeedGetPreferredMemType(ceed, &mem_type_backend); 177 178 ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr); 179 if (!vec_type) { // Not yet set by user -dm_vec_type 180 switch (mem_type_backend) { 181 case CEED_MEM_HOST: vec_type = VECSTANDARD; break; 182 case CEED_MEM_DEVICE: { 183 const char *resolved; 184 CeedGetResource(ceed, &resolved); 185 if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA; 186 else if (strstr(resolved, "/gpu/hip/occa")) 187 vec_type = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678 188 else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP; 189 else vec_type = VECSTANDARD; 190 } 191 } 192 ierr = DMSetVecType(dm, vec_type); CHKERRQ(ierr); 193 } 194 195 // Print summary 196 if (!test_mode) { 197 PetscInt P = degree + 1, Q = P + q_extra; 198 const char *used_resource; 199 CeedGetResource(ceed, &used_resource); 200 ierr = PetscPrintf(comm, 201 "\n-- CEED Benchmark Problem %d on the Sphere -- libCEED + PETSc --\n" 202 " libCEED:\n" 203 " libCEED Backend : %s\n" 204 " libCEED Backend MemType : %s\n" 205 " Mesh:\n" 206 " Number of 1D Basis Nodes (p) : %d\n" 207 " Number of 1D Quadrature Points (q) : %d\n" 208 " Global nodes : %D\n", 209 bp_choice+1, ceed_resource, CeedMemTypes[mem_type_backend], P, Q, 210 g_size/num_comp_u); CHKERRQ(ierr); 211 } 212 213 // Create RHS vector 214 ierr = VecDuplicate(X_loc, &rhs_loc); CHKERRQ(ierr); 215 ierr = VecZeroEntries(rhs_loc); CHKERRQ(ierr); 216 ierr = VecGetArrayAndMemType(rhs_loc, &r, &mem_type); CHKERRQ(ierr); 217 CeedVectorCreate(ceed, xl_size, &rhs_ceed); 218 CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r); 219 220 // Setup libCEED's objects 221 ierr = PetscMalloc1(1, &ceed_data); CHKERRQ(ierr); 222 ierr = SetupLibceedByDegree(dm, ceed, degree, topo_dim, q_extra, num_comp_x, 223 num_comp_u, g_size, xl_size, bp_options[bp_choice], 224 ceed_data, true, rhs_ceed, &target); CHKERRQ(ierr); 225 226 // Gather RHS 227 CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL); 228 ierr = VecRestoreArrayAndMemType(rhs_loc, &r); CHKERRQ(ierr); 229 ierr = VecZeroEntries(rhs); CHKERRQ(ierr); 230 ierr = DMLocalToGlobal(dm, rhs_loc, ADD_VALUES, rhs); CHKERRQ(ierr); 231 CeedVectorDestroy(&rhs_ceed); 232 233 // Create the error Q-function 234 CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].error, 235 bp_options[bp_choice].error_loc, &qf_error); 236 CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP); 237 CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE); 238 CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_NONE); 239 240 // Create the error operator 241 CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error); 242 CeedOperatorSetField(op_error, "u", ceed_data->elem_restr_u, 243 ceed_data->basis_u, CEED_VECTOR_ACTIVE); 244 CeedOperatorSetField(op_error, "true_soln", ceed_data->elem_restr_u_i, 245 CEED_BASIS_COLLOCATED, target); 246 CeedOperatorSetField(op_error, "error", ceed_data->elem_restr_u_i, 247 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 248 249 // Set up Mat 250 user_O->comm = comm; 251 user_O->dm = dm; 252 user_O->X_loc = X_loc; 253 ierr = VecDuplicate(X_loc, &user_O->Y_loc); CHKERRQ(ierr); 254 user_O->x_ceed = ceed_data->x_ceed; 255 user_O->y_ceed = ceed_data->y_ceed; 256 user_O->op = ceed_data->op_apply; 257 user_O->ceed = ceed; 258 259 // Setup solver 260 ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr); 261 { 262 PC pc; 263 ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 264 if (bp_choice == CEED_BP1 || bp_choice == CEED_BP2) { 265 ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr); 266 ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr); 267 } else { 268 ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr); 269 MatNullSpace nullspace; 270 271 ierr = MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, 0, &nullspace); 272 CHKERRQ(ierr); 273 ierr = MatSetNullSpace(mat_O, nullspace); CHKERRQ(ierr); 274 ierr = MatNullSpaceDestroy(&nullspace); CHKERRQ(ierr); 275 } 276 ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 277 ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 278 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 279 PETSC_DEFAULT); CHKERRQ(ierr); 280 } 281 ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr); 282 ierr = KSPSetOperators(ksp, mat_O, mat_O); CHKERRQ(ierr); 283 284 // First run, if benchmarking 285 if (benchmark_mode) { 286 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1); 287 CHKERRQ(ierr); 288 my_rt_start = MPI_Wtime(); 289 ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 290 my_rt = MPI_Wtime() - my_rt_start; 291 ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm); 292 CHKERRQ(ierr); 293 // Set maxits based on first iteration timing 294 if (my_rt > 0.02) { 295 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5); 296 CHKERRQ(ierr); 297 } else { 298 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20); 299 CHKERRQ(ierr); 300 } 301 } 302 303 // Timed solve 304 ierr = VecZeroEntries(X); CHKERRQ(ierr); 305 ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr); 306 307 // -- Performance logging 308 ierr = PetscLogStageRegister("Solve Stage", &solve_stage); CHKERRQ(ierr); 309 ierr = PetscLogStagePush(solve_stage); CHKERRQ(ierr); 310 311 // -- Solve 312 my_rt_start = MPI_Wtime(); 313 ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 314 my_rt = MPI_Wtime() - my_rt_start; 315 316 // -- Performance logging 317 ierr = PetscLogStagePop(); 318 319 // Output results 320 { 321 KSPType ksp_type; 322 KSPConvergedReason reason; 323 PetscReal rnorm; 324 PetscInt its; 325 ierr = KSPGetType(ksp, &ksp_type); CHKERRQ(ierr); 326 ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr); 327 ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr); 328 ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr); 329 if (!test_mode || reason < 0 || rnorm > 1e-8) { 330 ierr = PetscPrintf(comm, 331 " KSP:\n" 332 " KSP Type : %s\n" 333 " KSP Convergence : %s\n" 334 " Total KSP Iterations : %D\n" 335 " Final rnorm : %e\n", 336 ksp_type, KSPConvergedReasons[reason], its, 337 (double)rnorm); CHKERRQ(ierr); 338 } 339 if (!test_mode) { 340 ierr = PetscPrintf(comm," Performance:\n"); CHKERRQ(ierr); 341 } 342 { 343 PetscReal max_error; 344 ierr = ComputeErrorMax(user_O, op_error, X, target, &max_error); 345 CHKERRQ(ierr); 346 PetscReal tol = 5e-4; 347 if (!test_mode || max_error > tol) { 348 ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm); 349 CHKERRQ(ierr); 350 ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm); 351 CHKERRQ(ierr); 352 ierr = PetscPrintf(comm, 353 " Pointwise Error (max) : %e\n" 354 " CG Solve Time : %g (%g) sec\n", 355 (double)max_error, rt_max, rt_min); CHKERRQ(ierr); 356 } 357 } 358 if (benchmark_mode && (!test_mode)) { 359 ierr = PetscPrintf(comm, 360 " DoFs/Sec in CG : %g (%g) million\n", 361 1e-6*g_size*its/rt_max, 1e-6*g_size*its/rt_min); CHKERRQ(ierr); 362 } 363 } 364 365 // Output solution 366 if (write_solution) { 367 PetscViewer vtk_viewer_soln; 368 369 ierr = PetscViewerCreate(comm, &vtk_viewer_soln); CHKERRQ(ierr); 370 ierr = PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK); CHKERRQ(ierr); 371 ierr = PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"); CHKERRQ(ierr); 372 ierr = VecView(X, vtk_viewer_soln); CHKERRQ(ierr); 373 ierr = PetscViewerDestroy(&vtk_viewer_soln); CHKERRQ(ierr); 374 } 375 376 // Cleanup 377 ierr = VecDestroy(&X); CHKERRQ(ierr); 378 ierr = VecDestroy(&X_loc); CHKERRQ(ierr); 379 ierr = VecDestroy(&user_O->Y_loc); CHKERRQ(ierr); 380 ierr = MatDestroy(&mat_O); CHKERRQ(ierr); 381 ierr = PetscFree(user_O); CHKERRQ(ierr); 382 ierr = CeedDataDestroy(0, ceed_data); CHKERRQ(ierr); 383 ierr = DMDestroy(&dm); CHKERRQ(ierr); 384 385 ierr = VecDestroy(&rhs); CHKERRQ(ierr); 386 ierr = VecDestroy(&rhs_loc); CHKERRQ(ierr); 387 ierr = KSPDestroy(&ksp); CHKERRQ(ierr); 388 CeedVectorDestroy(&target); 389 CeedQFunctionDestroy(&qf_error); 390 CeedOperatorDestroy(&op_error); 391 CeedDestroy(&ceed); 392 return PetscFinalize(); 393 } 394