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