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 // ----------------------------------------------------------------------------- 74 // Parameter structure for running problems 75 // ----------------------------------------------------------------------------- 76 typedef struct RunParams_ *RunParams; 77 struct RunParams_ { 78 MPI_Comm comm; 79 PetscBool test_mode, benchmark_mode, read_mesh, userlnodes, setmemtyperequest, 80 petschavecuda, write_solution; 81 char *filename, *ceedresource, *hostname; 82 PetscInt localnodes, degree, qextra, dim, ncompu, *melem; 83 PetscInt ksp_max_it_clip[2]; 84 PetscMPIInt rankspernode; 85 bpType bpchoice; 86 CeedMemType memtyperequested; 87 PetscLogStage solvestage; 88 }; 89 90 // ----------------------------------------------------------------------------- 91 // Main body of program, called in a loop for performance benchmarking purposes 92 // ----------------------------------------------------------------------------- 93 94 static PetscErrorCode Run(RunParams rp) { 95 PetscInt ierr; 96 double my_rt_start, my_rt, rt_min, rt_max; 97 PetscInt xlsize, lsize, gsize; 98 PetscScalar *r; 99 Vec X, Xloc, rhs, rhsloc; 100 Mat matO; 101 KSP ksp; 102 DM dm; 103 UserO userO; 104 Ceed ceed; 105 CeedData ceeddata; 106 CeedQFunction qferror; 107 CeedOperator operror; 108 CeedVector rhsceed, target; 109 110 PetscFunctionBeginUser; 111 // Setup DM 112 if (rp->read_mesh) { 113 ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, rp->filename, PETSC_TRUE, &dm); 114 CHKERRQ(ierr); 115 } else { 116 if (rp->userlnodes) { 117 // Find a nicely composite number of elements no less than global nodes 118 PetscMPIInt size; 119 ierr = MPI_Comm_size(rp->comm, &size); CHKERRQ(ierr); 120 for (PetscInt gelem = 121 PetscMax(1, size * rp->localnodes / PetscPowInt(rp->degree, rp->dim)); 122 ; 123 gelem++) { 124 Split3(gelem, rp->melem, true); 125 if (Max3(rp->melem) / Min3(rp->melem) <= 2) break; 126 } 127 } 128 ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, rp->dim, PETSC_FALSE, rp->melem, 129 NULL, NULL, NULL, PETSC_TRUE, &dm); CHKERRQ(ierr); 130 } 131 132 { 133 DM dmDist = NULL; 134 PetscPartitioner part; 135 136 ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr); 137 ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 138 ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr); 139 if (dmDist) { 140 ierr = DMDestroy(&dm); CHKERRQ(ierr); 141 dm = dmDist; 142 } 143 } 144 145 // Set up libCEED 146 CeedInit(rp->ceedresource, &ceed); 147 CeedMemType memtypebackend; 148 CeedGetPreferredMemType(ceed, &memtypebackend); 149 150 // Check memtype compatibility 151 if (!rp->setmemtyperequest) 152 rp->memtyperequested = memtypebackend; 153 else if (!rp->petschavecuda && rp->memtyperequested == CEED_MEM_DEVICE) 154 SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, 155 "PETSc was not built with CUDA. " 156 "Requested MemType CEED_MEM_DEVICE is not supported.", NULL); 157 158 // Create DM 159 ierr = SetupDMByDegree(dm, rp->degree, rp->ncompu, rp->bpchoice); 160 CHKERRQ(ierr); 161 162 // Create vectors 163 if (rp->memtyperequested == CEED_MEM_DEVICE) { 164 ierr = DMSetVecType(dm, VECCUDA); CHKERRQ(ierr); 165 } 166 ierr = DMCreateGlobalVector(dm, &X); CHKERRQ(ierr); 167 ierr = VecGetLocalSize(X, &lsize); CHKERRQ(ierr); 168 ierr = VecGetSize(X, &gsize); CHKERRQ(ierr); 169 ierr = DMCreateLocalVector(dm, &Xloc); CHKERRQ(ierr); 170 ierr = VecGetSize(Xloc, &xlsize); CHKERRQ(ierr); 171 ierr = VecDuplicate(X, &rhs); CHKERRQ(ierr); 172 173 // Operator 174 ierr = PetscMalloc1(1, &userO); CHKERRQ(ierr); 175 ierr = MatCreateShell(rp->comm, lsize, lsize, gsize, gsize, 176 userO, &matO); CHKERRQ(ierr); 177 ierr = MatShellSetOperation(matO, MATOP_MULT, 178 (void(*)(void))MatMult_Ceed); CHKERRQ(ierr); 179 ierr = MatShellSetOperation(matO, MATOP_GET_DIAGONAL, 180 (void(*)(void))MatGetDiag); CHKERRQ(ierr); 181 if (rp->memtyperequested == CEED_MEM_DEVICE) { 182 ierr = MatShellSetVecType(matO, VECCUDA); CHKERRQ(ierr); 183 } 184 185 // Print summary 186 if (!rp->test_mode) { 187 PetscInt P = rp->degree + 1, Q = P + rp->qextra; 188 189 const char *usedresource; 190 CeedGetResource(ceed, &usedresource); 191 192 VecType vectype; 193 ierr = VecGetType(X, &vectype); CHKERRQ(ierr); 194 195 PetscInt cStart, cEnd; 196 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); CHKERRQ(ierr); 197 PetscMPIInt comm_size; 198 ierr = MPI_Comm_size(rp->comm, &comm_size); CHKERRQ(ierr); 199 ierr = PetscPrintf(rp->comm, 200 "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n" 201 " MPI:\n" 202 " Hostname : %s\n" 203 " Total ranks : %d\n" 204 " Ranks per compute node : %d\n" 205 " PETSc:\n" 206 " PETSc Vec Type : %s\n" 207 " libCEED:\n" 208 " libCEED Backend : %s\n" 209 " libCEED Backend MemType : %s\n" 210 " libCEED User Requested MemType : %s\n" 211 " Mesh:\n" 212 " Number of 1D Basis Nodes (P) : %d\n" 213 " Number of 1D Quadrature Points (Q) : %d\n" 214 " Global nodes : %D\n" 215 " Local Elements : %D\n" 216 " Owned nodes : %D\n" 217 " DoF per node : %D\n", 218 rp->bpchoice+1, 219 rp->hostname, 220 comm_size, rp->rankspernode, 221 vectype, usedresource, 222 CeedMemTypes[memtypebackend], 223 (rp->setmemtyperequest) ? 224 CeedMemTypes[rp->memtyperequested] : "none", 225 P, Q, gsize/rp->ncompu, cEnd - cStart, lsize/rp->ncompu, 226 rp->ncompu); 227 CHKERRQ(ierr); 228 } 229 230 // Create RHS vector 231 ierr = VecDuplicate(Xloc, &rhsloc); CHKERRQ(ierr); 232 ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr); 233 if (rp->memtyperequested == CEED_MEM_HOST) { 234 ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr); 235 } else { 236 ierr = VecCUDAGetArray(rhsloc, &r); CHKERRQ(ierr); 237 } 238 CeedVectorCreate(ceed, xlsize, &rhsceed); 239 CeedVectorSetArray(rhsceed, rp->memtyperequested, CEED_USE_POINTER, r); 240 241 ierr = PetscMalloc1(1, &ceeddata); CHKERRQ(ierr); 242 ierr = SetupLibceedByDegree(dm, ceed, rp->degree, rp->dim, rp->qextra, 243 rp->ncompu, gsize, xlsize, rp->bpchoice, ceeddata, 244 true, rhsceed, &target); CHKERRQ(ierr); 245 246 // Gather RHS 247 CeedVectorSyncArray(rhsceed, rp->memtyperequested); 248 if (rp->memtyperequested == CEED_MEM_HOST) { 249 ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr); 250 } else { 251 ierr = VecCUDARestoreArray(rhsloc, &r); CHKERRQ(ierr); 252 } 253 ierr = VecZeroEntries(rhs); CHKERRQ(ierr); 254 ierr = DMLocalToGlobalBegin(dm, rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr); 255 ierr = DMLocalToGlobalEnd(dm, rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr); 256 CeedVectorDestroy(&rhsceed); 257 258 // Create the error Q-function 259 CeedQFunctionCreateInterior(ceed, 1, bpOptions[rp->bpchoice].error, 260 bpOptions[rp->bpchoice].errorfname, &qferror); 261 CeedQFunctionAddInput(qferror, "u", rp->ncompu, CEED_EVAL_INTERP); 262 CeedQFunctionAddInput(qferror, "true_soln", rp->ncompu, CEED_EVAL_NONE); 263 CeedQFunctionAddOutput(qferror, "error", rp->ncompu, CEED_EVAL_NONE); 264 265 // Create the error operator 266 CeedOperatorCreate(ceed, qferror, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, 267 &operror); 268 CeedOperatorSetField(operror, "u", ceeddata->Erestrictu, 269 ceeddata->basisu, CEED_VECTOR_ACTIVE); 270 CeedOperatorSetField(operror, "true_soln", ceeddata->Erestrictui, 271 CEED_BASIS_COLLOCATED, target); 272 CeedOperatorSetField(operror, "error", ceeddata->Erestrictui, 273 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 274 275 // Set up Mat 276 userO->comm = rp->comm; 277 userO->dm = dm; 278 userO->Xloc = Xloc; 279 ierr = VecDuplicate(Xloc, &userO->Yloc); CHKERRQ(ierr); 280 userO->xceed = ceeddata->xceed; 281 userO->yceed = ceeddata->yceed; 282 userO->op = ceeddata->opapply; 283 userO->ceed = ceed; 284 userO->memtype = rp->memtyperequested; 285 if (rp->memtyperequested == CEED_MEM_HOST) { 286 userO->VecGetArray = VecGetArray; 287 userO->VecGetArrayRead = VecGetArrayRead; 288 userO->VecRestoreArray = VecRestoreArray; 289 userO->VecRestoreArrayRead = VecRestoreArrayRead; 290 } else { 291 userO->VecGetArray = VecCUDAGetArray; 292 userO->VecGetArrayRead = VecCUDAGetArrayRead; 293 userO->VecRestoreArray = VecCUDARestoreArray; 294 userO->VecRestoreArrayRead = VecCUDARestoreArrayRead; 295 } 296 297 ierr = KSPCreate(rp->comm, &ksp); CHKERRQ(ierr); 298 { 299 PC pc; 300 ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 301 if (rp->bpchoice == CEED_BP1 || rp->bpchoice == CEED_BP2) { 302 ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr); 303 ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr); 304 } else { 305 ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr); 306 } 307 ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 308 ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 309 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 310 PETSC_DEFAULT); CHKERRQ(ierr); 311 } 312 ierr = KSPSetOperators(ksp, matO, matO); CHKERRQ(ierr); 313 314 // First run, if benchmarking 315 if (rp->benchmark_mode) { 316 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1); 317 CHKERRQ(ierr); 318 my_rt_start = MPI_Wtime(); 319 ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 320 my_rt = MPI_Wtime() - my_rt_start; 321 ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, rp->comm); 322 CHKERRQ(ierr); 323 // Set maxits based on first iteration timing 324 if (my_rt > 0.02) { 325 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, rp->ksp_max_it_clip[0]); 326 CHKERRQ(ierr); 327 } else { 328 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, rp->ksp_max_it_clip[1]); 329 CHKERRQ(ierr); 330 } 331 } 332 ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr); 333 334 // Timed solve 335 ierr = VecZeroEntries(X); CHKERRQ(ierr); 336 ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr); 337 338 // -- Performance logging 339 ierr = PetscLogStagePush(rp->solvestage); CHKERRQ(ierr); 340 341 // -- Solve 342 my_rt_start = MPI_Wtime(); 343 ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 344 my_rt = MPI_Wtime() - my_rt_start; 345 346 // -- Performance logging 347 ierr = PetscLogStagePop(); 348 349 // Output results 350 { 351 KSPType ksptype; 352 KSPConvergedReason reason; 353 PetscReal rnorm; 354 PetscInt its; 355 ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr); 356 ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr); 357 ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr); 358 ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr); 359 if (!rp->test_mode || reason < 0 || rnorm > 1e-8) { 360 ierr = PetscPrintf(rp->comm, 361 " KSP:\n" 362 " KSP Type : %s\n" 363 " KSP Convergence : %s\n" 364 " Total KSP Iterations : %D\n" 365 " Final rnorm : %e\n", 366 ksptype, KSPConvergedReasons[reason], its, 367 (double)rnorm); CHKERRQ(ierr); 368 } 369 if (!rp->test_mode) { 370 ierr = PetscPrintf(rp->comm," Performance:\n"); CHKERRQ(ierr); 371 } 372 { 373 PetscReal maxerror; 374 ierr = ComputeErrorMax(userO, operror, X, target, &maxerror); 375 CHKERRQ(ierr); 376 PetscReal tol = 5e-2; 377 if (!rp->test_mode || maxerror > tol) { 378 ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, rp->comm); 379 CHKERRQ(ierr); 380 ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, rp->comm); 381 CHKERRQ(ierr); 382 ierr = PetscPrintf(rp->comm, 383 " Pointwise Error (max) : %e\n" 384 " CG Solve Time : %g (%g) sec\n", 385 (double)maxerror, rt_max, rt_min); CHKERRQ(ierr); 386 } 387 } 388 if (rp->benchmark_mode && (!rp->test_mode)) { 389 ierr = PetscPrintf(rp->comm, 390 " DoFs/Sec in CG : %g (%g) million\n", 391 1e-6*gsize*its/rt_max, 392 1e-6*gsize*its/rt_min); CHKERRQ(ierr); 393 } 394 } 395 396 if (rp->write_solution) { 397 PetscViewer vtkviewersoln; 398 399 ierr = PetscViewerCreate(rp->comm, &vtkviewersoln); CHKERRQ(ierr); 400 ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr); 401 ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr); 402 ierr = VecView(X, vtkviewersoln); CHKERRQ(ierr); 403 ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr); 404 } 405 406 // Cleanup 407 ierr = VecDestroy(&X); CHKERRQ(ierr); 408 ierr = VecDestroy(&Xloc); CHKERRQ(ierr); 409 ierr = VecDestroy(&userO->Yloc); CHKERRQ(ierr); 410 ierr = MatDestroy(&matO); CHKERRQ(ierr); 411 ierr = PetscFree(userO); CHKERRQ(ierr); 412 ierr = CeedDataDestroy(0, ceeddata); CHKERRQ(ierr); 413 ierr = DMDestroy(&dm); CHKERRQ(ierr); 414 415 ierr = VecDestroy(&rhs); CHKERRQ(ierr); 416 ierr = VecDestroy(&rhsloc); CHKERRQ(ierr); 417 ierr = KSPDestroy(&ksp); CHKERRQ(ierr); 418 CeedVectorDestroy(&target); 419 CeedQFunctionDestroy(&qferror); 420 CeedOperatorDestroy(&operror); 421 CeedDestroy(&ceed); 422 PetscFunctionReturn(0); 423 } 424 425 int main(int argc, char **argv) { 426 PetscInt ierr, commsize; 427 RunParams rp; 428 MPI_Comm comm; 429 char filename[PETSC_MAX_PATH_LEN]; 430 char ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self"; 431 char hostname[PETSC_MAX_PATH_LEN]; 432 433 PetscInt dim = 3, melem[3] = {3, 3, 3}; 434 PetscInt num_degrees = 30, degree[30] = {}, num_localnodes = 2, localnodes[2] = {}; 435 PetscMPIInt rankspernode; 436 PetscBool degree_set; 437 bpType bpchoice; 438 439 // Check PETSc CUDA support 440 PetscBool petschavecuda; 441 // *INDENT-OFF* 442 #ifdef PETSC_HAVE_CUDA 443 petschavecuda = PETSC_TRUE; 444 #else 445 petschavecuda = PETSC_FALSE; 446 #endif 447 // *INDENT-ON* 448 449 // Initialize PETSc 450 ierr = PetscInitialize(&argc, &argv, NULL, help); 451 if (ierr) return ierr; 452 comm = PETSC_COMM_WORLD; 453 ierr = MPI_Comm_size(comm, &commsize); 454 if (ierr != MPI_SUCCESS) return ierr; 455 #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) 456 { 457 MPI_Comm splitcomm; 458 ierr = MPI_Comm_split_type(comm, MPI_COMM_TYPE_SHARED, 0, MPI_INFO_NULL, 459 &splitcomm); 460 CHKERRQ(ierr); 461 ierr = MPI_Comm_size(splitcomm, &rankspernode); CHKERRQ(ierr); 462 ierr = MPI_Comm_free(&splitcomm); CHKERRQ(ierr); 463 } 464 #else 465 rankspernode = -1; // Unknown 466 #endif 467 468 // Setup all parameters needed in Run() 469 ierr = PetscMalloc1(1, &rp); CHKERRQ(ierr); 470 rp->comm = comm; 471 472 // Read command line options 473 ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); 474 CHKERRQ(ierr); 475 bpchoice = CEED_BP1; 476 ierr = PetscOptionsEnum("-problem", "CEED benchmark problem to solve", NULL, 477 bpTypes, (PetscEnum)bpchoice, (PetscEnum *)&bpchoice, 478 NULL); CHKERRQ(ierr); 479 rp->ncompu = bpOptions[bpchoice].ncompu; 480 rp->test_mode = PETSC_FALSE; 481 ierr = PetscOptionsBool("-test", 482 "Testing mode (do not print unless error is large)", 483 NULL, rp->test_mode, &rp->test_mode, NULL); CHKERRQ(ierr); 484 rp->benchmark_mode = PETSC_FALSE; 485 ierr = PetscOptionsBool("-benchmark", 486 "Benchmarking mode (prints benchmark statistics)", 487 NULL, rp->benchmark_mode, &rp->benchmark_mode, NULL); 488 CHKERRQ(ierr); 489 rp->write_solution = PETSC_FALSE; 490 ierr = PetscOptionsBool("-write_solution", "Write solution for visualization", 491 NULL, rp->write_solution, &rp->write_solution, NULL); 492 CHKERRQ(ierr); 493 degree[0] = rp->test_mode ? 3 : 2; 494 ierr = PetscOptionsIntArray("-degree", 495 "Polynomial degree of tensor product basis", NULL, 496 degree, &num_degrees, °ree_set); CHKERRQ(ierr); 497 if (!degree_set) 498 num_degrees = 1; 499 rp->qextra = bpOptions[bpchoice].qextra; 500 ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", NULL, 501 rp->qextra, &rp->qextra, NULL); CHKERRQ(ierr); 502 ierr = PetscOptionsString("-ceed", "CEED resource specifier", NULL, 503 ceedresource, 504 ceedresource, sizeof(ceedresource), NULL); CHKERRQ(ierr); 505 ierr = PetscGetHostName(hostname, sizeof hostname); CHKERRQ(ierr); 506 ierr = PetscOptionsString("-hostname", "Hostname for output", NULL, hostname, 507 hostname, sizeof(hostname), NULL); CHKERRQ(ierr); 508 rp->read_mesh = PETSC_FALSE; 509 ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL, filename, 510 filename, sizeof(filename), &rp->read_mesh); 511 CHKERRQ(ierr); 512 rp->filename = filename; 513 if (!rp->read_mesh) { 514 PetscInt tmp = dim; 515 ierr = PetscOptionsIntArray("-cells", "Number of cells per dimension", NULL, 516 melem, &tmp, NULL); CHKERRQ(ierr); 517 } 518 rp->memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST; 519 ierr = PetscOptionsEnum("-memtype", "CEED MemType requested", NULL, memTypes, 520 (PetscEnum)rp->memtyperequested, 521 (PetscEnum *)&rp->memtyperequested, &rp->setmemtyperequest); 522 CHKERRQ(ierr); 523 524 localnodes[0] = 1000; 525 ierr = PetscOptionsIntArray("-local_nodes", 526 "Target number of locally owned nodes per " 527 "process (single value or min,max)", 528 NULL, localnodes, &num_localnodes, &rp->userlnodes); 529 CHKERRQ(ierr); 530 if (num_localnodes < 2) 531 localnodes[1] = 2 * localnodes[0]; 532 { 533 PetscInt two = 2; 534 rp->ksp_max_it_clip[0] = 5; 535 rp->ksp_max_it_clip[1] = 20; 536 ierr = PetscOptionsIntArray("-ksp_max_it_clip", 537 "Min and max number of iterations to use during benchmarking", 538 NULL, rp->ksp_max_it_clip, &two, NULL); CHKERRQ(ierr); 539 } 540 541 if (rp->benchmark_mode && !degree_set) { 542 PetscInt maxdegree = 8; 543 ierr = PetscOptionsInt("-max_degree", "Max degree to run in benchmark mode", 544 NULL, maxdegree, &maxdegree, NULL); 545 CHKERRQ(ierr); 546 for (PetscInt i = 0; i < maxdegree; i++) 547 degree[i] = i + 1; 548 num_degrees = maxdegree; 549 } 550 { 551 PetscBool flg; 552 PetscInt p = rankspernode; 553 ierr = PetscOptionsInt("-p", "Number of MPI ranks per node", NULL, 554 p, &p, &flg); 555 CHKERRQ(ierr); 556 if (flg) rankspernode = p; 557 } 558 559 ierr = PetscOptionsEnd(); 560 CHKERRQ(ierr); 561 562 // Register PETSc logging stage 563 ierr = PetscLogStageRegister("Solve Stage", &rp->solvestage); 564 CHKERRQ(ierr); 565 566 rp->petschavecuda = petschavecuda; 567 rp->ceedresource = ceedresource; 568 rp->hostname = hostname; 569 rp->dim = dim; 570 rp->melem = melem; 571 rp->bpchoice = bpchoice; 572 rp->rankspernode = rankspernode; 573 574 for (PetscInt d = 0; d < num_degrees; d++) { 575 PetscInt deg = degree[d]; 576 for (PetscInt n = localnodes[0]; n < localnodes[1]; n *= 2) { 577 rp->degree = deg; 578 rp->localnodes = n; 579 ierr = Run(rp); CHKERRQ(ierr); 580 } 581 } 582 583 // Clear memory 584 ierr = PetscFree(rp); CHKERRQ(ierr); 585 586 return PetscFinalize(); 587 } 588