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")) vectype = VECHIP; 160 else vectype = VECSTANDARD; 161 } 162 } 163 ierr = DMSetVecType(dmorig, vectype); CHKERRQ(ierr); 164 ierr = DMSetFromOptions(dmorig); CHKERRQ(ierr); 165 166 // Allocate arrays for PETSc objects for each level 167 switch (coarsen) { 168 case COARSEN_UNIFORM: 169 numlevels = degree; 170 break; 171 case COARSEN_LOGARITHMIC: 172 numlevels = ceil(log(degree)/log(2)) + 1; 173 break; 174 } 175 ierr = PetscMalloc1(numlevels, &leveldegrees); CHKERRQ(ierr); 176 fineLevel = numlevels - 1; 177 178 switch (coarsen) { 179 case COARSEN_UNIFORM: 180 for (int i=0; i<numlevels; i++) leveldegrees[i] = i + 1; 181 break; 182 case COARSEN_LOGARITHMIC: 183 for (int i=0; i<numlevels - 1; i++) leveldegrees[i] = pow(2,i); 184 leveldegrees[fineLevel] = degree; 185 break; 186 } 187 ierr = PetscMalloc1(numlevels, &dm); CHKERRQ(ierr); 188 ierr = PetscMalloc1(numlevels, &X); CHKERRQ(ierr); 189 ierr = PetscMalloc1(numlevels, &Xloc); CHKERRQ(ierr); 190 ierr = PetscMalloc1(numlevels, &mult); CHKERRQ(ierr); 191 ierr = PetscMalloc1(numlevels, &userO); CHKERRQ(ierr); 192 ierr = PetscMalloc1(numlevels, &userPR); CHKERRQ(ierr); 193 ierr = PetscMalloc1(numlevels, &matO); CHKERRQ(ierr); 194 ierr = PetscMalloc1(numlevels, &matPR); CHKERRQ(ierr); 195 ierr = PetscMalloc1(numlevels, &lsize); CHKERRQ(ierr); 196 ierr = PetscMalloc1(numlevels, &xlsize); CHKERRQ(ierr); 197 ierr = PetscMalloc1(numlevels, &gsize); CHKERRQ(ierr); 198 199 // Setup DM and Operator Mat Shells for each level 200 for (CeedInt i=0; i<numlevels; i++) { 201 // Create DM 202 ierr = DMClone(dmorig, &dm[i]); CHKERRQ(ierr); 203 ierr = DMGetVecType(dmorig, &vectype); CHKERRQ(ierr); 204 ierr = DMSetVecType(dm[i], vectype); CHKERRQ(ierr); 205 ierr = SetupDMByDegree(dm[i], leveldegrees[i], ncompu, bpchoice); 206 CHKERRQ(ierr); 207 208 // Create vectors 209 ierr = DMCreateGlobalVector(dm[i], &X[i]); CHKERRQ(ierr); 210 ierr = VecGetLocalSize(X[i], &lsize[i]); CHKERRQ(ierr); 211 ierr = VecGetSize(X[i], &gsize[i]); CHKERRQ(ierr); 212 ierr = DMCreateLocalVector(dm[i], &Xloc[i]); CHKERRQ(ierr); 213 ierr = VecGetSize(Xloc[i], &xlsize[i]); CHKERRQ(ierr); 214 215 // Operator 216 ierr = PetscMalloc1(1, &userO[i]); CHKERRQ(ierr); 217 ierr = MatCreateShell(comm, lsize[i], lsize[i], gsize[i], gsize[i], 218 userO[i], &matO[i]); CHKERRQ(ierr); 219 ierr = MatShellSetOperation(matO[i], MATOP_MULT, 220 (void(*)(void))MatMult_Ceed); CHKERRQ(ierr); 221 ierr = MatShellSetOperation(matO[i], MATOP_GET_DIAGONAL, 222 (void(*)(void))MatGetDiag); CHKERRQ(ierr); 223 ierr = MatShellSetVecType(matO[i], vectype); CHKERRQ(ierr); 224 225 // Level transfers 226 if (i > 0) { 227 // Interp 228 ierr = PetscMalloc1(1, &userPR[i]); CHKERRQ(ierr); 229 ierr = MatCreateShell(comm, lsize[i], lsize[i-1], gsize[i], gsize[i-1], 230 userPR[i], &matPR[i]); CHKERRQ(ierr); 231 ierr = MatShellSetOperation(matPR[i], MATOP_MULT, 232 (void(*)(void))MatMult_Prolong); 233 CHKERRQ(ierr); 234 ierr = MatShellSetOperation(matPR[i], MATOP_MULT_TRANSPOSE, 235 (void(*)(void))MatMult_Restrict); 236 CHKERRQ(ierr); 237 ierr = MatShellSetVecType(matPR[i], vectype); CHKERRQ(ierr); 238 } 239 } 240 ierr = VecDuplicate(X[fineLevel], &rhs); CHKERRQ(ierr); 241 242 // Print global grid information 243 if (!test_mode) { 244 PetscInt P = degree + 1, Q = P + qextra; 245 246 const char *usedresource; 247 CeedGetResource(ceed, &usedresource); 248 249 ierr = VecGetType(X[0], &vectype); CHKERRQ(ierr); 250 251 ierr = PetscPrintf(comm, 252 "\n-- CEED Benchmark Problem %d -- libCEED + PETSc + PCMG --\n" 253 " PETSc:\n" 254 " PETSc Vec Type : %s\n" 255 " libCEED:\n" 256 " libCEED Backend : %s\n" 257 " libCEED Backend MemType : %s\n" 258 " Mesh:\n" 259 " Number of 1D Basis Nodes (p) : %d\n" 260 " Number of 1D Quadrature Points (q) : %d\n" 261 " Global Nodes : %D\n" 262 " Owned Nodes : %D\n" 263 " DoF per node : %D\n" 264 " Multigrid:\n" 265 " Number of Levels : %d\n", 266 bpchoice+1, vectype, usedresource, 267 CeedMemTypes[memtypebackend], 268 P, Q, gsize[fineLevel]/ncompu, lsize[fineLevel]/ncompu, 269 ncompu, numlevels); CHKERRQ(ierr); 270 } 271 272 // Create RHS vector 273 ierr = VecDuplicate(Xloc[fineLevel], &rhsloc); CHKERRQ(ierr); 274 ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr); 275 ierr = VecGetArrayAndMemType(rhsloc, &r, &memtype); CHKERRQ(ierr); 276 CeedVectorCreate(ceed, xlsize[fineLevel], &rhsceed); 277 CeedVectorSetArray(rhsceed, MemTypeP2C(memtype), CEED_USE_POINTER, r); 278 279 // Set up libCEED operators on each level 280 ierr = PetscMalloc1(numlevels, &ceeddata); CHKERRQ(ierr); 281 for (int i=0; i<numlevels; i++) { 282 // Print level information 283 if (!test_mode && (i == 0 || i == fineLevel)) { 284 ierr = PetscPrintf(comm," Level %D (%s):\n" 285 " Number of 1D Basis Nodes (p) : %d\n" 286 " Global Nodes : %D\n" 287 " Owned Nodes : %D\n", 288 i, (i? "fine" : "coarse"), leveldegrees[i] + 1, 289 gsize[i]/ncompu, lsize[i]/ncompu); CHKERRQ(ierr); 290 } 291 ierr = PetscMalloc1(1, &ceeddata[i]); CHKERRQ(ierr); 292 ierr = SetupLibceedByDegree(dm[i], ceed, leveldegrees[i], dim, qextra, 293 ncompu, gsize[i], xlsize[i], bpchoice, 294 ceeddata[i], i==(fineLevel), rhsceed, &target); 295 CHKERRQ(ierr); 296 } 297 298 // Gather RHS 299 CeedVectorTakeArray(rhsceed, MemTypeP2C(memtype), NULL); 300 ierr = VecRestoreArrayAndMemType(rhsloc, &r); CHKERRQ(ierr); 301 ierr = VecZeroEntries(rhs); CHKERRQ(ierr); 302 ierr = DMLocalToGlobal(dm[fineLevel], rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr); 303 CeedVectorDestroy(&rhsceed); 304 305 // Create the restriction/interpolation QFunction 306 CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_NONE, CEED_EVAL_INTERP, 307 &qfrestrict); 308 CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_INTERP, CEED_EVAL_NONE, 309 &qfprolong); 310 311 // Set up libCEED level transfer operators 312 ierr = CeedLevelTransferSetup(ceed, numlevels, ncompu, bpchoice, ceeddata, 313 leveldegrees, qfrestrict, qfprolong); 314 CHKERRQ(ierr); 315 316 // Create the error QFunction 317 CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].error, 318 bpOptions[bpchoice].errorfname, &qferror); 319 CeedQFunctionAddInput(qferror, "u", ncompu, CEED_EVAL_INTERP); 320 CeedQFunctionAddInput(qferror, "true_soln", ncompu, CEED_EVAL_NONE); 321 CeedQFunctionAddOutput(qferror, "error", ncompu, CEED_EVAL_NONE); 322 323 // Create the error operator 324 CeedOperatorCreate(ceed, qferror, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, 325 &operror); 326 CeedOperatorSetField(operror, "u", ceeddata[fineLevel]->Erestrictu, 327 ceeddata[fineLevel]->basisu, CEED_VECTOR_ACTIVE); 328 CeedOperatorSetField(operror, "true_soln", ceeddata[fineLevel]->Erestrictui, 329 CEED_BASIS_COLLOCATED, target); 330 CeedOperatorSetField(operror, "error", ceeddata[fineLevel]->Erestrictui, 331 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 332 333 // Calculate multiplicity 334 for (int i=0; i<numlevels; i++) { 335 PetscScalar *x; 336 337 // CEED vector 338 ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr); 339 ierr = VecGetArray(Xloc[i], &x); CHKERRQ(ierr); 340 CeedVectorSetArray(ceeddata[i]->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x); 341 342 // Multiplicity 343 CeedElemRestrictionGetMultiplicity(ceeddata[i]->Erestrictu, 344 ceeddata[i]->xceed); 345 CeedVectorSyncArray(ceeddata[i]->xceed, CEED_MEM_HOST); 346 347 // Restore vector 348 ierr = VecRestoreArray(Xloc[i], &x); CHKERRQ(ierr); 349 350 // Creat mult vector 351 ierr = VecDuplicate(Xloc[i], &mult[i]); CHKERRQ(ierr); 352 353 // Local-to-global 354 ierr = VecZeroEntries(X[i]); CHKERRQ(ierr); 355 ierr = DMLocalToGlobal(dm[i], Xloc[i], ADD_VALUES, X[i]); 356 CHKERRQ(ierr); 357 ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr); 358 359 // Global-to-local 360 ierr = DMGlobalToLocal(dm[i], X[i], INSERT_VALUES, mult[i]); 361 CHKERRQ(ierr); 362 ierr = VecZeroEntries(X[i]); CHKERRQ(ierr); 363 364 // Multiplicity scaling 365 ierr = VecReciprocal(mult[i]); 366 } 367 368 // Set up Mat 369 for (int i=0; i<numlevels; i++) { 370 // User Operator 371 userO[i]->comm = comm; 372 userO[i]->dm = dm[i]; 373 userO[i]->Xloc = Xloc[i]; 374 ierr = VecDuplicate(Xloc[i], &userO[i]->Yloc); CHKERRQ(ierr); 375 userO[i]->xceed = ceeddata[i]->xceed; 376 userO[i]->yceed = ceeddata[i]->yceed; 377 userO[i]->op = ceeddata[i]->opapply; 378 userO[i]->ceed = ceed; 379 380 if (i > 0) { 381 // Prolongation/Restriction Operator 382 userPR[i]->comm = comm; 383 userPR[i]->dmf = dm[i]; 384 userPR[i]->dmc = dm[i-1]; 385 userPR[i]->locvecc = Xloc[i-1]; 386 userPR[i]->locvecf = userO[i]->Yloc; 387 userPR[i]->multvec = mult[i]; 388 userPR[i]->ceedvecc = userO[i-1]->xceed; 389 userPR[i]->ceedvecf = userO[i]->yceed; 390 userPR[i]->opprolong = ceeddata[i]->opprolong; 391 userPR[i]->oprestrict = ceeddata[i]->oprestrict; 392 userPR[i]->ceed = ceed; 393 } 394 } 395 396 // Setup dummy SNES for AMG coarse solve 397 ierr = SNESCreate(comm, &snesdummy); CHKERRQ(ierr); 398 ierr = SNESSetDM(snesdummy, dm[0]); CHKERRQ(ierr); 399 ierr = SNESSetSolution(snesdummy, X[0]); CHKERRQ(ierr); 400 401 // -- Jacobian matrix 402 ierr = DMSetMatType(dm[0], MATAIJ); CHKERRQ(ierr); 403 ierr = DMCreateMatrix(dm[0], &matcoarse); CHKERRQ(ierr); 404 ierr = SNESSetJacobian(snesdummy, matcoarse, matcoarse, NULL, 405 NULL); CHKERRQ(ierr); 406 407 // -- Residual evaluation function 408 ierr = SNESSetFunction(snesdummy, X[0], FormResidual_Ceed, 409 userO[0]); CHKERRQ(ierr); 410 411 // -- Form Jacobian 412 ierr = SNESComputeJacobianDefaultColor(snesdummy, X[0], matO[0], 413 matcoarse, NULL); CHKERRQ(ierr); 414 415 // Set up KSP 416 ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr); 417 { 418 ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 419 ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 420 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 421 PETSC_DEFAULT); CHKERRQ(ierr); 422 } 423 ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr); 424 ierr = KSPSetOperators(ksp, matO[fineLevel], matO[fineLevel]); 425 CHKERRQ(ierr); 426 427 // Set up PCMG 428 ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 429 PCMGCycleType pcgmcycletype = PC_MG_CYCLE_V; 430 { 431 ierr = PCSetType(pc, PCMG); CHKERRQ(ierr); 432 433 // PCMG levels 434 ierr = PCMGSetLevels(pc, numlevels, NULL); CHKERRQ(ierr); 435 for (int i=0; i<numlevels; i++) { 436 // Smoother 437 KSP smoother; 438 PC smoother_pc; 439 ierr = PCMGGetSmoother(pc, i, &smoother); CHKERRQ(ierr); 440 ierr = KSPSetType(smoother, KSPCHEBYSHEV); CHKERRQ(ierr); 441 ierr = KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1); CHKERRQ(ierr); 442 ierr = KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE); CHKERRQ(ierr); 443 ierr = KSPSetOperators(smoother, matO[i], matO[i]); CHKERRQ(ierr); 444 ierr = KSPGetPC(smoother, &smoother_pc); CHKERRQ(ierr); 445 ierr = PCSetType(smoother_pc, PCJACOBI); CHKERRQ(ierr); 446 ierr = PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 447 448 // Work vector 449 if (i < numlevels - 1) { 450 ierr = PCMGSetX(pc, i, X[i]); CHKERRQ(ierr); 451 } 452 453 // Level transfers 454 if (i > 0) { 455 // Interpolation 456 ierr = PCMGSetInterpolation(pc, i, matPR[i]); CHKERRQ(ierr); 457 } 458 459 // Coarse solve 460 KSP coarse; 461 PC coarse_pc; 462 ierr = PCMGGetCoarseSolve(pc, &coarse); CHKERRQ(ierr); 463 ierr = KSPSetType(coarse, KSPPREONLY); CHKERRQ(ierr); 464 ierr = KSPSetOperators(coarse, matcoarse, matcoarse); CHKERRQ(ierr); 465 466 ierr = KSPGetPC(coarse, &coarse_pc); CHKERRQ(ierr); 467 ierr = PCSetType(coarse_pc, PCGAMG); CHKERRQ(ierr); 468 469 ierr = KSPSetOptionsPrefix(coarse, "coarse_"); CHKERRQ(ierr); 470 ierr = PCSetOptionsPrefix(coarse_pc, "coarse_"); CHKERRQ(ierr); 471 ierr = KSPSetFromOptions(coarse); CHKERRQ(ierr); 472 ierr = PCSetFromOptions(coarse_pc); CHKERRQ(ierr); 473 } 474 475 // PCMG options 476 ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr); 477 ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr); 478 ierr = PCMGSetCycleType(pc, pcgmcycletype); CHKERRQ(ierr); 479 } 480 481 // First run, if benchmarking 482 if (benchmark_mode) { 483 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1); 484 CHKERRQ(ierr); 485 ierr = VecZeroEntries(X[fineLevel]); CHKERRQ(ierr); 486 my_rt_start = MPI_Wtime(); 487 ierr = KSPSolve(ksp, rhs, X[fineLevel]); CHKERRQ(ierr); 488 my_rt = MPI_Wtime() - my_rt_start; 489 ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm); 490 CHKERRQ(ierr); 491 // Set maxits based on first iteration timing 492 if (my_rt > 0.02) { 493 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5); 494 CHKERRQ(ierr); 495 } else { 496 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20); 497 CHKERRQ(ierr); 498 } 499 } 500 501 // Timed solve 502 ierr = VecZeroEntries(X[fineLevel]); CHKERRQ(ierr); 503 ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr); 504 505 // -- Performance logging 506 ierr = PetscLogStageRegister("Solve Stage", &solvestage); CHKERRQ(ierr); 507 ierr = PetscLogStagePush(solvestage); CHKERRQ(ierr); 508 509 // -- Solve 510 my_rt_start = MPI_Wtime(); 511 ierr = KSPSolve(ksp, rhs, X[fineLevel]); CHKERRQ(ierr); 512 my_rt = MPI_Wtime() - my_rt_start; 513 514 515 // -- Performance logging 516 ierr = PetscLogStagePop(); 517 518 // Output results 519 { 520 KSPType ksptype; 521 PCMGType pcmgtype; 522 KSPConvergedReason reason; 523 PetscReal rnorm; 524 PetscInt its; 525 ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr); 526 ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr); 527 ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr); 528 ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr); 529 ierr = PCMGGetType(pc, &pcmgtype); CHKERRQ(ierr); 530 if (!test_mode || reason < 0 || rnorm > 1e-8) { 531 ierr = PetscPrintf(comm, 532 " KSP:\n" 533 " KSP Type : %s\n" 534 " KSP Convergence : %s\n" 535 " Total KSP Iterations : %D\n" 536 " Final rnorm : %e\n", 537 ksptype, KSPConvergedReasons[reason], its, 538 (double)rnorm); CHKERRQ(ierr); 539 ierr = PetscPrintf(comm, 540 " PCMG:\n" 541 " PCMG Type : %s\n" 542 " PCMG Cycle Type : %s\n", 543 PCMGTypes[pcmgtype], 544 PCMGCycleTypes[pcgmcycletype]); CHKERRQ(ierr); 545 } 546 if (!test_mode) { 547 ierr = PetscPrintf(comm," Performance:\n"); CHKERRQ(ierr); 548 } 549 { 550 PetscReal maxerror; 551 ierr = ComputeErrorMax(userO[fineLevel], operror, X[fineLevel], target, 552 &maxerror); CHKERRQ(ierr); 553 PetscReal tol = 5e-2; 554 if (!test_mode || maxerror > tol) { 555 ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm); 556 CHKERRQ(ierr); 557 ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm); 558 CHKERRQ(ierr); 559 ierr = PetscPrintf(comm, 560 " Pointwise Error (max) : %e\n" 561 " CG Solve Time : %g (%g) sec\n", 562 (double)maxerror, rt_max, rt_min); CHKERRQ(ierr); 563 } 564 } 565 if (benchmark_mode && (!test_mode)) { 566 ierr = PetscPrintf(comm, 567 " DoFs/Sec in CG : %g (%g) million\n", 568 1e-6*gsize[fineLevel]*its/rt_max, 569 1e-6*gsize[fineLevel]*its/rt_min); 570 CHKERRQ(ierr); 571 } 572 } 573 574 if (write_solution) { 575 PetscViewer vtkviewersoln; 576 577 ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr); 578 ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr); 579 ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr); 580 ierr = VecView(X[fineLevel], vtkviewersoln); CHKERRQ(ierr); 581 ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr); 582 } 583 584 // Cleanup 585 for (int i=0; i<numlevels; i++) { 586 ierr = VecDestroy(&X[i]); CHKERRQ(ierr); 587 ierr = VecDestroy(&Xloc[i]); CHKERRQ(ierr); 588 ierr = VecDestroy(&mult[i]); CHKERRQ(ierr); 589 ierr = VecDestroy(&userO[i]->Yloc); CHKERRQ(ierr); 590 ierr = MatDestroy(&matO[i]); CHKERRQ(ierr); 591 ierr = PetscFree(userO[i]); CHKERRQ(ierr); 592 if (i > 0) { 593 ierr = MatDestroy(&matPR[i]); CHKERRQ(ierr); 594 ierr = PetscFree(userPR[i]); CHKERRQ(ierr); 595 } 596 ierr = CeedDataDestroy(i, ceeddata[i]); CHKERRQ(ierr); 597 ierr = DMDestroy(&dm[i]); CHKERRQ(ierr); 598 } 599 ierr = PetscFree(leveldegrees); CHKERRQ(ierr); 600 ierr = PetscFree(dm); CHKERRQ(ierr); 601 ierr = PetscFree(X); CHKERRQ(ierr); 602 ierr = PetscFree(Xloc); CHKERRQ(ierr); 603 ierr = PetscFree(mult); CHKERRQ(ierr); 604 ierr = PetscFree(matO); CHKERRQ(ierr); 605 ierr = PetscFree(matPR); CHKERRQ(ierr); 606 ierr = PetscFree(ceeddata); CHKERRQ(ierr); 607 ierr = PetscFree(userO); CHKERRQ(ierr); 608 ierr = PetscFree(userPR); CHKERRQ(ierr); 609 ierr = PetscFree(lsize); CHKERRQ(ierr); 610 ierr = PetscFree(xlsize); CHKERRQ(ierr); 611 ierr = PetscFree(gsize); CHKERRQ(ierr); 612 ierr = VecDestroy(&rhs); CHKERRQ(ierr); 613 ierr = VecDestroy(&rhsloc); CHKERRQ(ierr); 614 ierr = MatDestroy(&matcoarse); CHKERRQ(ierr); 615 ierr = KSPDestroy(&ksp); CHKERRQ(ierr); 616 ierr = SNESDestroy(&snesdummy); CHKERRQ(ierr); 617 ierr = DMDestroy(&dmorig); CHKERRQ(ierr); 618 CeedVectorDestroy(&target); 619 CeedQFunctionDestroy(&qferror); 620 CeedQFunctionDestroy(&qfrestrict); 621 CeedQFunctionDestroy(&qfprolong); 622 CeedOperatorDestroy(&operror); 623 CeedDestroy(&ceed); 624 return PetscFinalize(); 625 } 626