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 ceedresource[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, qextra, lsize, gsize, topodim = 2, ncompx = 3, 60 ncompu = 1, xlsize; 61 PetscScalar *r; 62 PetscBool test_mode, benchmark_mode, read_mesh, write_solution, simplex; 63 PetscLogStage solvestage; 64 Vec X, Xloc, rhs, rhsloc; 65 Mat matO; 66 KSP ksp; 67 DM dm; 68 UserO userO; 69 Ceed ceed; 70 CeedData ceeddata; 71 CeedQFunction qfError; 72 CeedOperator opError; 73 CeedVector rhsceed, target; 74 bpType bpChoice; 75 VecType vectype; 76 PetscMemType memtype; 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 bpChoice = CEED_BP1; 85 ierr = PetscOptionsEnum("-problem", 86 "CEED benchmark problem to solve", NULL, 87 bpTypes, (PetscEnum)bpChoice, (PetscEnum *)&bpChoice, 88 NULL); CHKERRQ(ierr); 89 ncompu = bpOptions[bpChoice].ncompu; 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 qextra = bpOptions[bpChoice].qextra; 108 ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", 109 NULL, qextra, &qextra, NULL); CHKERRQ(ierr); 110 ierr = PetscOptionsString("-ceed", "CEED resource specifier", 111 NULL, ceedresource, ceedresource, 112 sizeof(ceedresource), 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, topodim, 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 dmDist = NULL; 135 PetscPartitioner part; 136 137 ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr); 138 ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 139 ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr); 140 if (dmDist) { 141 ierr = DMDestroy(&dm); CHKERRQ(ierr); 142 dm = dmDist; 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, ncompu, topodim, false, (BCFunction)NULL); 155 CHKERRQ(ierr); 156 157 // Create vectors 158 ierr = DMCreateGlobalVector(dm, &X); CHKERRQ(ierr); 159 ierr = VecGetLocalSize(X, &lsize); CHKERRQ(ierr); 160 ierr = VecGetSize(X, &gsize); CHKERRQ(ierr); 161 ierr = DMCreateLocalVector(dm, &Xloc); CHKERRQ(ierr); 162 ierr = VecGetSize(Xloc, &xlsize); CHKERRQ(ierr); 163 ierr = VecDuplicate(X, &rhs); CHKERRQ(ierr); 164 165 // Operator 166 ierr = PetscMalloc1(1, &userO); CHKERRQ(ierr); 167 ierr = MatCreateShell(comm, lsize, lsize, gsize, gsize, 168 userO, &matO); CHKERRQ(ierr); 169 ierr = MatShellSetOperation(matO, MATOP_MULT, 170 (void(*)(void))MatMult_Ceed); CHKERRQ(ierr); 171 172 // Set up libCEED 173 CeedInit(ceedresource, &ceed); 174 CeedMemType memtypebackend; 175 CeedGetPreferredMemType(ceed, &memtypebackend); 176 177 ierr = DMGetVecType(dm, &vectype); CHKERRQ(ierr); 178 if (!vectype) { // Not yet set by user -dm_vec_type 179 switch (memtypebackend) { 180 case CEED_MEM_HOST: vectype = VECSTANDARD; break; 181 case CEED_MEM_DEVICE: { 182 const char *resolved; 183 CeedGetResource(ceed, &resolved); 184 if (strstr(resolved, "/gpu/cuda")) vectype = VECCUDA; 185 else if (strstr(resolved, "/gpu/hip/occa")) 186 vectype = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678 187 else if (strstr(resolved, "/gpu/hip")) vectype = VECHIP; 188 else vectype = VECSTANDARD; 189 } 190 } 191 ierr = DMSetVecType(dm, vectype); CHKERRQ(ierr); 192 } 193 194 // Print summary 195 if (!test_mode) { 196 PetscInt P = degree + 1, Q = P + qextra; 197 const char *usedresource; 198 CeedGetResource(ceed, &usedresource); 199 ierr = PetscPrintf(comm, 200 "\n-- CEED Benchmark Problem %d on the Sphere -- libCEED + PETSc --\n" 201 " libCEED:\n" 202 " libCEED Backend : %s\n" 203 " libCEED Backend MemType : %s\n" 204 " Mesh:\n" 205 " Number of 1D Basis Nodes (p) : %d\n" 206 " Number of 1D Quadrature Points (q) : %d\n" 207 " Global nodes : %D\n", 208 bpChoice+1, ceedresource, CeedMemTypes[memtypebackend], P, Q, 209 gsize/ncompu); CHKERRQ(ierr); 210 } 211 212 // Create RHS vector 213 ierr = VecDuplicate(Xloc, &rhsloc); CHKERRQ(ierr); 214 ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr); 215 ierr = VecGetArrayAndMemType(rhsloc, &r, &memtype); CHKERRQ(ierr); 216 CeedVectorCreate(ceed, xlsize, &rhsceed); 217 CeedVectorSetArray(rhsceed, MemTypeP2C(memtype), CEED_USE_POINTER, r); 218 219 // Setup libCEED's objects 220 ierr = PetscMalloc1(1, &ceeddata); CHKERRQ(ierr); 221 ierr = SetupLibceedByDegree(dm, ceed, degree, topodim, qextra, 222 ncompx, ncompu, gsize, xlsize, bpOptions[bpChoice], 223 ceeddata, true, rhsceed, &target); CHKERRQ(ierr); 224 225 // Gather RHS 226 CeedVectorTakeArray(rhsceed, MemTypeP2C(memtype), NULL); 227 ierr = VecRestoreArrayAndMemType(rhsloc, &r); CHKERRQ(ierr); 228 ierr = VecZeroEntries(rhs); CHKERRQ(ierr); 229 ierr = DMLocalToGlobal(dm, rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr); 230 CeedVectorDestroy(&rhsceed); 231 232 // Create the error Q-function 233 CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].error, 234 bpOptions[bpChoice].errorfname, &qfError); 235 CeedQFunctionAddInput(qfError, "u", ncompu, CEED_EVAL_INTERP); 236 CeedQFunctionAddInput(qfError, "true_soln", ncompu, CEED_EVAL_NONE); 237 CeedQFunctionAddOutput(qfError, "error", ncompu, CEED_EVAL_NONE); 238 239 // Create the error operator 240 CeedOperatorCreate(ceed, qfError, NULL, NULL, &opError); 241 CeedOperatorSetField(opError, "u", ceeddata->Erestrictu, 242 ceeddata->basisu, CEED_VECTOR_ACTIVE); 243 CeedOperatorSetField(opError, "true_soln", ceeddata->Erestrictui, 244 CEED_BASIS_COLLOCATED, target); 245 CeedOperatorSetField(opError, "error", ceeddata->Erestrictui, 246 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 247 248 // Set up Mat 249 userO->comm = comm; 250 userO->dm = dm; 251 userO->Xloc = Xloc; 252 ierr = VecDuplicate(Xloc, &userO->Yloc); CHKERRQ(ierr); 253 userO->Xceed = ceeddata->Xceed; 254 userO->Yceed = ceeddata->Yceed; 255 userO->op = ceeddata->opApply; 256 userO->ceed = ceed; 257 258 // Setup solver 259 ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr); 260 { 261 PC pc; 262 ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 263 if (bpChoice == CEED_BP1 || bpChoice == CEED_BP2) { 264 ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr); 265 ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr); 266 } else { 267 ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr); 268 MatNullSpace nullspace; 269 270 ierr = MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, 0, &nullspace); 271 CHKERRQ(ierr); 272 ierr = MatSetNullSpace(matO, nullspace); CHKERRQ(ierr); 273 ierr = MatNullSpaceDestroy(&nullspace); CHKERRQ(ierr); 274 } 275 ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 276 ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 277 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 278 PETSC_DEFAULT); CHKERRQ(ierr); 279 } 280 ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr); 281 ierr = KSPSetOperators(ksp, matO, matO); CHKERRQ(ierr); 282 283 // First run, if benchmarking 284 if (benchmark_mode) { 285 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1); 286 CHKERRQ(ierr); 287 my_rt_start = MPI_Wtime(); 288 ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 289 my_rt = MPI_Wtime() - my_rt_start; 290 ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm); 291 CHKERRQ(ierr); 292 // Set maxits based on first iteration timing 293 if (my_rt > 0.02) { 294 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5); 295 CHKERRQ(ierr); 296 } else { 297 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20); 298 CHKERRQ(ierr); 299 } 300 } 301 302 // Timed solve 303 ierr = VecZeroEntries(X); CHKERRQ(ierr); 304 ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr); 305 306 // -- Performance logging 307 ierr = PetscLogStageRegister("Solve Stage", &solvestage); CHKERRQ(ierr); 308 ierr = PetscLogStagePush(solvestage); CHKERRQ(ierr); 309 310 // -- Solve 311 my_rt_start = MPI_Wtime(); 312 ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 313 my_rt = MPI_Wtime() - my_rt_start; 314 315 // -- Performance logging 316 ierr = PetscLogStagePop(); 317 318 // Output results 319 { 320 KSPType ksptype; 321 KSPConvergedReason reason; 322 PetscReal rnorm; 323 PetscInt its; 324 ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr); 325 ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr); 326 ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr); 327 ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr); 328 if (!test_mode || reason < 0 || rnorm > 1e-8) { 329 ierr = PetscPrintf(comm, 330 " KSP:\n" 331 " KSP Type : %s\n" 332 " KSP Convergence : %s\n" 333 " Total KSP Iterations : %D\n" 334 " Final rnorm : %e\n", 335 ksptype, KSPConvergedReasons[reason], its, 336 (double)rnorm); CHKERRQ(ierr); 337 } 338 if (!test_mode) { 339 ierr = PetscPrintf(comm," Performance:\n"); CHKERRQ(ierr); 340 } 341 { 342 PetscReal maxerror; 343 ierr = ComputeErrorMax(userO, opError, X, target, &maxerror); 344 CHKERRQ(ierr); 345 PetscReal tol = 5e-4; 346 if (!test_mode || maxerror > tol) { 347 ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm); 348 CHKERRQ(ierr); 349 ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm); 350 CHKERRQ(ierr); 351 ierr = PetscPrintf(comm, 352 " Pointwise Error (max) : %e\n" 353 " CG Solve Time : %g (%g) sec\n", 354 (double)maxerror, rt_max, rt_min); CHKERRQ(ierr); 355 } 356 } 357 if (benchmark_mode && (!test_mode)) { 358 ierr = PetscPrintf(comm, 359 " DoFs/Sec in CG : %g (%g) million\n", 360 1e-6*gsize*its/rt_max, 1e-6*gsize*its/rt_min); CHKERRQ(ierr); 361 } 362 } 363 364 // Output solution 365 if (write_solution) { 366 PetscViewer vtkviewersoln; 367 368 ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr); 369 ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr); 370 ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtu"); CHKERRQ(ierr); 371 ierr = VecView(X, vtkviewersoln); CHKERRQ(ierr); 372 ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr); 373 } 374 375 // Cleanup 376 ierr = VecDestroy(&X); CHKERRQ(ierr); 377 ierr = VecDestroy(&Xloc); CHKERRQ(ierr); 378 ierr = VecDestroy(&userO->Yloc); CHKERRQ(ierr); 379 ierr = MatDestroy(&matO); CHKERRQ(ierr); 380 ierr = PetscFree(userO); CHKERRQ(ierr); 381 ierr = CeedDataDestroy(0, ceeddata); CHKERRQ(ierr); 382 ierr = DMDestroy(&dm); CHKERRQ(ierr); 383 384 ierr = VecDestroy(&rhs); CHKERRQ(ierr); 385 ierr = VecDestroy(&rhsloc); CHKERRQ(ierr); 386 ierr = KSPDestroy(&ksp); CHKERRQ(ierr); 387 CeedVectorDestroy(&target); 388 CeedQFunctionDestroy(&qfError); 389 CeedOperatorDestroy(&opError); 390 CeedDestroy(&ceed); 391 return PetscFinalize(); 392 } 393