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