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 3-6 with Multigrid 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 multigrid [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] 27 // 28 // Sample runs: 29 // 30 // multigrid -problem bp3 31 // multigrid -problem bp4 32 // multigrid -problem bp5 -ceed /cpu/self 33 // multigrid -problem bp6 -ceed /gpu/cuda 34 // 35 //TESTARGS -ceed {ceed_resource} -test -problem bp3 -degree 3 36 37 /// @file 38 /// CEED BPs 1-6 multigrid example using PETSc 39 const char help[] = "Solve CEED BPs using p-multigrid with PETSc and DMPlex\n"; 40 41 #define multigrid 42 #include "setup.h" 43 44 int main(int argc, char **argv) { 45 PetscInt ierr; 46 MPI_Comm comm; 47 char filename[PETSC_MAX_PATH_LEN], 48 ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self"; 49 double my_rt_start, my_rt, rt_min, rt_max; 50 PetscInt degree = 3, qextra, *lsize, *xlsize, *gsize, dim = 3, fineLevel, 51 melem[3] = {3, 3, 3}, ncompu = 1, numlevels = degree, *leveldegrees; 52 PetscScalar *r; 53 PetscBool test_mode, benchmark_mode, read_mesh, write_solution; 54 PetscLogStage solvestage; 55 DM *dm, dmorig; 56 SNES snesdummy; 57 KSP ksp; 58 PC pc; 59 Mat *matO, *matPR, matcoarse; 60 Vec *X, *Xloc, *mult, rhs, rhsloc; 61 PetscMemType memtype; 62 UserO *userO; 63 UserProlongRestr *userPR; 64 Ceed ceed; 65 CeedData *ceeddata; 66 CeedVector rhsceed, target; 67 CeedQFunction qferror, qfrestrict, qfprolong; 68 CeedOperator operror; 69 bpType bpchoice; 70 coarsenType coarsen; 71 72 ierr = PetscInitialize(&argc, &argv, NULL, help); 73 if (ierr) return ierr; 74 comm = PETSC_COMM_WORLD; 75 76 // Parse command line options 77 ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr); 78 bpchoice = CEED_BP3; 79 ierr = PetscOptionsEnum("-problem", 80 "CEED benchmark problem to solve", NULL, 81 bpTypes, (PetscEnum)bpchoice, (PetscEnum *)&bpchoice, 82 NULL); CHKERRQ(ierr); 83 ncompu = bpOptions[bpchoice].ncompu; 84 test_mode = PETSC_FALSE; 85 ierr = PetscOptionsBool("-test", 86 "Testing mode (do not print unless error is large)", 87 NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr); 88 benchmark_mode = PETSC_FALSE; 89 ierr = PetscOptionsBool("-benchmark", 90 "Benchmarking mode (prints benchmark statistics)", 91 NULL, benchmark_mode, &benchmark_mode, NULL); 92 CHKERRQ(ierr); 93 write_solution = PETSC_FALSE; 94 ierr = PetscOptionsBool("-write_solution", 95 "Write solution for visualization", 96 NULL, write_solution, &write_solution, NULL); 97 CHKERRQ(ierr); 98 degree = test_mode ? 3 : 2; 99 ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", 100 NULL, degree, °ree, NULL); CHKERRQ(ierr); 101 if (degree < 1) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, 102 "-degree %D must be at least 1", degree); 103 qextra = bpOptions[bpchoice].qextra; 104 ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", 105 NULL, qextra, &qextra, NULL); CHKERRQ(ierr); 106 ierr = PetscOptionsString("-ceed", "CEED resource specifier", 107 NULL, ceedresource, ceedresource, 108 sizeof(ceedresource), NULL); CHKERRQ(ierr); 109 coarsen = COARSEN_UNIFORM; 110 ierr = PetscOptionsEnum("-coarsen", 111 "Coarsening strategy to use", NULL, 112 coarsenTypes, (PetscEnum)coarsen, 113 (PetscEnum *)&coarsen, NULL); CHKERRQ(ierr); 114 read_mesh = PETSC_FALSE; 115 ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL, 116 filename, filename, sizeof(filename), &read_mesh); 117 CHKERRQ(ierr); 118 if (!read_mesh) { 119 PetscInt tmp = dim; 120 ierr = PetscOptionsIntArray("-cells","Number of cells per dimension", NULL, 121 melem, &tmp, NULL); CHKERRQ(ierr); 122 } 123 ierr = PetscOptionsEnd(); CHKERRQ(ierr); 124 125 // Set up libCEED 126 CeedInit(ceedresource, &ceed); 127 CeedMemType memtypebackend; 128 CeedGetPreferredMemType(ceed, &memtypebackend); 129 130 // Setup DM 131 if (read_mesh) { 132 ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dmorig); 133 CHKERRQ(ierr); 134 } else { 135 ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_FALSE, melem, NULL, 136 NULL, NULL, PETSC_TRUE,&dmorig); CHKERRQ(ierr); 137 } 138 139 { 140 DM dmDist = NULL; 141 PetscPartitioner part; 142 143 ierr = DMPlexGetPartitioner(dmorig, &part); CHKERRQ(ierr); 144 ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 145 ierr = DMPlexDistribute(dmorig, 0, NULL, &dmDist); CHKERRQ(ierr); 146 if (dmDist) { 147 ierr = DMDestroy(&dmorig); CHKERRQ(ierr); 148 dmorig = dmDist; 149 } 150 } 151 152 VecType vectype; 153 switch (memtypebackend) { 154 case CEED_MEM_HOST: vectype = VECSTANDARD; break; 155 case CEED_MEM_DEVICE: { 156 const char *resolved; 157 CeedGetResource(ceed, &resolved); 158 if (strstr(resolved, "/gpu/cuda")) vectype = VECCUDA; 159 else if (strstr(resolved, "/gpu/hip/occa")) 160 vectype = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678 161 else if (strstr(resolved, "/gpu/hip")) vectype = VECHIP; 162 else vectype = VECSTANDARD; 163 } 164 } 165 ierr = DMSetVecType(dmorig, vectype); CHKERRQ(ierr); 166 ierr = DMSetFromOptions(dmorig); CHKERRQ(ierr); 167 168 // Allocate arrays for PETSc objects for each level 169 switch (coarsen) { 170 case COARSEN_UNIFORM: 171 numlevels = degree; 172 break; 173 case COARSEN_LOGARITHMIC: 174 numlevels = ceil(log(degree)/log(2)) + 1; 175 break; 176 } 177 ierr = PetscMalloc1(numlevels, &leveldegrees); CHKERRQ(ierr); 178 fineLevel = numlevels - 1; 179 180 switch (coarsen) { 181 case COARSEN_UNIFORM: 182 for (int i=0; i<numlevels; i++) leveldegrees[i] = i + 1; 183 break; 184 case COARSEN_LOGARITHMIC: 185 for (int i=0; i<numlevels - 1; i++) leveldegrees[i] = pow(2,i); 186 leveldegrees[fineLevel] = degree; 187 break; 188 } 189 ierr = PetscMalloc1(numlevels, &dm); CHKERRQ(ierr); 190 ierr = PetscMalloc1(numlevels, &X); CHKERRQ(ierr); 191 ierr = PetscMalloc1(numlevels, &Xloc); CHKERRQ(ierr); 192 ierr = PetscMalloc1(numlevels, &mult); CHKERRQ(ierr); 193 ierr = PetscMalloc1(numlevels, &userO); CHKERRQ(ierr); 194 ierr = PetscMalloc1(numlevels, &userPR); CHKERRQ(ierr); 195 ierr = PetscMalloc1(numlevels, &matO); CHKERRQ(ierr); 196 ierr = PetscMalloc1(numlevels, &matPR); CHKERRQ(ierr); 197 ierr = PetscMalloc1(numlevels, &lsize); CHKERRQ(ierr); 198 ierr = PetscMalloc1(numlevels, &xlsize); CHKERRQ(ierr); 199 ierr = PetscMalloc1(numlevels, &gsize); CHKERRQ(ierr); 200 201 // Setup DM and Operator Mat Shells for each level 202 for (CeedInt i=0; i<numlevels; i++) { 203 // Create DM 204 ierr = DMClone(dmorig, &dm[i]); CHKERRQ(ierr); 205 ierr = DMGetVecType(dmorig, &vectype); CHKERRQ(ierr); 206 ierr = DMSetVecType(dm[i], vectype); CHKERRQ(ierr); 207 ierr = SetupDMByDegree(dm[i], leveldegrees[i], ncompu, bpchoice); 208 CHKERRQ(ierr); 209 210 // Create vectors 211 ierr = DMCreateGlobalVector(dm[i], &X[i]); CHKERRQ(ierr); 212 ierr = VecGetLocalSize(X[i], &lsize[i]); CHKERRQ(ierr); 213 ierr = VecGetSize(X[i], &gsize[i]); CHKERRQ(ierr); 214 ierr = DMCreateLocalVector(dm[i], &Xloc[i]); CHKERRQ(ierr); 215 ierr = VecGetSize(Xloc[i], &xlsize[i]); CHKERRQ(ierr); 216 217 // Operator 218 ierr = PetscMalloc1(1, &userO[i]); CHKERRQ(ierr); 219 ierr = MatCreateShell(comm, lsize[i], lsize[i], gsize[i], gsize[i], 220 userO[i], &matO[i]); CHKERRQ(ierr); 221 ierr = MatShellSetOperation(matO[i], MATOP_MULT, 222 (void(*)(void))MatMult_Ceed); CHKERRQ(ierr); 223 ierr = MatShellSetOperation(matO[i], MATOP_GET_DIAGONAL, 224 (void(*)(void))MatGetDiag); CHKERRQ(ierr); 225 ierr = MatShellSetVecType(matO[i], vectype); CHKERRQ(ierr); 226 227 // Level transfers 228 if (i > 0) { 229 // Interp 230 ierr = PetscMalloc1(1, &userPR[i]); CHKERRQ(ierr); 231 ierr = MatCreateShell(comm, lsize[i], lsize[i-1], gsize[i], gsize[i-1], 232 userPR[i], &matPR[i]); CHKERRQ(ierr); 233 ierr = MatShellSetOperation(matPR[i], MATOP_MULT, 234 (void(*)(void))MatMult_Prolong); 235 CHKERRQ(ierr); 236 ierr = MatShellSetOperation(matPR[i], MATOP_MULT_TRANSPOSE, 237 (void(*)(void))MatMult_Restrict); 238 CHKERRQ(ierr); 239 ierr = MatShellSetVecType(matPR[i], vectype); CHKERRQ(ierr); 240 } 241 } 242 ierr = VecDuplicate(X[fineLevel], &rhs); CHKERRQ(ierr); 243 244 // Print global grid information 245 if (!test_mode) { 246 PetscInt P = degree + 1, Q = P + qextra; 247 248 const char *usedresource; 249 CeedGetResource(ceed, &usedresource); 250 251 ierr = VecGetType(X[0], &vectype); CHKERRQ(ierr); 252 253 ierr = PetscPrintf(comm, 254 "\n-- CEED Benchmark Problem %d -- libCEED + PETSc + PCMG --\n" 255 " PETSc:\n" 256 " PETSc Vec Type : %s\n" 257 " libCEED:\n" 258 " libCEED Backend : %s\n" 259 " libCEED Backend MemType : %s\n" 260 " Mesh:\n" 261 " Number of 1D Basis Nodes (p) : %d\n" 262 " Number of 1D Quadrature Points (q) : %d\n" 263 " Global Nodes : %D\n" 264 " Owned Nodes : %D\n" 265 " DoF per node : %D\n" 266 " Multigrid:\n" 267 " Number of Levels : %d\n", 268 bpchoice+1, vectype, usedresource, 269 CeedMemTypes[memtypebackend], 270 P, Q, gsize[fineLevel]/ncompu, lsize[fineLevel]/ncompu, 271 ncompu, numlevels); CHKERRQ(ierr); 272 } 273 274 // Create RHS vector 275 ierr = VecDuplicate(Xloc[fineLevel], &rhsloc); CHKERRQ(ierr); 276 ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr); 277 ierr = VecGetArrayAndMemType(rhsloc, &r, &memtype); CHKERRQ(ierr); 278 CeedVectorCreate(ceed, xlsize[fineLevel], &rhsceed); 279 CeedVectorSetArray(rhsceed, MemTypeP2C(memtype), CEED_USE_POINTER, r); 280 281 // Set up libCEED operators on each level 282 ierr = PetscMalloc1(numlevels, &ceeddata); CHKERRQ(ierr); 283 for (int i=0; i<numlevels; i++) { 284 // Print level information 285 if (!test_mode && (i == 0 || i == fineLevel)) { 286 ierr = PetscPrintf(comm," Level %D (%s):\n" 287 " Number of 1D Basis Nodes (p) : %d\n" 288 " Global Nodes : %D\n" 289 " Owned Nodes : %D\n", 290 i, (i? "fine" : "coarse"), leveldegrees[i] + 1, 291 gsize[i]/ncompu, lsize[i]/ncompu); CHKERRQ(ierr); 292 } 293 ierr = PetscMalloc1(1, &ceeddata[i]); CHKERRQ(ierr); 294 ierr = SetupLibceedByDegree(dm[i], ceed, leveldegrees[i], dim, qextra, 295 ncompu, gsize[i], xlsize[i], bpchoice, 296 ceeddata[i], i==(fineLevel), rhsceed, &target); 297 CHKERRQ(ierr); 298 } 299 300 // Gather RHS 301 CeedVectorTakeArray(rhsceed, MemTypeP2C(memtype), NULL); 302 ierr = VecRestoreArrayAndMemType(rhsloc, &r); CHKERRQ(ierr); 303 ierr = VecZeroEntries(rhs); CHKERRQ(ierr); 304 ierr = DMLocalToGlobal(dm[fineLevel], rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr); 305 CeedVectorDestroy(&rhsceed); 306 307 // Create the restriction/interpolation QFunction 308 CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_NONE, CEED_EVAL_INTERP, 309 &qfrestrict); 310 CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_INTERP, CEED_EVAL_NONE, 311 &qfprolong); 312 313 // Set up libCEED level transfer operators 314 ierr = CeedLevelTransferSetup(ceed, numlevels, ncompu, bpchoice, ceeddata, 315 leveldegrees, qfrestrict, qfprolong); 316 CHKERRQ(ierr); 317 318 // Create the error QFunction 319 CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].error, 320 bpOptions[bpchoice].errorfname, &qferror); 321 CeedQFunctionAddInput(qferror, "u", ncompu, CEED_EVAL_INTERP); 322 CeedQFunctionAddInput(qferror, "true_soln", ncompu, CEED_EVAL_NONE); 323 CeedQFunctionAddOutput(qferror, "error", ncompu, CEED_EVAL_NONE); 324 325 // Create the error operator 326 CeedOperatorCreate(ceed, qferror, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, 327 &operror); 328 CeedOperatorSetField(operror, "u", ceeddata[fineLevel]->Erestrictu, 329 ceeddata[fineLevel]->basisu, CEED_VECTOR_ACTIVE); 330 CeedOperatorSetField(operror, "true_soln", ceeddata[fineLevel]->Erestrictui, 331 CEED_BASIS_COLLOCATED, target); 332 CeedOperatorSetField(operror, "error", ceeddata[fineLevel]->Erestrictui, 333 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 334 335 // Calculate multiplicity 336 for (int i=0; i<numlevels; i++) { 337 PetscScalar *x; 338 339 // CEED vector 340 ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr); 341 ierr = VecGetArray(Xloc[i], &x); CHKERRQ(ierr); 342 CeedVectorSetArray(ceeddata[i]->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x); 343 344 // Multiplicity 345 CeedElemRestrictionGetMultiplicity(ceeddata[i]->Erestrictu, 346 ceeddata[i]->xceed); 347 CeedVectorSyncArray(ceeddata[i]->xceed, CEED_MEM_HOST); 348 349 // Restore vector 350 ierr = VecRestoreArray(Xloc[i], &x); CHKERRQ(ierr); 351 352 // Creat mult vector 353 ierr = VecDuplicate(Xloc[i], &mult[i]); CHKERRQ(ierr); 354 355 // Local-to-global 356 ierr = VecZeroEntries(X[i]); CHKERRQ(ierr); 357 ierr = DMLocalToGlobal(dm[i], Xloc[i], ADD_VALUES, X[i]); 358 CHKERRQ(ierr); 359 ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr); 360 361 // Global-to-local 362 ierr = DMGlobalToLocal(dm[i], X[i], INSERT_VALUES, mult[i]); 363 CHKERRQ(ierr); 364 ierr = VecZeroEntries(X[i]); CHKERRQ(ierr); 365 366 // Multiplicity scaling 367 ierr = VecReciprocal(mult[i]); 368 } 369 370 // Set up Mat 371 for (int i=0; i<numlevels; i++) { 372 // User Operator 373 userO[i]->comm = comm; 374 userO[i]->dm = dm[i]; 375 userO[i]->Xloc = Xloc[i]; 376 ierr = VecDuplicate(Xloc[i], &userO[i]->Yloc); CHKERRQ(ierr); 377 userO[i]->xceed = ceeddata[i]->xceed; 378 userO[i]->yceed = ceeddata[i]->yceed; 379 userO[i]->op = ceeddata[i]->opapply; 380 userO[i]->ceed = ceed; 381 382 if (i > 0) { 383 // Prolongation/Restriction Operator 384 userPR[i]->comm = comm; 385 userPR[i]->dmf = dm[i]; 386 userPR[i]->dmc = dm[i-1]; 387 userPR[i]->locvecc = Xloc[i-1]; 388 userPR[i]->locvecf = userO[i]->Yloc; 389 userPR[i]->multvec = mult[i]; 390 userPR[i]->ceedvecc = userO[i-1]->xceed; 391 userPR[i]->ceedvecf = userO[i]->yceed; 392 userPR[i]->opprolong = ceeddata[i]->opprolong; 393 userPR[i]->oprestrict = ceeddata[i]->oprestrict; 394 userPR[i]->ceed = ceed; 395 } 396 } 397 398 // Setup dummy SNES for AMG coarse solve 399 ierr = SNESCreate(comm, &snesdummy); CHKERRQ(ierr); 400 ierr = SNESSetDM(snesdummy, dm[0]); CHKERRQ(ierr); 401 ierr = SNESSetSolution(snesdummy, X[0]); CHKERRQ(ierr); 402 403 // -- Jacobian matrix 404 ierr = DMSetMatType(dm[0], MATAIJ); CHKERRQ(ierr); 405 ierr = DMCreateMatrix(dm[0], &matcoarse); CHKERRQ(ierr); 406 ierr = SNESSetJacobian(snesdummy, matcoarse, matcoarse, NULL, 407 NULL); CHKERRQ(ierr); 408 409 // -- Residual evaluation function 410 ierr = SNESSetFunction(snesdummy, X[0], FormResidual_Ceed, 411 userO[0]); CHKERRQ(ierr); 412 413 // -- Form Jacobian 414 ierr = SNESComputeJacobianDefaultColor(snesdummy, X[0], matO[0], 415 matcoarse, NULL); CHKERRQ(ierr); 416 417 // Set up KSP 418 ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr); 419 { 420 ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 421 ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 422 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 423 PETSC_DEFAULT); CHKERRQ(ierr); 424 } 425 ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr); 426 ierr = KSPSetOperators(ksp, matO[fineLevel], matO[fineLevel]); 427 CHKERRQ(ierr); 428 429 // Set up PCMG 430 ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 431 PCMGCycleType pcgmcycletype = PC_MG_CYCLE_V; 432 { 433 ierr = PCSetType(pc, PCMG); CHKERRQ(ierr); 434 435 // PCMG levels 436 ierr = PCMGSetLevels(pc, numlevels, NULL); CHKERRQ(ierr); 437 for (int i=0; i<numlevels; i++) { 438 // Smoother 439 KSP smoother; 440 PC smoother_pc; 441 ierr = PCMGGetSmoother(pc, i, &smoother); CHKERRQ(ierr); 442 ierr = KSPSetType(smoother, KSPCHEBYSHEV); CHKERRQ(ierr); 443 ierr = KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1); CHKERRQ(ierr); 444 ierr = KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE); CHKERRQ(ierr); 445 ierr = KSPSetOperators(smoother, matO[i], matO[i]); CHKERRQ(ierr); 446 ierr = KSPGetPC(smoother, &smoother_pc); CHKERRQ(ierr); 447 ierr = PCSetType(smoother_pc, PCJACOBI); CHKERRQ(ierr); 448 ierr = PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 449 450 // Work vector 451 if (i < numlevels - 1) { 452 ierr = PCMGSetX(pc, i, X[i]); CHKERRQ(ierr); 453 } 454 455 // Level transfers 456 if (i > 0) { 457 // Interpolation 458 ierr = PCMGSetInterpolation(pc, i, matPR[i]); CHKERRQ(ierr); 459 } 460 461 // Coarse solve 462 KSP coarse; 463 PC coarse_pc; 464 ierr = PCMGGetCoarseSolve(pc, &coarse); CHKERRQ(ierr); 465 ierr = KSPSetType(coarse, KSPPREONLY); CHKERRQ(ierr); 466 ierr = KSPSetOperators(coarse, matcoarse, matcoarse); CHKERRQ(ierr); 467 468 ierr = KSPGetPC(coarse, &coarse_pc); CHKERRQ(ierr); 469 ierr = PCSetType(coarse_pc, PCGAMG); CHKERRQ(ierr); 470 471 ierr = KSPSetOptionsPrefix(coarse, "coarse_"); CHKERRQ(ierr); 472 ierr = PCSetOptionsPrefix(coarse_pc, "coarse_"); CHKERRQ(ierr); 473 ierr = KSPSetFromOptions(coarse); CHKERRQ(ierr); 474 ierr = PCSetFromOptions(coarse_pc); CHKERRQ(ierr); 475 } 476 477 // PCMG options 478 ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr); 479 ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr); 480 ierr = PCMGSetCycleType(pc, pcgmcycletype); CHKERRQ(ierr); 481 } 482 483 // First run, if benchmarking 484 if (benchmark_mode) { 485 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1); 486 CHKERRQ(ierr); 487 ierr = VecZeroEntries(X[fineLevel]); CHKERRQ(ierr); 488 my_rt_start = MPI_Wtime(); 489 ierr = KSPSolve(ksp, rhs, X[fineLevel]); CHKERRQ(ierr); 490 my_rt = MPI_Wtime() - my_rt_start; 491 ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm); 492 CHKERRQ(ierr); 493 // Set maxits based on first iteration timing 494 if (my_rt > 0.02) { 495 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5); 496 CHKERRQ(ierr); 497 } else { 498 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20); 499 CHKERRQ(ierr); 500 } 501 } 502 503 // Timed solve 504 ierr = VecZeroEntries(X[fineLevel]); CHKERRQ(ierr); 505 ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr); 506 507 // -- Performance logging 508 ierr = PetscLogStageRegister("Solve Stage", &solvestage); CHKERRQ(ierr); 509 ierr = PetscLogStagePush(solvestage); CHKERRQ(ierr); 510 511 // -- Solve 512 my_rt_start = MPI_Wtime(); 513 ierr = KSPSolve(ksp, rhs, X[fineLevel]); CHKERRQ(ierr); 514 my_rt = MPI_Wtime() - my_rt_start; 515 516 517 // -- Performance logging 518 ierr = PetscLogStagePop(); 519 520 // Output results 521 { 522 KSPType ksptype; 523 PCMGType pcmgtype; 524 KSPConvergedReason reason; 525 PetscReal rnorm; 526 PetscInt its; 527 ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr); 528 ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr); 529 ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr); 530 ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr); 531 ierr = PCMGGetType(pc, &pcmgtype); CHKERRQ(ierr); 532 if (!test_mode || reason < 0 || rnorm > 1e-8) { 533 ierr = PetscPrintf(comm, 534 " KSP:\n" 535 " KSP Type : %s\n" 536 " KSP Convergence : %s\n" 537 " Total KSP Iterations : %D\n" 538 " Final rnorm : %e\n", 539 ksptype, KSPConvergedReasons[reason], its, 540 (double)rnorm); CHKERRQ(ierr); 541 ierr = PetscPrintf(comm, 542 " PCMG:\n" 543 " PCMG Type : %s\n" 544 " PCMG Cycle Type : %s\n", 545 PCMGTypes[pcmgtype], 546 PCMGCycleTypes[pcgmcycletype]); CHKERRQ(ierr); 547 } 548 if (!test_mode) { 549 ierr = PetscPrintf(comm," Performance:\n"); CHKERRQ(ierr); 550 } 551 { 552 PetscReal maxerror; 553 ierr = ComputeErrorMax(userO[fineLevel], operror, X[fineLevel], target, 554 &maxerror); CHKERRQ(ierr); 555 PetscReal tol = 5e-2; 556 if (!test_mode || maxerror > tol) { 557 ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm); 558 CHKERRQ(ierr); 559 ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm); 560 CHKERRQ(ierr); 561 ierr = PetscPrintf(comm, 562 " Pointwise Error (max) : %e\n" 563 " CG Solve Time : %g (%g) sec\n", 564 (double)maxerror, rt_max, rt_min); CHKERRQ(ierr); 565 } 566 } 567 if (benchmark_mode && (!test_mode)) { 568 ierr = PetscPrintf(comm, 569 " DoFs/Sec in CG : %g (%g) million\n", 570 1e-6*gsize[fineLevel]*its/rt_max, 571 1e-6*gsize[fineLevel]*its/rt_min); 572 CHKERRQ(ierr); 573 } 574 } 575 576 if (write_solution) { 577 PetscViewer vtkviewersoln; 578 579 ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr); 580 ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr); 581 ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr); 582 ierr = VecView(X[fineLevel], vtkviewersoln); CHKERRQ(ierr); 583 ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr); 584 } 585 586 // Cleanup 587 for (int i=0; i<numlevels; i++) { 588 ierr = VecDestroy(&X[i]); CHKERRQ(ierr); 589 ierr = VecDestroy(&Xloc[i]); CHKERRQ(ierr); 590 ierr = VecDestroy(&mult[i]); CHKERRQ(ierr); 591 ierr = VecDestroy(&userO[i]->Yloc); CHKERRQ(ierr); 592 ierr = MatDestroy(&matO[i]); CHKERRQ(ierr); 593 ierr = PetscFree(userO[i]); CHKERRQ(ierr); 594 if (i > 0) { 595 ierr = MatDestroy(&matPR[i]); CHKERRQ(ierr); 596 ierr = PetscFree(userPR[i]); CHKERRQ(ierr); 597 } 598 ierr = CeedDataDestroy(i, ceeddata[i]); CHKERRQ(ierr); 599 ierr = DMDestroy(&dm[i]); CHKERRQ(ierr); 600 } 601 ierr = PetscFree(leveldegrees); CHKERRQ(ierr); 602 ierr = PetscFree(dm); CHKERRQ(ierr); 603 ierr = PetscFree(X); CHKERRQ(ierr); 604 ierr = PetscFree(Xloc); CHKERRQ(ierr); 605 ierr = PetscFree(mult); CHKERRQ(ierr); 606 ierr = PetscFree(matO); CHKERRQ(ierr); 607 ierr = PetscFree(matPR); CHKERRQ(ierr); 608 ierr = PetscFree(ceeddata); CHKERRQ(ierr); 609 ierr = PetscFree(userO); CHKERRQ(ierr); 610 ierr = PetscFree(userPR); CHKERRQ(ierr); 611 ierr = PetscFree(lsize); CHKERRQ(ierr); 612 ierr = PetscFree(xlsize); CHKERRQ(ierr); 613 ierr = PetscFree(gsize); CHKERRQ(ierr); 614 ierr = VecDestroy(&rhs); CHKERRQ(ierr); 615 ierr = VecDestroy(&rhsloc); CHKERRQ(ierr); 616 ierr = MatDestroy(&matcoarse); CHKERRQ(ierr); 617 ierr = KSPDestroy(&ksp); CHKERRQ(ierr); 618 ierr = SNESDestroy(&snesdummy); CHKERRQ(ierr); 619 ierr = DMDestroy(&dmorig); CHKERRQ(ierr); 620 CeedVectorDestroy(&target); 621 CeedQFunctionDestroy(&qferror); 622 CeedQFunctionDestroy(&qfrestrict); 623 CeedQFunctionDestroy(&qfprolong); 624 CeedOperatorDestroy(&operror); 625 CeedDestroy(&ceed); 626 return PetscFinalize(); 627 } 628