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