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