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 // 22 // The code uses higher level communication protocols in DMPlex. 23 // 24 // Build with: 25 // 26 // make bps [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] 27 // 28 // Sample runs: 29 // 30 // ./bps -problem bp1 -degree 3 31 // ./bps -problem bp2 -ceed /cpu/self -degree 3 32 // ./bps -problem bp3 -ceed /gpu/occa -degree 3 33 // ./bps -problem bp4 -ceed /cpu/occa -degree 3 34 // ./bps -problem bp5 -ceed /omp/occa -degree 3 35 // ./bps -problem bp6 -ceed /ocl/occa -degree 3 36 // 37 //TESTARGS -ceed {ceed_resource} -test -problem bp5 -degree 3 38 39 /// @file 40 /// CEED BPs example using PETSc with DMPlex 41 /// See bpsraw.c for a "raw" implementation using a structured grid. 42 const char help[] = "Solve CEED BPs using PETSc with DMPlex\n"; 43 44 #include <stdbool.h> 45 #include <string.h> 46 #include <petscksp.h> 47 #include <petscdmplex.h> 48 #include <ceed.h> 49 #include "setup.h" 50 51 // ----------------------------------------------------------------------------- 52 // Utilities 53 // ----------------------------------------------------------------------------- 54 55 // Utility function, compute three factors of an integer 56 static void Split3(PetscInt size, PetscInt m[3], bool reverse) { 57 for (PetscInt d=0,sizeleft=size; d<3; d++) { 58 PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(sizeleft, 1./(3 - d))); 59 while (try * (sizeleft / try) != sizeleft) try++; 60 m[reverse ? 2-d : d] = try; 61 sizeleft /= try; 62 } 63 } 64 65 static int Max3(const PetscInt a[3]) { 66 return PetscMax(a[0], PetscMax(a[1], a[2])); 67 } 68 69 static int Min3(const PetscInt a[3]) { 70 return PetscMin(a[0], PetscMin(a[1], a[2])); 71 } 72 73 int main(int argc, char **argv) { 74 PetscInt ierr; 75 MPI_Comm comm; 76 char filename[PETSC_MAX_PATH_LEN], 77 ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self"; 78 double my_rt_start, my_rt, rt_min, rt_max; 79 PetscInt degree = 3, qextra, lsize, gsize, dim = 3, melem[3] = {3, 3, 3}, 80 ncompu = 1, xlsize, localnodes; 81 PetscScalar *r; 82 PetscBool test_mode, benchmark_mode, read_mesh, write_solution, 83 userlnodes = PETSC_FALSE; 84 PetscLogStage solvestage; 85 Vec X, Xloc, rhs, rhsloc; 86 Mat matO; 87 KSP ksp; 88 DM dm; 89 UserO userO; 90 Ceed ceed; 91 CeedData ceeddata; 92 CeedQFunction qferror; 93 CeedOperator operror; 94 CeedVector rhsceed, target; 95 CeedMemType memtyperequested; 96 bpType bpchoice; 97 98 // Check PETSc CUDA support 99 PetscBool petschavecuda, setmemtyperequest = PETSC_FALSE; 100 // *INDENT-OFF* 101 #ifdef PETSC_HAVE_CUDA 102 petschavecuda = PETSC_TRUE; 103 #else 104 petschavecuda = PETSC_FALSE; 105 #endif 106 // *INDENT-ON* 107 108 ierr = PetscInitialize(&argc, &argv, NULL, help); 109 if (ierr) return ierr; 110 comm = PETSC_COMM_WORLD; 111 112 // Read command line options 113 ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr); 114 bpchoice = CEED_BP1; 115 ierr = PetscOptionsEnum("-problem", 116 "CEED benchmark problem to solve", NULL, 117 bpTypes, (PetscEnum)bpchoice, (PetscEnum *)&bpchoice, 118 NULL); CHKERRQ(ierr); 119 ncompu = bpOptions[bpchoice].ncompu; 120 test_mode = PETSC_FALSE; 121 ierr = PetscOptionsBool("-test", 122 "Testing mode (do not print unless error is large)", 123 NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr); 124 benchmark_mode = PETSC_FALSE; 125 ierr = PetscOptionsBool("-benchmark", 126 "Benchmarking mode (prints benchmark statistics)", 127 NULL, benchmark_mode, &benchmark_mode, NULL); 128 CHKERRQ(ierr); 129 write_solution = PETSC_FALSE; 130 ierr = PetscOptionsBool("-write_solution", 131 "Write solution for visualization", 132 NULL, write_solution, &write_solution, NULL); 133 CHKERRQ(ierr); 134 degree = test_mode ? 3 : 2; 135 ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", 136 NULL, degree, °ree, NULL); CHKERRQ(ierr); 137 qextra = bpOptions[bpchoice].qextra; 138 ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", 139 NULL, qextra, &qextra, NULL); CHKERRQ(ierr); 140 ierr = PetscOptionsString("-ceed", "CEED resource specifier", 141 NULL, ceedresource, ceedresource, 142 sizeof(ceedresource), NULL); CHKERRQ(ierr); 143 read_mesh = PETSC_FALSE; 144 ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL, 145 filename, filename, sizeof(filename), &read_mesh); 146 CHKERRQ(ierr); 147 if (!read_mesh) { 148 PetscInt tmp = dim; 149 ierr = PetscOptionsIntArray("-cells", "Number of cells per dimension", NULL, 150 melem, &tmp, NULL); CHKERRQ(ierr); 151 } 152 memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST; 153 ierr = PetscOptionsEnum("-memtype", 154 "CEED MemType requested", NULL, 155 memTypes, (PetscEnum)memtyperequested, 156 (PetscEnum *)&memtyperequested, &setmemtyperequest); 157 CHKERRQ(ierr); 158 localnodes = 1000; 159 ierr = PetscOptionsInt("-local_nodes", 160 "Target number of locally owned nodes per process", 161 NULL, localnodes, &localnodes, &userlnodes); 162 CHKERRQ(ierr); 163 164 ierr = PetscOptionsEnd(); CHKERRQ(ierr); 165 166 // Setup DM 167 if (read_mesh) { 168 ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dm); 169 CHKERRQ(ierr); 170 } else { 171 if (userlnodes) { 172 // Find a nicely composite number of elements no less than global nodes 173 PetscMPIInt size; 174 ierr = MPI_Comm_size(comm, &size); CHKERRQ(ierr); 175 for (PetscInt gelem = PetscMax(1, size * localnodes / PetscPowInt(degree, dim)); 176 ; 177 gelem++) { 178 Split3(gelem, melem, true); 179 if (Max3(melem) / Min3(melem) <= 2) break; 180 } 181 } 182 ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_FALSE, melem, NULL, 183 NULL, NULL, PETSC_TRUE, &dm); CHKERRQ(ierr); 184 } 185 186 { 187 DM dmDist = NULL; 188 PetscPartitioner part; 189 190 ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr); 191 ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 192 ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr); 193 if (dmDist) { 194 ierr = DMDestroy(&dm); CHKERRQ(ierr); 195 dm = dmDist; 196 } 197 } 198 199 // Set up libCEED 200 CeedInit(ceedresource, &ceed); 201 CeedMemType memtypebackend; 202 CeedGetPreferredMemType(ceed, &memtypebackend); 203 204 // Check memtype compatibility 205 if (!setmemtyperequest) 206 memtyperequested = memtypebackend; 207 else if (!petschavecuda && memtyperequested == CEED_MEM_DEVICE) 208 SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, 209 "PETSc was not built with CUDA. " 210 "Requested MemType CEED_MEM_DEVICE is not supported.", NULL); 211 212 // Create DM 213 ierr = SetupDMByDegree(dm, degree, ncompu, bpchoice); 214 CHKERRQ(ierr); 215 216 // Create vectors 217 if (memtyperequested == CEED_MEM_DEVICE) { 218 ierr = DMSetVecType(dm, VECCUDA); CHKERRQ(ierr); 219 } 220 ierr = DMCreateGlobalVector(dm, &X); CHKERRQ(ierr); 221 ierr = VecGetLocalSize(X, &lsize); CHKERRQ(ierr); 222 ierr = VecGetSize(X, &gsize); CHKERRQ(ierr); 223 ierr = DMCreateLocalVector(dm, &Xloc); CHKERRQ(ierr); 224 ierr = VecGetSize(Xloc, &xlsize); CHKERRQ(ierr); 225 ierr = VecDuplicate(X, &rhs); CHKERRQ(ierr); 226 227 // Operator 228 ierr = PetscMalloc1(1, &userO); CHKERRQ(ierr); 229 ierr = MatCreateShell(comm, lsize, lsize, gsize, gsize, 230 userO, &matO); CHKERRQ(ierr); 231 ierr = MatShellSetOperation(matO, MATOP_MULT, 232 (void(*)(void))MatMult_Ceed); CHKERRQ(ierr); 233 ierr = MatShellSetOperation(matO, MATOP_GET_DIAGONAL, 234 (void(*)(void))MatGetDiag); CHKERRQ(ierr); 235 if (memtyperequested == CEED_MEM_DEVICE) { 236 ierr = MatShellSetVecType(matO, VECCUDA); CHKERRQ(ierr); 237 } 238 239 // Print summary 240 if (!test_mode) { 241 PetscInt P = degree + 1, Q = P + qextra; 242 243 const char *usedresource; 244 CeedGetResource(ceed, &usedresource); 245 246 VecType vectype; 247 ierr = VecGetType(X, &vectype); CHKERRQ(ierr); 248 249 PetscInt cStart, cEnd; 250 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); CHKERRQ(ierr); 251 ierr = PetscPrintf(comm, 252 "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n" 253 " PETSc:\n" 254 " PETSc Vec Type : %s\n" 255 " libCEED:\n" 256 " libCEED Backend : %s\n" 257 " libCEED Backend MemType : %s\n" 258 " libCEED User Requested MemType : %s\n" 259 " Mesh:\n" 260 " Number of 1D Basis Nodes (p) : %d\n" 261 " Number of 1D Quadrature Points (q) : %d\n" 262 " Global nodes : %D\n" 263 " Local Elements : %D\n" 264 " Owned nodes : %D\n" 265 " DoF per node : %D\n", 266 bpchoice+1, vectype, usedresource, 267 CeedMemTypes[memtypebackend], 268 (setmemtyperequest) ? 269 CeedMemTypes[memtyperequested] : "none", 270 P, Q, gsize/ncompu, cEnd - cStart, lsize/ncompu, 271 ncompu); 272 CHKERRQ(ierr); 273 } 274 275 // Create RHS vector 276 ierr = VecDuplicate(Xloc, &rhsloc); CHKERRQ(ierr); 277 ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr); 278 if (memtyperequested == CEED_MEM_HOST) { 279 ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr); 280 } else { 281 ierr = VecCUDAGetArray(rhsloc, &r); CHKERRQ(ierr); 282 } 283 CeedVectorCreate(ceed, xlsize, &rhsceed); 284 CeedVectorSetArray(rhsceed, memtyperequested, CEED_USE_POINTER, r); 285 286 ierr = PetscMalloc1(1, &ceeddata); CHKERRQ(ierr); 287 ierr = SetupLibceedByDegree(dm, ceed, degree, dim, qextra, 288 ncompu, gsize, xlsize, bpchoice, ceeddata, 289 true, rhsceed, &target); CHKERRQ(ierr); 290 291 // Gather RHS 292 CeedVectorSyncArray(rhsceed, memtyperequested); 293 if (memtyperequested == CEED_MEM_HOST) { 294 ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr); 295 } else { 296 ierr = VecCUDARestoreArray(rhsloc, &r); CHKERRQ(ierr); 297 } 298 ierr = VecZeroEntries(rhs); CHKERRQ(ierr); 299 ierr = DMLocalToGlobalBegin(dm, rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr); 300 ierr = DMLocalToGlobalEnd(dm, rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr); 301 CeedVectorDestroy(&rhsceed); 302 303 // Create the error Q-function 304 CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].error, 305 bpOptions[bpchoice].errorfname, &qferror); 306 CeedQFunctionAddInput(qferror, "u", ncompu, CEED_EVAL_INTERP); 307 CeedQFunctionAddInput(qferror, "true_soln", ncompu, CEED_EVAL_NONE); 308 CeedQFunctionAddOutput(qferror, "error", ncompu, CEED_EVAL_NONE); 309 310 // Create the error operator 311 CeedOperatorCreate(ceed, qferror, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, 312 &operror); 313 CeedOperatorSetField(operror, "u", ceeddata->Erestrictu, 314 ceeddata->basisu, CEED_VECTOR_ACTIVE); 315 CeedOperatorSetField(operror, "true_soln", ceeddata->Erestrictui, 316 CEED_BASIS_COLLOCATED, target); 317 CeedOperatorSetField(operror, "error", ceeddata->Erestrictui, 318 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 319 320 // Set up Mat 321 userO->comm = comm; 322 userO->dm = dm; 323 userO->Xloc = Xloc; 324 ierr = VecDuplicate(Xloc, &userO->Yloc); CHKERRQ(ierr); 325 userO->xceed = ceeddata->xceed; 326 userO->yceed = ceeddata->yceed; 327 userO->op = ceeddata->opapply; 328 userO->ceed = ceed; 329 userO->memtype = memtyperequested; 330 if (memtyperequested == CEED_MEM_HOST) { 331 userO->VecGetArray = VecGetArray; 332 userO->VecGetArrayRead = VecGetArrayRead; 333 userO->VecRestoreArray = VecRestoreArray; 334 userO->VecRestoreArrayRead = VecRestoreArrayRead; 335 } else { 336 userO->VecGetArray = VecCUDAGetArray; 337 userO->VecGetArrayRead = VecCUDAGetArrayRead; 338 userO->VecRestoreArray = VecCUDARestoreArray; 339 userO->VecRestoreArrayRead = VecCUDARestoreArrayRead; 340 } 341 342 ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr); 343 { 344 PC pc; 345 ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 346 if (bpchoice == CEED_BP1 || bpchoice == CEED_BP2) { 347 ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr); 348 ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr); 349 } else { 350 ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr); 351 } 352 ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 353 ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 354 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 355 PETSC_DEFAULT); CHKERRQ(ierr); 356 } 357 ierr = KSPSetOperators(ksp, matO, matO); CHKERRQ(ierr); 358 359 // First run, if benchmarking 360 if (benchmark_mode) { 361 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1); 362 CHKERRQ(ierr); 363 my_rt_start = MPI_Wtime(); 364 ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 365 my_rt = MPI_Wtime() - my_rt_start; 366 ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm); 367 CHKERRQ(ierr); 368 // Set maxits based on first iteration timing 369 if (my_rt > 0.02) { 370 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5); 371 CHKERRQ(ierr); 372 } else { 373 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20); 374 CHKERRQ(ierr); 375 } 376 } 377 ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr); 378 379 // Timed solve 380 ierr = VecZeroEntries(X); CHKERRQ(ierr); 381 ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr); 382 383 // -- Performance logging 384 ierr = PetscLogStageRegister("Solve Stage", &solvestage); CHKERRQ(ierr); 385 ierr = PetscLogStagePush(solvestage); CHKERRQ(ierr); 386 387 // -- Solve 388 my_rt_start = MPI_Wtime(); 389 ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 390 my_rt = MPI_Wtime() - my_rt_start; 391 392 // -- Performance logging 393 ierr = PetscLogStagePop(); 394 395 // Output results 396 { 397 KSPType ksptype; 398 KSPConvergedReason reason; 399 PetscReal rnorm; 400 PetscInt its; 401 ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr); 402 ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr); 403 ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr); 404 ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr); 405 if (!test_mode || reason < 0 || rnorm > 1e-8) { 406 ierr = PetscPrintf(comm, 407 " KSP:\n" 408 " KSP Type : %s\n" 409 " KSP Convergence : %s\n" 410 " Total KSP Iterations : %D\n" 411 " Final rnorm : %e\n", 412 ksptype, KSPConvergedReasons[reason], its, 413 (double)rnorm); CHKERRQ(ierr); 414 } 415 if (!test_mode) { 416 ierr = PetscPrintf(comm," Performance:\n"); CHKERRQ(ierr); 417 } 418 { 419 PetscReal maxerror; 420 ierr = ComputeErrorMax(userO, operror, X, target, &maxerror); 421 CHKERRQ(ierr); 422 PetscReal tol = 5e-2; 423 if (!test_mode || maxerror > tol) { 424 ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm); 425 CHKERRQ(ierr); 426 ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm); 427 CHKERRQ(ierr); 428 ierr = PetscPrintf(comm, 429 " Pointwise Error (max) : %e\n" 430 " CG Solve Time : %g (%g) sec\n", 431 (double)maxerror, rt_max, rt_min); CHKERRQ(ierr); 432 } 433 } 434 if (benchmark_mode && (!test_mode)) { 435 ierr = PetscPrintf(comm, 436 " DoFs/Sec in CG : %g (%g) million\n", 437 1e-6*gsize*its/rt_max, 438 1e-6*gsize*its/rt_min); CHKERRQ(ierr); 439 } 440 } 441 442 if (write_solution) { 443 PetscViewer vtkviewersoln; 444 445 ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr); 446 ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr); 447 ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr); 448 ierr = VecView(X, vtkviewersoln); CHKERRQ(ierr); 449 ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr); 450 } 451 452 // Cleanup 453 ierr = VecDestroy(&X); CHKERRQ(ierr); 454 ierr = VecDestroy(&Xloc); CHKERRQ(ierr); 455 ierr = VecDestroy(&userO->Yloc); CHKERRQ(ierr); 456 ierr = MatDestroy(&matO); CHKERRQ(ierr); 457 ierr = PetscFree(userO); CHKERRQ(ierr); 458 ierr = CeedDataDestroy(0, ceeddata); CHKERRQ(ierr); 459 ierr = DMDestroy(&dm); CHKERRQ(ierr); 460 461 ierr = VecDestroy(&rhs); CHKERRQ(ierr); 462 ierr = VecDestroy(&rhsloc); CHKERRQ(ierr); 463 ierr = KSPDestroy(&ksp); CHKERRQ(ierr); 464 CeedVectorDestroy(&target); 465 CeedQFunctionDestroy(&qferror); 466 CeedOperatorDestroy(&operror); 467 CeedDestroy(&ceed); 468 return PetscFinalize(); 469 } 470