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