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 #include "bps.h" 42 43 int main(int argc, char **argv) { 44 PetscInt ierr; 45 MPI_Comm comm; 46 char filename[PETSC_MAX_PATH_LEN], 47 ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self"; 48 double my_rt_start, my_rt, rt_min, rt_max; 49 PetscInt degree = 3, qextra, *lsize, *xlsize, *gsize, dim = 3, fineLevel, 50 melem[3] = {3, 3, 3}, ncompu = 1, numlevels = degree, *leveldegrees; 51 PetscScalar *r; 52 PetscScalar eps = 1.0; 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 ierr = PetscOptionsScalar("-eps", 99 "Epsilon parameter for Kershaw mesh transformation", 100 NULL, eps, &eps, NULL); 101 if (eps > 1 || eps <= 0) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, 102 "-eps %D must be (0,1]", eps); 103 degree = test_mode ? 3 : 2; 104 ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", 105 NULL, degree, °ree, NULL); CHKERRQ(ierr); 106 if (degree < 1) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, 107 "-degree %D must be at least 1", degree); 108 qextra = bpOptions[bpchoice].qextra; 109 ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", 110 NULL, qextra, &qextra, NULL); CHKERRQ(ierr); 111 ierr = PetscOptionsString("-ceed", "CEED resource specifier", 112 NULL, ceedresource, ceedresource, 113 sizeof(ceedresource), NULL); CHKERRQ(ierr); 114 coarsen = COARSEN_UNIFORM; 115 ierr = PetscOptionsEnum("-coarsen", 116 "Coarsening strategy to use", NULL, 117 coarsenTypes, (PetscEnum)coarsen, 118 (PetscEnum *)&coarsen, 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 ierr = PetscOptionsEnd(); CHKERRQ(ierr); 129 130 // Set up libCEED 131 CeedInit(ceedresource, &ceed); 132 CeedMemType memtypebackend; 133 CeedGetPreferredMemType(ceed, &memtypebackend); 134 135 // Setup DM 136 if (read_mesh) { 137 ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dmorig); 138 CHKERRQ(ierr); 139 } else { 140 ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_FALSE, melem, NULL, 141 NULL, NULL, PETSC_TRUE,&dmorig); CHKERRQ(ierr); 142 } 143 144 { 145 DM dmDist = NULL; 146 PetscPartitioner part; 147 148 ierr = DMPlexGetPartitioner(dmorig, &part); CHKERRQ(ierr); 149 ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 150 ierr = DMPlexDistribute(dmorig, 0, NULL, &dmDist); CHKERRQ(ierr); 151 if (dmDist) { 152 ierr = DMDestroy(&dmorig); CHKERRQ(ierr); 153 dmorig = dmDist; 154 } 155 } 156 157 // apply Kershaw mesh transformation 158 ierr = kershaw(dmorig, eps); CHKERRQ(ierr); 159 160 VecType vectype; 161 switch (memtypebackend) { 162 case CEED_MEM_HOST: vectype = VECSTANDARD; break; 163 case CEED_MEM_DEVICE: { 164 const char *resolved; 165 CeedGetResource(ceed, &resolved); 166 if (strstr(resolved, "/gpu/cuda")) vectype = VECCUDA; 167 else if (strstr(resolved, "/gpu/hip/occa")) 168 vectype = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678 169 else if (strstr(resolved, "/gpu/hip")) vectype = VECHIP; 170 else vectype = VECSTANDARD; 171 } 172 } 173 ierr = DMSetVecType(dmorig, vectype); CHKERRQ(ierr); 174 ierr = DMSetFromOptions(dmorig); CHKERRQ(ierr); 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 = DMGetVecType(dmorig, &vectype); CHKERRQ(ierr); 214 ierr = DMSetVecType(dm[i], vectype); CHKERRQ(ierr); 215 PetscInt dim; 216 ierr = DMGetDimension(dm[i], &dim); CHKERRQ(ierr); 217 ierr = SetupDMByDegree(dm[i], leveldegrees[i], ncompu, dim, 218 bpOptions[bpchoice].enforcebc, bpOptions[bpchoice].bcsfunc); 219 CHKERRQ(ierr); 220 221 // Create vectors 222 ierr = DMCreateGlobalVector(dm[i], &X[i]); CHKERRQ(ierr); 223 ierr = VecGetLocalSize(X[i], &lsize[i]); CHKERRQ(ierr); 224 ierr = VecGetSize(X[i], &gsize[i]); CHKERRQ(ierr); 225 ierr = DMCreateLocalVector(dm[i], &Xloc[i]); CHKERRQ(ierr); 226 ierr = VecGetSize(Xloc[i], &xlsize[i]); CHKERRQ(ierr); 227 228 // Operator 229 ierr = PetscMalloc1(1, &userO[i]); CHKERRQ(ierr); 230 ierr = MatCreateShell(comm, lsize[i], lsize[i], gsize[i], gsize[i], 231 userO[i], &matO[i]); CHKERRQ(ierr); 232 ierr = MatShellSetOperation(matO[i], MATOP_MULT, 233 (void(*)(void))MatMult_Ceed); CHKERRQ(ierr); 234 ierr = MatShellSetOperation(matO[i], MATOP_GET_DIAGONAL, 235 (void(*)(void))MatGetDiag); CHKERRQ(ierr); 236 ierr = MatShellSetVecType(matO[i], vectype); CHKERRQ(ierr); 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 ierr = MatShellSetVecType(matPR[i], vectype); CHKERRQ(ierr); 251 } 252 } 253 ierr = VecDuplicate(X[fineLevel], &rhs); CHKERRQ(ierr); 254 255 // Print global grid information 256 if (!test_mode) { 257 PetscInt P = degree + 1, Q = P + qextra; 258 259 const char *usedresource; 260 CeedGetResource(ceed, &usedresource); 261 262 ierr = VecGetType(X[0], &vectype); CHKERRQ(ierr); 263 264 ierr = PetscPrintf(comm, 265 "\n-- CEED Benchmark Problem %d -- libCEED + PETSc + PCMG --\n" 266 " PETSc:\n" 267 " PETSc Vec Type : %s\n" 268 " libCEED:\n" 269 " libCEED Backend : %s\n" 270 " libCEED Backend MemType : %s\n" 271 " Mesh:\n" 272 " Number of 1D Basis Nodes (p) : %d\n" 273 " Number of 1D Quadrature Points (q) : %d\n" 274 " Global Nodes : %D\n" 275 " Owned Nodes : %D\n" 276 " DoF per node : %D\n" 277 " Multigrid:\n" 278 " Number of Levels : %d\n", 279 bpchoice+1, vectype, usedresource, 280 CeedMemTypes[memtypebackend], 281 P, Q, gsize[fineLevel]/ncompu, lsize[fineLevel]/ncompu, 282 ncompu, numlevels); CHKERRQ(ierr); 283 } 284 285 // Create RHS vector 286 ierr = VecDuplicate(Xloc[fineLevel], &rhsloc); CHKERRQ(ierr); 287 ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr); 288 ierr = VecGetArrayAndMemType(rhsloc, &r, &memtype); CHKERRQ(ierr); 289 CeedVectorCreate(ceed, xlsize[fineLevel], &rhsceed); 290 CeedVectorSetArray(rhsceed, MemTypeP2C(memtype), CEED_USE_POINTER, r); 291 292 // Set up libCEED operators on each level 293 ierr = PetscMalloc1(numlevels, &ceeddata); CHKERRQ(ierr); 294 for (int i=0; i<numlevels; i++) { 295 // Print level information 296 if (!test_mode && (i == 0 || i == fineLevel)) { 297 ierr = PetscPrintf(comm," Level %D (%s):\n" 298 " Number of 1D Basis Nodes (p) : %d\n" 299 " Global Nodes : %D\n" 300 " Owned Nodes : %D\n", 301 i, (i? "fine" : "coarse"), leveldegrees[i] + 1, 302 gsize[i]/ncompu, lsize[i]/ncompu); CHKERRQ(ierr); 303 } 304 ierr = PetscMalloc1(1, &ceeddata[i]); CHKERRQ(ierr); 305 ierr = SetupLibceedByDegree(dm[i], ceed, leveldegrees[i], dim, qextra, 306 dim, ncompu, gsize[i], xlsize[i], bpOptions[bpchoice], 307 ceeddata[i], i==(fineLevel), rhsceed, &target); 308 CHKERRQ(ierr); 309 } 310 311 // Gather RHS 312 CeedVectorTakeArray(rhsceed, MemTypeP2C(memtype), NULL); 313 ierr = VecRestoreArrayAndMemType(rhsloc, &r); CHKERRQ(ierr); 314 ierr = VecZeroEntries(rhs); CHKERRQ(ierr); 315 ierr = DMLocalToGlobal(dm[fineLevel], rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr); 316 CeedVectorDestroy(&rhsceed); 317 318 // Create the restriction/interpolation QFunction 319 CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_NONE, CEED_EVAL_INTERP, 320 &qfrestrict); 321 CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_INTERP, CEED_EVAL_NONE, 322 &qfprolong); 323 324 // Set up libCEED level transfer operators 325 ierr = CeedLevelTransferSetup(ceed, numlevels, ncompu, ceeddata, leveldegrees, 326 qfrestrict, qfprolong); CHKERRQ(ierr); 327 328 // Create the error QFunction 329 CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].error, 330 bpOptions[bpchoice].errorfname, &qferror); 331 CeedQFunctionAddInput(qferror, "u", ncompu, CEED_EVAL_INTERP); 332 CeedQFunctionAddInput(qferror, "true_soln", ncompu, CEED_EVAL_NONE); 333 CeedQFunctionAddOutput(qferror, "error", ncompu, CEED_EVAL_NONE); 334 335 // Create the error operator 336 CeedOperatorCreate(ceed, qferror, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, 337 &operror); 338 CeedOperatorSetField(operror, "u", ceeddata[fineLevel]->Erestrictu, 339 ceeddata[fineLevel]->basisu, CEED_VECTOR_ACTIVE); 340 CeedOperatorSetField(operror, "true_soln", ceeddata[fineLevel]->Erestrictui, 341 CEED_BASIS_COLLOCATED, target); 342 CeedOperatorSetField(operror, "error", ceeddata[fineLevel]->Erestrictui, 343 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 344 345 // Calculate multiplicity 346 for (int i=0; i<numlevels; i++) { 347 PetscScalar *x; 348 349 // CEED vector 350 ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr); 351 ierr = VecGetArray(Xloc[i], &x); CHKERRQ(ierr); 352 CeedVectorSetArray(ceeddata[i]->Xceed, CEED_MEM_HOST, CEED_USE_POINTER, x); 353 354 // Multiplicity 355 CeedElemRestrictionGetMultiplicity(ceeddata[i]->Erestrictu, 356 ceeddata[i]->Xceed); 357 CeedVectorSyncArray(ceeddata[i]->Xceed, CEED_MEM_HOST); 358 359 // Restore vector 360 ierr = VecRestoreArray(Xloc[i], &x); CHKERRQ(ierr); 361 362 // Creat mult vector 363 ierr = VecDuplicate(Xloc[i], &mult[i]); CHKERRQ(ierr); 364 365 // Local-to-global 366 ierr = VecZeroEntries(X[i]); CHKERRQ(ierr); 367 ierr = DMLocalToGlobal(dm[i], Xloc[i], ADD_VALUES, X[i]); 368 CHKERRQ(ierr); 369 ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr); 370 371 // Global-to-local 372 ierr = DMGlobalToLocal(dm[i], X[i], INSERT_VALUES, mult[i]); 373 CHKERRQ(ierr); 374 ierr = VecZeroEntries(X[i]); CHKERRQ(ierr); 375 376 // Multiplicity scaling 377 ierr = VecReciprocal(mult[i]); 378 } 379 380 // Set up Mat 381 for (int i=0; i<numlevels; i++) { 382 // User Operator 383 userO[i]->comm = comm; 384 userO[i]->dm = dm[i]; 385 userO[i]->Xloc = Xloc[i]; 386 ierr = VecDuplicate(Xloc[i], &userO[i]->Yloc); CHKERRQ(ierr); 387 userO[i]->Xceed = ceeddata[i]->Xceed; 388 userO[i]->Yceed = ceeddata[i]->Yceed; 389 userO[i]->op = ceeddata[i]->opApply; 390 userO[i]->ceed = ceed; 391 392 if (i > 0) { 393 // Prolongation/Restriction Operator 394 userPR[i]->comm = comm; 395 userPR[i]->dmf = dm[i]; 396 userPR[i]->dmc = dm[i-1]; 397 userPR[i]->locvecc = Xloc[i-1]; 398 userPR[i]->locvecf = userO[i]->Yloc; 399 userPR[i]->multvec = mult[i]; 400 userPR[i]->ceedvecc = userO[i-1]->Xceed; 401 userPR[i]->ceedvecf = userO[i]->Yceed; 402 userPR[i]->opProlong = ceeddata[i]->opProlong; 403 userPR[i]->opRestrict = ceeddata[i]->opRestrict; 404 userPR[i]->ceed = ceed; 405 } 406 } 407 408 // Setup dummy SNES for AMG coarse solve 409 ierr = SNESCreate(comm, &snesdummy); CHKERRQ(ierr); 410 ierr = SNESSetDM(snesdummy, dm[0]); CHKERRQ(ierr); 411 ierr = SNESSetSolution(snesdummy, X[0]); CHKERRQ(ierr); 412 413 // -- Jacobian matrix 414 ierr = DMSetMatType(dm[0], MATAIJ); CHKERRQ(ierr); 415 ierr = DMCreateMatrix(dm[0], &matcoarse); CHKERRQ(ierr); 416 ierr = SNESSetJacobian(snesdummy, matcoarse, matcoarse, NULL, 417 NULL); CHKERRQ(ierr); 418 419 // -- Residual evaluation function 420 ierr = SNESSetFunction(snesdummy, X[0], FormResidual_Ceed, 421 userO[0]); CHKERRQ(ierr); 422 423 // -- Form Jacobian 424 ierr = SNESComputeJacobianDefaultColor(snesdummy, X[0], matO[0], 425 matcoarse, NULL); CHKERRQ(ierr); 426 427 // Set up KSP 428 ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr); 429 { 430 ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 431 ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 432 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 433 PETSC_DEFAULT); CHKERRQ(ierr); 434 } 435 ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr); 436 ierr = KSPSetOperators(ksp, matO[fineLevel], matO[fineLevel]); 437 CHKERRQ(ierr); 438 439 // Set up PCMG 440 ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 441 PCMGCycleType pcgmcycletype = PC_MG_CYCLE_V; 442 { 443 ierr = PCSetType(pc, PCMG); CHKERRQ(ierr); 444 445 // PCMG levels 446 ierr = PCMGSetLevels(pc, numlevels, NULL); CHKERRQ(ierr); 447 for (int i=0; i<numlevels; i++) { 448 // Smoother 449 KSP smoother; 450 PC smoother_pc; 451 ierr = PCMGGetSmoother(pc, i, &smoother); CHKERRQ(ierr); 452 ierr = KSPSetType(smoother, KSPCHEBYSHEV); CHKERRQ(ierr); 453 ierr = KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1); CHKERRQ(ierr); 454 ierr = KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE); CHKERRQ(ierr); 455 ierr = KSPSetOperators(smoother, matO[i], matO[i]); CHKERRQ(ierr); 456 ierr = KSPGetPC(smoother, &smoother_pc); CHKERRQ(ierr); 457 ierr = PCSetType(smoother_pc, PCJACOBI); CHKERRQ(ierr); 458 ierr = PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 459 460 // Work vector 461 if (i < numlevels - 1) { 462 ierr = PCMGSetX(pc, i, X[i]); CHKERRQ(ierr); 463 } 464 465 // Level transfers 466 if (i > 0) { 467 // Interpolation 468 ierr = PCMGSetInterpolation(pc, i, matPR[i]); CHKERRQ(ierr); 469 } 470 471 // Coarse solve 472 KSP coarse; 473 PC coarse_pc; 474 ierr = PCMGGetCoarseSolve(pc, &coarse); CHKERRQ(ierr); 475 ierr = KSPSetType(coarse, KSPPREONLY); CHKERRQ(ierr); 476 ierr = KSPSetOperators(coarse, matcoarse, matcoarse); CHKERRQ(ierr); 477 478 ierr = KSPGetPC(coarse, &coarse_pc); CHKERRQ(ierr); 479 ierr = PCSetType(coarse_pc, PCGAMG); CHKERRQ(ierr); 480 481 ierr = KSPSetOptionsPrefix(coarse, "coarse_"); CHKERRQ(ierr); 482 ierr = PCSetOptionsPrefix(coarse_pc, "coarse_"); CHKERRQ(ierr); 483 ierr = KSPSetFromOptions(coarse); CHKERRQ(ierr); 484 ierr = PCSetFromOptions(coarse_pc); CHKERRQ(ierr); 485 } 486 487 // PCMG options 488 ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr); 489 ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr); 490 ierr = PCMGSetCycleType(pc, pcgmcycletype); CHKERRQ(ierr); 491 } 492 493 // First run, if benchmarking 494 if (benchmark_mode) { 495 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1); 496 CHKERRQ(ierr); 497 ierr = VecZeroEntries(X[fineLevel]); CHKERRQ(ierr); 498 my_rt_start = MPI_Wtime(); 499 ierr = KSPSolve(ksp, rhs, X[fineLevel]); CHKERRQ(ierr); 500 my_rt = MPI_Wtime() - my_rt_start; 501 ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm); 502 CHKERRQ(ierr); 503 // Set maxits based on first iteration timing 504 if (my_rt > 0.02) { 505 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5); 506 CHKERRQ(ierr); 507 } else { 508 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20); 509 CHKERRQ(ierr); 510 } 511 } 512 513 // Timed solve 514 ierr = VecZeroEntries(X[fineLevel]); CHKERRQ(ierr); 515 ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr); 516 517 // -- Performance logging 518 ierr = PetscLogStageRegister("Solve Stage", &solvestage); CHKERRQ(ierr); 519 ierr = PetscLogStagePush(solvestage); CHKERRQ(ierr); 520 521 // -- Solve 522 my_rt_start = MPI_Wtime(); 523 ierr = KSPSolve(ksp, rhs, X[fineLevel]); CHKERRQ(ierr); 524 my_rt = MPI_Wtime() - my_rt_start; 525 526 527 // -- Performance logging 528 ierr = PetscLogStagePop(); 529 530 // Output results 531 { 532 KSPType ksptype; 533 PCMGType pcmgtype; 534 KSPConvergedReason reason; 535 PetscReal rnorm; 536 PetscInt its; 537 ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr); 538 ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr); 539 ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr); 540 ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr); 541 ierr = PCMGGetType(pc, &pcmgtype); CHKERRQ(ierr); 542 if (!test_mode || reason < 0 || rnorm > 1e-8) { 543 ierr = PetscPrintf(comm, 544 " KSP:\n" 545 " KSP Type : %s\n" 546 " KSP Convergence : %s\n" 547 " Total KSP Iterations : %D\n" 548 " Final rnorm : %e\n", 549 ksptype, KSPConvergedReasons[reason], its, 550 (double)rnorm); CHKERRQ(ierr); 551 ierr = PetscPrintf(comm, 552 " PCMG:\n" 553 " PCMG Type : %s\n" 554 " PCMG Cycle Type : %s\n", 555 PCMGTypes[pcmgtype], 556 PCMGCycleTypes[pcgmcycletype]); CHKERRQ(ierr); 557 } 558 if (!test_mode) { 559 ierr = PetscPrintf(comm," Performance:\n"); CHKERRQ(ierr); 560 } 561 { 562 PetscReal maxerror; 563 ierr = ComputeErrorMax(userO[fineLevel], operror, X[fineLevel], target, 564 &maxerror); CHKERRQ(ierr); 565 PetscReal tol = 5e-2; 566 if (!test_mode || maxerror > tol) { 567 ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm); 568 CHKERRQ(ierr); 569 ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm); 570 CHKERRQ(ierr); 571 ierr = PetscPrintf(comm, 572 " Pointwise Error (max) : %e\n" 573 " CG Solve Time : %g (%g) sec\n", 574 (double)maxerror, rt_max, rt_min); CHKERRQ(ierr); 575 } 576 } 577 if (benchmark_mode && (!test_mode)) { 578 ierr = PetscPrintf(comm, 579 " DoFs/Sec in CG : %g (%g) million\n", 580 1e-6*gsize[fineLevel]*its/rt_max, 581 1e-6*gsize[fineLevel]*its/rt_min); 582 CHKERRQ(ierr); 583 } 584 } 585 586 if (write_solution) { 587 PetscViewer vtkviewersoln; 588 589 ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr); 590 ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr); 591 ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtu"); CHKERRQ(ierr); 592 ierr = VecView(X[fineLevel], vtkviewersoln); CHKERRQ(ierr); 593 ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr); 594 } 595 596 // Cleanup 597 for (int i=0; i<numlevels; i++) { 598 ierr = VecDestroy(&X[i]); CHKERRQ(ierr); 599 ierr = VecDestroy(&Xloc[i]); CHKERRQ(ierr); 600 ierr = VecDestroy(&mult[i]); CHKERRQ(ierr); 601 ierr = VecDestroy(&userO[i]->Yloc); CHKERRQ(ierr); 602 ierr = MatDestroy(&matO[i]); CHKERRQ(ierr); 603 ierr = PetscFree(userO[i]); CHKERRQ(ierr); 604 if (i > 0) { 605 ierr = MatDestroy(&matPR[i]); CHKERRQ(ierr); 606 ierr = PetscFree(userPR[i]); CHKERRQ(ierr); 607 } 608 ierr = CeedDataDestroy(i, ceeddata[i]); CHKERRQ(ierr); 609 ierr = DMDestroy(&dm[i]); CHKERRQ(ierr); 610 } 611 ierr = PetscFree(leveldegrees); CHKERRQ(ierr); 612 ierr = PetscFree(dm); CHKERRQ(ierr); 613 ierr = PetscFree(X); CHKERRQ(ierr); 614 ierr = PetscFree(Xloc); CHKERRQ(ierr); 615 ierr = PetscFree(mult); CHKERRQ(ierr); 616 ierr = PetscFree(matO); CHKERRQ(ierr); 617 ierr = PetscFree(matPR); CHKERRQ(ierr); 618 ierr = PetscFree(ceeddata); CHKERRQ(ierr); 619 ierr = PetscFree(userO); CHKERRQ(ierr); 620 ierr = PetscFree(userPR); CHKERRQ(ierr); 621 ierr = PetscFree(lsize); CHKERRQ(ierr); 622 ierr = PetscFree(xlsize); CHKERRQ(ierr); 623 ierr = PetscFree(gsize); CHKERRQ(ierr); 624 ierr = VecDestroy(&rhs); CHKERRQ(ierr); 625 ierr = VecDestroy(&rhsloc); CHKERRQ(ierr); 626 ierr = MatDestroy(&matcoarse); CHKERRQ(ierr); 627 ierr = KSPDestroy(&ksp); CHKERRQ(ierr); 628 ierr = SNESDestroy(&snesdummy); CHKERRQ(ierr); 629 ierr = DMDestroy(&dmorig); CHKERRQ(ierr); 630 CeedVectorDestroy(&target); 631 CeedQFunctionDestroy(&qferror); 632 CeedQFunctionDestroy(&qfrestrict); 633 CeedQFunctionDestroy(&qfprolong); 634 CeedOperatorDestroy(&operror); 635 CeedDestroy(&ceed); 636 return PetscFinalize(); 637 } 638