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