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, diagloc; 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 " Multigrid:\n" 240 " Number of Levels : %d\n", 241 bpChoice+1, usedresource, P, Q, 242 gsize[numlevels-1]/ncompu, lsize[numlevels-1]/ncompu, 243 numlevels); CHKERRQ(ierr); 244 } 245 246 // Create RHS vector 247 ierr = VecDuplicate(Xloc[numlevels-1], &rhsloc); CHKERRQ(ierr); 248 ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr); 249 ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr); 250 CeedVectorCreate(ceed, xlsize[numlevels-1], &rhsceed); 251 CeedVectorSetArray(rhsceed, CEED_MEM_HOST, CEED_USE_POINTER, r); 252 253 // Set up libCEED operators on each level 254 ierr = PetscMalloc1(numlevels, &ceeddata); CHKERRQ(ierr); 255 for (int i=0; i<numlevels; i++) { 256 // Print level information 257 if (!test_mode && (i == 0 || i == numlevels-1)) { 258 ierr = PetscPrintf(comm," Level %D (%s):\n" 259 " Number of 1D Basis Nodes (p) : %d\n" 260 " Global Nodes : %D\n" 261 " Owned Nodes : %D\n", 262 i, (i? "fine" : "coarse"), leveldegrees[i] + 1, 263 gsize[i]/ncompu, lsize[i]/ncompu); CHKERRQ(ierr); 264 } 265 ierr = PetscMalloc1(1, &ceeddata[i]); CHKERRQ(ierr); 266 ierr = SetupLibceedByDegree(dm[i], ceed, leveldegrees[i], dim, qextra, 267 ncompu, gsize[i], xlsize[i], bpChoice, 268 ceeddata[i], i==(numlevels-1), rhsceed, 269 &target); CHKERRQ(ierr); 270 } 271 272 // Gather RHS 273 ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr); 274 ierr = VecZeroEntries(rhs); CHKERRQ(ierr); 275 ierr = DMLocalToGlobalBegin(dm[numlevels-1], rhsloc, ADD_VALUES, rhs); 276 CHKERRQ(ierr); 277 ierr = DMLocalToGlobalEnd(dm[numlevels-1], rhsloc, ADD_VALUES, rhs); 278 CHKERRQ(ierr); 279 CeedVectorDestroy(&rhsceed); 280 281 // Create the restriction/interpolation Q-function 282 CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_NONE, CEED_EVAL_INTERP, 283 &qf_restrict); 284 CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_INTERP, CEED_EVAL_NONE, 285 &qf_prolong); 286 287 // Set up libCEED level transfer operators 288 ierr = CeedLevelTransferSetup(ceed, numlevels, ncompu, bpChoice, ceeddata, 289 leveldegrees, qf_restrict, qf_prolong); 290 CHKERRQ(ierr); 291 292 // Create the error Q-function 293 CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].error, 294 bpOptions[bpChoice].errorfname, &qf_error); 295 CeedQFunctionAddInput(qf_error, "u", ncompu, CEED_EVAL_INTERP); 296 CeedQFunctionAddInput(qf_error, "true_soln", ncompu, CEED_EVAL_NONE); 297 CeedQFunctionAddOutput(qf_error, "error", ncompu, CEED_EVAL_NONE); 298 299 // Create the error operator 300 CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error); 301 CeedOperatorSetField(op_error, "u", ceeddata[numlevels-1]->Erestrictu, 302 CEED_TRANSPOSE, ceeddata[numlevels-1]->basisu, 303 CEED_VECTOR_ACTIVE); 304 CeedOperatorSetField(op_error, "true_soln", ceeddata[numlevels-1]->Erestrictui, 305 CEED_NOTRANSPOSE, CEED_BASIS_COLLOCATED, target); 306 CeedOperatorSetField(op_error, "error", ceeddata[numlevels-1]->Erestrictui, 307 CEED_NOTRANSPOSE, CEED_BASIS_COLLOCATED, 308 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 PetscScalar *petscarray; 362 CeedInt length; 363 364 ierr = VecDuplicate(X[i], &userO[i]->diag); CHKERRQ(ierr); 365 ierr = VecDuplicate(Xloc[i], &diagloc); CHKERRQ(ierr); 366 367 // -- Local diagonal 368 CeedOperatorAssembleLinearDiagonal(userO[i]->op, &diagceed, 369 CEED_REQUEST_IMMEDIATE); 370 371 // -- Copy values 372 CeedVectorGetArrayRead(diagceed, CEED_MEM_HOST, &ceedarray); 373 ierr = VecGetArray(diagloc, &petscarray); CHKERRQ(ierr); 374 CeedVectorGetLength(diagceed, &length); 375 for (CeedInt i=0; i<length; i++) 376 petscarray[i] = ceedarray[i]; 377 CeedVectorRestoreArrayRead(diagceed, &ceedarray); 378 ierr = VecRestoreArray(diagloc, &petscarray); CHKERRQ(ierr); 379 380 // -- Global diagonal 381 ierr = VecZeroEntries(userO[i]->diag); CHKERRQ(ierr); 382 ierr = DMLocalToGlobalBegin(userO[i]->dm, diagloc, ADD_VALUES, 383 userO[i]->diag); CHKERRQ(ierr); 384 ierr = DMLocalToGlobalEnd(userO[i]->dm, diagloc, ADD_VALUES, 385 userO[i]->diag); CHKERRQ(ierr); 386 387 // -- Cleanup 388 ierr = VecDestroy(&diagloc); CHKERRQ(ierr); 389 CeedVectorDestroy(&diagceed); 390 391 if (i > 0) { 392 // Interp Operator 393 userI[i]->comm = comm; 394 userI[i]->dmc = dm[i-1]; 395 userI[i]->dmf = dm[i]; 396 userI[i]->Xloc = Xloc[i-1]; 397 userI[i]->Yloc = userO[i]->Yloc; 398 userI[i]->mult = mult[i]; 399 userI[i]->ceedvecc = userO[i-1]->xceed; 400 userI[i]->ceedvecf = userO[i]->yceed; 401 userI[i]->op = ceeddata[i]->op_interp; 402 userI[i]->ceed = ceed; 403 404 // Restrict Operator 405 userR[i]->comm = comm; 406 userR[i]->dmc = dm[i-1]; 407 userR[i]->dmf = dm[i]; 408 userR[i]->Xloc = Xloc[i]; 409 userR[i]->Yloc = userO[i-1]->Yloc; 410 userR[i]->mult = mult[i]; 411 userR[i]->ceedvecf = userO[i]->xceed; 412 userR[i]->ceedvecc = userO[i-1]->yceed; 413 userR[i]->op = ceeddata[i]->op_restrict; 414 userR[i]->ceed = ceed; 415 } 416 } 417 418 // Set up KSP 419 ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr); 420 { 421 ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 422 ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 423 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 424 PETSC_DEFAULT); CHKERRQ(ierr); 425 } 426 ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr); 427 ierr = KSPSetOperators(ksp, matO[numlevels-1], matO[numlevels-1]); 428 CHKERRQ(ierr); 429 430 // Set up PCMG 431 ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 432 PCMGCycleType pcgmcycletype = PC_MG_CYCLE_V; 433 { 434 ierr = PCSetType(pc, PCMG); CHKERRQ(ierr); 435 436 // PCMG levels 437 ierr = PCMGSetLevels(pc, numlevels, NULL); CHKERRQ(ierr); 438 for (int i=0; i<numlevels; i++) { 439 // Smoother 440 KSP smoother; 441 PC smoother_pc; 442 ierr = PCMGGetSmoother(pc, i, &smoother); CHKERRQ(ierr); 443 ierr = KSPSetType(smoother, KSPCHEBYSHEV); CHKERRQ(ierr); 444 ierr = KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1); CHKERRQ(ierr); 445 ierr = KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE); CHKERRQ(ierr); 446 ierr = KSPSetOperators(smoother, matO[i], matO[i]); CHKERRQ(ierr); 447 ierr = KSPGetPC(smoother, &smoother_pc); CHKERRQ(ierr); 448 ierr = PCSetType(smoother_pc, PCJACOBI); CHKERRQ(ierr); 449 ierr = PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 450 451 // Work vector 452 if (i < numlevels-1) { 453 ierr = PCMGSetX(pc, i, X[i]); CHKERRQ(ierr); 454 } 455 456 // Level transfers 457 if (i > 0) { 458 // Interpolation 459 ierr = PCMGSetInterpolation(pc, i, matI[i]); CHKERRQ(ierr); 460 461 // Restriction 462 ierr = PCMGSetRestriction(pc, i, matR[i]); CHKERRQ(ierr); 463 } 464 465 // Coarse solve 466 KSP coarse; 467 PC coarse_pc; 468 ierr = PCMGGetCoarseSolve(pc, &coarse); CHKERRQ(ierr); 469 ierr = KSPSetType(coarse, KSPCG); CHKERRQ(ierr); 470 ierr = KSPSetOperators(coarse, matO[0], matO[0]); CHKERRQ(ierr); 471 ierr = KSPSetTolerances(coarse, 1e-10, 1e-10, PETSC_DEFAULT, 472 PETSC_DEFAULT); CHKERRQ(ierr); 473 ierr = KSPGetPC(coarse, &coarse_pc); CHKERRQ(ierr); 474 ierr = PCSetType(coarse_pc, PCJACOBI); CHKERRQ(ierr); 475 ierr = PCJacobiSetType(coarse_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 476 } 477 478 // PCMG options 479 ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr); 480 ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr); 481 ierr = PCMGSetCycleType(pc, pcgmcycletype); CHKERRQ(ierr); 482 } 483 484 // First run, if benchmarking 485 if (benchmark_mode) { 486 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1); 487 CHKERRQ(ierr); 488 ierr = VecZeroEntries(X[numlevels-1]); CHKERRQ(ierr); 489 my_rt_start = MPI_Wtime(); 490 ierr = KSPSolve(ksp, rhs, X[numlevels-1]); CHKERRQ(ierr); 491 my_rt = MPI_Wtime() - my_rt_start; 492 ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm); 493 CHKERRQ(ierr); 494 // Set maxits based on first iteration timing 495 if (my_rt > 0.02) { 496 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5); 497 CHKERRQ(ierr); 498 } else { 499 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20); 500 CHKERRQ(ierr); 501 } 502 } 503 504 // Timed solve 505 ierr = VecZeroEntries(X[numlevels-1]); CHKERRQ(ierr); 506 ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr); 507 my_rt_start = MPI_Wtime(); 508 ierr = KSPSolve(ksp, rhs, X[numlevels-1]); CHKERRQ(ierr); 509 my_rt = MPI_Wtime() - my_rt_start; 510 511 // Output results 512 { 513 KSPType ksptype; 514 PCMGType pcmgtype; 515 KSPConvergedReason reason; 516 PetscReal rnorm; 517 PetscInt its; 518 ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr); 519 ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr); 520 ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr); 521 ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr); 522 ierr = PCMGGetType(pc, &pcmgtype); CHKERRQ(ierr); 523 if (!test_mode || reason < 0 || rnorm > 1e-8) { 524 ierr = PetscPrintf(comm, 525 " KSP:\n" 526 " KSP Type : %s\n" 527 " KSP Convergence : %s\n" 528 " Total KSP Iterations : %D\n" 529 " Final rnorm : %e\n", 530 ksptype, KSPConvergedReasons[reason], its, 531 (double)rnorm); CHKERRQ(ierr); 532 ierr = PetscPrintf(comm, 533 " PCMG:\n" 534 " PCMG Type : %s\n" 535 " PCMG Cycle Type : %s\n", 536 PCMGTypes[pcmgtype], 537 PCMGCycleTypes[pcgmcycletype]); CHKERRQ(ierr); 538 } 539 if (!test_mode) { 540 ierr = PetscPrintf(comm," Performance:\n"); CHKERRQ(ierr); 541 } 542 { 543 PetscReal maxerror; 544 ierr = ComputeErrorMax(userO[numlevels-1], op_error, X[numlevels-1], target, 545 &maxerror); CHKERRQ(ierr); 546 PetscReal tol = 5e-2; 547 if (!test_mode || maxerror > tol) { 548 ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm); 549 CHKERRQ(ierr); 550 ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm); 551 CHKERRQ(ierr); 552 ierr = PetscPrintf(comm, 553 " Pointwise Error (max) : %e\n" 554 " CG Solve Time : %g (%g) sec\n", 555 (double)maxerror, rt_max, rt_min); CHKERRQ(ierr); 556 } 557 } 558 if (benchmark_mode && (!test_mode)) { 559 ierr = PetscPrintf(comm, 560 " DoFs/Sec in CG : %g (%g) million\n", 561 1e-6*gsize[numlevels-1]*its/rt_max, 562 1e-6*gsize[numlevels-1]*its/rt_min); 563 CHKERRQ(ierr); 564 } 565 } 566 567 if (write_solution) { 568 PetscViewer vtkviewersoln; 569 570 ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr); 571 ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr); 572 ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr); 573 ierr = VecView(X[numlevels-1], vtkviewersoln); CHKERRQ(ierr); 574 ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr); 575 } 576 577 // Cleanup 578 for (int i=0; i<numlevels; i++) { 579 ierr = VecDestroy(&X[i]); CHKERRQ(ierr); 580 ierr = VecDestroy(&Xloc[i]); CHKERRQ(ierr); 581 ierr = VecDestroy(&mult[i]); CHKERRQ(ierr); 582 ierr = VecDestroy(&userO[i]->Yloc); CHKERRQ(ierr); 583 ierr = VecDestroy(&userO[i]->diag); CHKERRQ(ierr); 584 ierr = MatDestroy(&matO[i]); CHKERRQ(ierr); 585 ierr = PetscFree(userO[i]); CHKERRQ(ierr); 586 if (i > 0) { 587 ierr = MatDestroy(&matI[i]); CHKERRQ(ierr); 588 ierr = PetscFree(userI[i]); CHKERRQ(ierr); 589 ierr = MatDestroy(&matR[i]); CHKERRQ(ierr); 590 ierr = PetscFree(userR[i]); CHKERRQ(ierr); 591 } 592 ierr = CeedDataDestroy(i, ceeddata[i]); CHKERRQ(ierr); 593 ierr = DMDestroy(&dm[i]); CHKERRQ(ierr); 594 } 595 ierr = PetscFree(leveldegrees); CHKERRQ(ierr); 596 ierr = PetscFree(dm); CHKERRQ(ierr); 597 ierr = PetscFree(X); CHKERRQ(ierr); 598 ierr = PetscFree(Xloc); CHKERRQ(ierr); 599 ierr = PetscFree(mult); CHKERRQ(ierr); 600 ierr = PetscFree(matO); CHKERRQ(ierr); 601 ierr = PetscFree(matI); CHKERRQ(ierr); 602 ierr = PetscFree(matR); CHKERRQ(ierr); 603 ierr = PetscFree(ceeddata); CHKERRQ(ierr); 604 ierr = PetscFree(userO); CHKERRQ(ierr); 605 ierr = PetscFree(userI); CHKERRQ(ierr); 606 ierr = PetscFree(userR); CHKERRQ(ierr); 607 ierr = PetscFree(lsize); CHKERRQ(ierr); 608 ierr = PetscFree(xlsize); CHKERRQ(ierr); 609 ierr = PetscFree(gsize); CHKERRQ(ierr); 610 ierr = VecDestroy(&rhs); CHKERRQ(ierr); 611 ierr = VecDestroy(&rhsloc); CHKERRQ(ierr); 612 ierr = KSPDestroy(&ksp); CHKERRQ(ierr); 613 ierr = DMDestroy(&dmOrig); CHKERRQ(ierr); 614 CeedVectorDestroy(&target); 615 CeedQFunctionDestroy(&qf_error); 616 CeedQFunctionDestroy(&qf_restrict); 617 CeedQFunctionDestroy(&qf_prolong); 618 CeedOperatorDestroy(&op_error); 619 CeedDestroy(&ceed); 620 return PetscFinalize(); 621 } 622