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 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 is intentionally "raw", using only low-level communication 23 // primitives. 24 // 25 // Build with: 26 // 27 // make bpsraw [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] 28 // 29 // Sample runs: 30 // 31 // ./bpsraw -problem bp1 32 // ./bpsraw -problem bp2 -ceed /cpu/self 33 // ./bpsraw -problem bp3 -ceed /gpu/occa 34 // ./bpsraw -problem bp4 -ceed /cpu/occa 35 // ./bpsraw -problem bp5 -ceed /omp/occa 36 // ./bpsraw -problem bp6 -ceed /ocl/occa 37 // 38 //TESTARGS -ceed {ceed_resource} -test -problem bp2 -degree 5 -qextra 5 39 40 /// @file 41 /// CEED BPs example using PETSc 42 /// See bps.c for an implementation using DMPlex unstructured grids. 43 const char help[] = "Solve CEED BPs using PETSc\n"; 44 45 #include <stdbool.h> 46 #include <string.h> 47 #include <petscksp.h> 48 #include <ceed.h> 49 #include "qfunctions/bps/common.h" 50 #include "qfunctions/bps/bp1.h" 51 #include "qfunctions/bps/bp2.h" 52 #include "qfunctions/bps/bp3.h" 53 #include "qfunctions/bps/bp4.h" 54 55 #include <petscsys.h> 56 #if PETSC_VERSION_LT(3,12,0) 57 #ifdef PETSC_HAVE_CUDA 58 #include <petsccuda.h> 59 // Note: With PETSc prior to version 3.12.0, providing the source path to 60 // include 'cublas_v2.h' will be needed to use 'petsccuda.h'. 61 #endif 62 #endif 63 64 static void Split3(PetscInt size, PetscInt m[3], bool reverse) { 65 for (PetscInt d=0,sizeleft=size; d<3; d++) { 66 PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(sizeleft, 1./(3 - d))); 67 while (try * (sizeleft / try) != sizeleft) try++; 68 m[reverse ? 2-d : d] = try; 69 sizeleft /= try; 70 } 71 } 72 73 static PetscInt Max3(const PetscInt a[3]) { 74 return PetscMax(a[0], PetscMax(a[1], a[2])); 75 } 76 static PetscInt Min3(const PetscInt a[3]) { 77 return PetscMin(a[0], PetscMin(a[1], a[2])); 78 } 79 static void GlobalNodes(const PetscInt p[3], const PetscInt irank[3], 80 PetscInt degree, const PetscInt melem[3], 81 PetscInt mnodes[3]) { 82 for (int d=0; d<3; d++) 83 mnodes[d] = degree*melem[d] + (irank[d] == p[d]-1); 84 } 85 static PetscInt GlobalStart(const PetscInt p[3], const PetscInt irank[3], 86 PetscInt degree, const PetscInt melem[3]) { 87 PetscInt start = 0; 88 // Dumb brute-force is easier to read 89 for (PetscInt i=0; i<p[0]; i++) { 90 for (PetscInt j=0; j<p[1]; j++) { 91 for (PetscInt k=0; k<p[2]; k++) { 92 PetscInt mnodes[3], ijkrank[] = {i,j,k}; 93 if (i == irank[0] && j == irank[1] && k == irank[2]) return start; 94 GlobalNodes(p, ijkrank, degree, melem, mnodes); 95 start += mnodes[0] * mnodes[1] * mnodes[2]; 96 } 97 } 98 } 99 return -1; 100 } 101 static int CreateRestriction(Ceed ceed, CeedInterlaceMode imode, 102 const CeedInt melem[3], CeedInt P, CeedInt ncomp, 103 CeedElemRestriction *Erestrict) { 104 const PetscInt nelem = melem[0]*melem[1]*melem[2]; 105 PetscInt mnodes[3], *idx, *idxp; 106 107 // Get indicies 108 for (int d=0; d<3; d++) mnodes[d] = melem[d]*(P-1) + 1; 109 idxp = idx = malloc(nelem*P*P*P*sizeof idx[0]); 110 for (CeedInt i=0; i<melem[0]; i++) 111 for (CeedInt j=0; j<melem[1]; j++) 112 for (CeedInt k=0; k<melem[2]; k++,idxp += P*P*P) 113 for (CeedInt ii=0; ii<P; ii++) 114 for (CeedInt jj=0; jj<P; jj++) 115 for (CeedInt kk=0; kk<P; kk++) { 116 if (0) { // This is the C-style (i,j,k) ordering that I prefer 117 idxp[(ii*P+jj)*P+kk] = (((i*(P-1)+ii)*mnodes[1] 118 + (j*(P-1)+jj))*mnodes[2] 119 + (k*(P-1)+kk)); 120 } else { // (k,j,i) ordering for consistency with MFEM example 121 idxp[ii+P*(jj+P*kk)] = (((i*(P-1)+ii)*mnodes[1] 122 + (j*(P-1)+jj))*mnodes[2] 123 + (k*(P-1)+kk)); 124 } 125 } 126 127 // Setup CEED restriction 128 CeedElemRestrictionCreate(ceed, imode, nelem, P*P*P, 129 mnodes[0]*mnodes[1]*mnodes[2], ncomp, 130 CEED_MEM_HOST, CEED_OWN_POINTER, idx, Erestrict); 131 132 PetscFunctionReturn(0); 133 } 134 135 // Data for PETSc 136 typedef struct User_ *User; 137 struct User_ { 138 MPI_Comm comm; 139 VecScatter ltog; // Scatter for all entries 140 VecScatter ltog0; // Skip Dirichlet values 141 VecScatter gtogD; // global-to-global; only Dirichlet values 142 Vec Xloc, Yloc; 143 CeedVector xceed, yceed; 144 CeedOperator op; 145 CeedVector qdata; 146 Ceed ceed; 147 CeedMemType memtype; 148 int (*VecGetArray)(Vec, PetscScalar **); 149 int (*VecGetArrayRead)(Vec, const PetscScalar **); 150 int (*VecRestoreArray)(Vec, PetscScalar **); 151 int (*VecRestoreArrayRead)(Vec, const PetscScalar **); 152 }; 153 154 // MemType Options 155 static const char *const memTypes[] = {"host","device","memType", 156 "CEED_MEM_",0 157 }; 158 159 // BP Options 160 typedef enum { 161 CEED_BP1 = 0, CEED_BP2 = 1, CEED_BP3 = 2, 162 CEED_BP4 = 3, CEED_BP5 = 4, CEED_BP6 = 5 163 } bpType; 164 static const char *const bpTypes[] = {"bp1","bp2","bp3","bp4","bp5","bp6", 165 "bpType","CEED_BP",0 166 }; 167 168 // BP specific data 169 typedef struct { 170 CeedInt ncompu, qdatasize, qextra; 171 CeedQFunctionUser setupgeo, setuprhs, apply, error; 172 const char *setupgeofname, *setuprhsfname, *applyfname, *errorfname; 173 CeedEvalMode inmode, outmode; 174 CeedQuadMode qmode; 175 } bpData; 176 177 bpData bpOptions[6] = { 178 [CEED_BP1] = { 179 .ncompu = 1, 180 .qdatasize = 1, 181 .qextra = 1, 182 .setupgeo = SetupMassGeo, 183 .setuprhs = SetupMassRhs, 184 .apply = Mass, 185 .error = Error, 186 .setupgeofname = SetupMassGeo_loc, 187 .setuprhsfname = SetupMassRhs_loc, 188 .applyfname = Mass_loc, 189 .errorfname = Error_loc, 190 .inmode = CEED_EVAL_INTERP, 191 .outmode = CEED_EVAL_INTERP, 192 .qmode = CEED_GAUSS 193 }, 194 [CEED_BP2] = { 195 .ncompu = 3, 196 .qdatasize = 1, 197 .qextra = 1, 198 .setupgeo = SetupMassGeo, 199 .setuprhs = SetupMassRhs3, 200 .apply = Mass3, 201 .error = Error3, 202 .setupgeofname = SetupMassGeo_loc, 203 .setuprhsfname = SetupMassRhs3_loc, 204 .applyfname = Mass3_loc, 205 .errorfname = Error3_loc, 206 .inmode = CEED_EVAL_INTERP, 207 .outmode = CEED_EVAL_INTERP, 208 .qmode = CEED_GAUSS 209 }, 210 [CEED_BP3] = { 211 .ncompu = 1, 212 .qdatasize = 6, 213 .qextra = 1, 214 .setupgeo = SetupDiffGeo, 215 .setuprhs = SetupDiffRhs, 216 .apply = Diff, 217 .error = Error, 218 .setupgeofname = SetupDiffGeo_loc, 219 .setuprhsfname = SetupDiffRhs_loc, 220 .applyfname = Diff_loc, 221 .errorfname = Error_loc, 222 .inmode = CEED_EVAL_GRAD, 223 .outmode = CEED_EVAL_GRAD, 224 .qmode = CEED_GAUSS 225 }, 226 [CEED_BP4] = { 227 .ncompu = 3, 228 .qdatasize = 6, 229 .qextra = 1, 230 .setupgeo = SetupDiffGeo, 231 .setuprhs = SetupDiffRhs3, 232 .apply = Diff3, 233 .error = Error3, 234 .setupgeofname = SetupDiffGeo_loc, 235 .setuprhsfname = SetupDiffRhs3_loc, 236 .applyfname = Diff_loc, 237 .errorfname = Error3_loc, 238 .inmode = CEED_EVAL_GRAD, 239 .outmode = CEED_EVAL_GRAD, 240 .qmode = CEED_GAUSS 241 }, 242 [CEED_BP5] = { 243 .ncompu = 1, 244 .qdatasize = 6, 245 .qextra = 0, 246 .setupgeo = SetupDiffGeo, 247 .setuprhs = SetupDiffRhs, 248 .apply = Diff, 249 .error = Error, 250 .setupgeofname = SetupDiffGeo_loc, 251 .setuprhsfname = SetupDiffRhs_loc, 252 .applyfname = Diff_loc, 253 .errorfname = Error_loc, 254 .inmode = CEED_EVAL_GRAD, 255 .outmode = CEED_EVAL_GRAD, 256 .qmode = CEED_GAUSS_LOBATTO 257 }, 258 [CEED_BP6] = { 259 .ncompu = 3, 260 .qdatasize = 6, 261 .qextra = 0, 262 .setupgeo = SetupDiffGeo, 263 .setuprhs = SetupDiffRhs3, 264 .apply = Diff3, 265 .error = Error3, 266 .setupgeofname = SetupDiffGeo_loc, 267 .setuprhsfname = SetupDiffRhs3_loc, 268 .applyfname = Diff_loc, 269 .errorfname = Error3_loc, 270 .inmode = CEED_EVAL_GRAD, 271 .outmode = CEED_EVAL_GRAD, 272 .qmode = CEED_GAUSS_LOBATTO 273 } 274 }; 275 276 // This function uses libCEED to compute the action of the mass matrix 277 static PetscErrorCode MatMult_Mass(Mat A, Vec X, Vec Y) { 278 PetscErrorCode ierr; 279 User user; 280 PetscScalar *x, *y; 281 282 PetscFunctionBeginUser; 283 284 ierr = MatShellGetContext(A, &user); CHKERRQ(ierr); 285 286 // Global-to-local 287 ierr = VecScatterBegin(user->ltog, X, user->Xloc, INSERT_VALUES, 288 SCATTER_REVERSE); CHKERRQ(ierr); 289 ierr = VecScatterEnd(user->ltog, X, user->Xloc, INSERT_VALUES, 290 SCATTER_REVERSE); CHKERRQ(ierr); 291 ierr = VecZeroEntries(user->Yloc); CHKERRQ(ierr); 292 293 // Setup libCEED vectors 294 ierr = user->VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); 295 CHKERRQ(ierr); 296 ierr = user->VecGetArray(user->Yloc, &y); CHKERRQ(ierr); 297 CeedVectorSetArray(user->xceed, user->memtype, CEED_USE_POINTER, x); 298 CeedVectorSetArray(user->yceed, user->memtype, CEED_USE_POINTER, y); 299 300 // Apply libCEED operator 301 CeedOperatorApply(user->op, user->xceed, user->yceed, 302 CEED_REQUEST_IMMEDIATE); 303 CeedVectorSyncArray(user->yceed, user->memtype); 304 305 // Restore PETSc vectors 306 ierr = user->VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x); 307 CHKERRQ(ierr); 308 ierr = user->VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr); 309 310 // Local-to-global 311 if (Y) { 312 ierr = VecZeroEntries(Y); CHKERRQ(ierr); 313 ierr = VecScatterBegin(user->ltog, user->Yloc, Y, ADD_VALUES, 314 SCATTER_FORWARD); CHKERRQ(ierr); 315 ierr = VecScatterEnd(user->ltog, user->Yloc, Y, ADD_VALUES, 316 SCATTER_FORWARD); CHKERRQ(ierr); 317 } 318 PetscFunctionReturn(0); 319 } 320 321 // This function uses libCEED to compute the action of the Laplacian with 322 // Dirichlet boundary conditions 323 static PetscErrorCode MatMult_Diff(Mat A, Vec X, Vec Y) { 324 PetscErrorCode ierr; 325 User user; 326 PetscScalar *x, *y; 327 328 PetscFunctionBeginUser; 329 330 ierr = MatShellGetContext(A, &user); CHKERRQ(ierr); 331 332 // Global-to-local 333 ierr = VecScatterBegin(user->ltog0, X, user->Xloc, INSERT_VALUES, 334 SCATTER_REVERSE); CHKERRQ(ierr); 335 ierr = VecScatterEnd(user->ltog0, X, user->Xloc, INSERT_VALUES, 336 SCATTER_REVERSE); 337 CHKERRQ(ierr); 338 ierr = VecZeroEntries(user->Yloc); CHKERRQ(ierr); 339 340 // Setup libCEED vectors 341 ierr = user->VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); 342 CHKERRQ(ierr); 343 ierr = user->VecGetArray(user->Yloc, &y); CHKERRQ(ierr); 344 CeedVectorSetArray(user->xceed, user->memtype, CEED_USE_POINTER, x); 345 CeedVectorSetArray(user->yceed, user->memtype, CEED_USE_POINTER, y); 346 347 // Apply libCEED operator 348 CeedOperatorApply(user->op, user->xceed, user->yceed, 349 CEED_REQUEST_IMMEDIATE); 350 CeedVectorSyncArray(user->yceed, user->memtype); 351 352 // Restore PETSc vectors 353 ierr = user->VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x); 354 CHKERRQ(ierr); 355 ierr = user->VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr); 356 357 // Local-to-global 358 ierr = VecZeroEntries(Y); CHKERRQ(ierr); 359 ierr = VecScatterBegin(user->gtogD, X, Y, INSERT_VALUES, SCATTER_FORWARD); 360 CHKERRQ(ierr); 361 ierr = VecScatterEnd(user->gtogD, X, Y, INSERT_VALUES, SCATTER_FORWARD); 362 CHKERRQ(ierr); 363 ierr = VecScatterBegin(user->ltog0, user->Yloc, Y, ADD_VALUES, 364 SCATTER_FORWARD); CHKERRQ(ierr); 365 ierr = VecScatterEnd(user->ltog0, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD); 366 CHKERRQ(ierr); 367 368 PetscFunctionReturn(0); 369 } 370 371 // This function calculates the error in the final solution 372 static PetscErrorCode ComputeErrorMax(User user, CeedOperator operror, Vec X, 373 CeedVector target, PetscReal *maxerror) { 374 PetscErrorCode ierr; 375 PetscScalar *x; 376 CeedVector collocated_error; 377 CeedInt length; 378 379 PetscFunctionBeginUser; 380 381 CeedVectorGetLength(target, &length); 382 CeedVectorCreate(user->ceed, length, &collocated_error); 383 384 // Global-to-local 385 ierr = VecScatterBegin(user->ltog, X, user->Xloc, INSERT_VALUES, 386 SCATTER_REVERSE); CHKERRQ(ierr); 387 ierr = VecScatterEnd(user->ltog, X, user->Xloc, INSERT_VALUES, 388 SCATTER_REVERSE); CHKERRQ(ierr); 389 390 // Setup libCEED vector 391 ierr = user->VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); 392 CHKERRQ(ierr); 393 CeedVectorSetArray(user->xceed, user->memtype, CEED_USE_POINTER, x); 394 395 // Apply libCEED operator 396 CeedOperatorApply(operror, user->xceed, collocated_error, 397 CEED_REQUEST_IMMEDIATE); 398 399 // Restore PETSc vector 400 ierr = user->VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x); 401 CHKERRQ(ierr); 402 403 // Reduce max error 404 *maxerror = 0; 405 const CeedScalar *e; 406 CeedVectorGetArrayRead(collocated_error, CEED_MEM_HOST, &e); 407 for (CeedInt i=0; i<length; i++) { 408 *maxerror = PetscMax(*maxerror, PetscAbsScalar(e[i])); 409 } 410 CeedVectorRestoreArrayRead(collocated_error, &e); 411 ierr = MPI_Allreduce(MPI_IN_PLACE, maxerror, 1, MPIU_REAL, MPIU_MAX, 412 user->comm); CHKERRQ(ierr); 413 414 // Cleanup 415 CeedVectorDestroy(&collocated_error); 416 417 PetscFunctionReturn(0); 418 } 419 420 int main(int argc, char **argv) { 421 PetscInt ierr; 422 MPI_Comm comm; 423 char ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self"; 424 double my_rt_start, my_rt, rt_min, rt_max; 425 PetscInt degree, qextra, localnodes, localelem, melem[3], mnodes[3], p[3], 426 irank[3], lnodes[3], lsize, ncompu = 1; 427 PetscScalar *r; 428 PetscBool test_mode, benchmark_mode, write_solution; 429 PetscMPIInt size, rank; 430 Vec X, Xloc, rhs, rhsloc; 431 Mat mat; 432 KSP ksp; 433 VecScatter ltog, ltog0, gtogD; 434 User user; 435 Ceed ceed; 436 CeedBasis basisx, basisu; 437 CeedElemRestriction Erestrictx, Erestrictu, Erestrictui, Erestrictqdi; 438 CeedQFunction qfsetupgeo, qfsetuprhs, qfapply, qferror; 439 CeedOperator opsetupgeo, opsetuprhs, opapply, operror; 440 CeedVector xcoord, qdata, rhsceed, target; 441 CeedMemType memtyperequested; 442 CeedInt P, Q; 443 const CeedInt dim = 3, ncompx = 3; 444 bpType bpchoice; 445 446 // Check for PETSc CUDA availability 447 PetscBool petschavecuda, setmemtyperequest = PETSC_FALSE; 448 // *INDENT-OFF* 449 #ifdef PETSC_HAVE_CUDA 450 petschavecuda = PETSC_TRUE; 451 #else 452 petschavecuda = PETSC_FALSE; 453 #endif 454 // *INDENT-ON* 455 456 ierr = PetscInitialize(&argc, &argv, NULL, help); 457 if (ierr) return ierr; 458 comm = PETSC_COMM_WORLD; 459 460 // Read command line options 461 ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr); 462 bpchoice = CEED_BP1; 463 ierr = PetscOptionsEnum("-problem", 464 "CEED benchmark problem to solve", NULL, 465 bpTypes, (PetscEnum)bpchoice, (PetscEnum *)&bpchoice, 466 NULL); CHKERRQ(ierr); 467 ncompu = bpOptions[bpchoice].ncompu; 468 test_mode = PETSC_FALSE; 469 ierr = PetscOptionsBool("-test", 470 "Testing mode (do not print unless error is large)", 471 NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr); 472 benchmark_mode = PETSC_FALSE; 473 ierr = PetscOptionsBool("-benchmark", 474 "Benchmarking mode (prints benchmark statistics)", 475 NULL, benchmark_mode, &benchmark_mode, NULL); 476 CHKERRQ(ierr); 477 write_solution = PETSC_FALSE; 478 ierr = PetscOptionsBool("-write_solution", 479 "Write solution for visualization", 480 NULL, write_solution, &write_solution, NULL); 481 CHKERRQ(ierr); 482 degree = test_mode ? 3 : 1; 483 ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", 484 NULL, degree, °ree, NULL); CHKERRQ(ierr); 485 qextra = bpOptions[bpchoice].qextra; 486 ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", 487 NULL, qextra, &qextra, NULL); CHKERRQ(ierr); 488 ierr = PetscOptionsString("-ceed", "CEED resource specifier", 489 NULL, ceedresource, ceedresource, 490 sizeof(ceedresource), NULL); CHKERRQ(ierr); 491 localnodes = 1000; 492 ierr = PetscOptionsInt("-local", 493 "Target number of locally owned nodes per process", 494 NULL, localnodes, &localnodes, NULL); CHKERRQ(ierr); 495 memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST; 496 ierr = PetscOptionsEnum("-memtype", 497 "CEED MemType requested", NULL, 498 memTypes, (PetscEnum)memtyperequested, 499 (PetscEnum *)&memtyperequested, &setmemtyperequest); 500 CHKERRQ(ierr); 501 ierr = PetscOptionsEnd(); CHKERRQ(ierr); 502 P = degree + 1; 503 Q = P + qextra; 504 505 // Set up libCEED 506 CeedInit(ceedresource, &ceed); 507 CeedMemType memtypebackend; 508 CeedGetPreferredMemType(ceed, &memtypebackend); 509 510 // Check memtype compatibility 511 if (!setmemtyperequest) 512 memtyperequested = memtypebackend; 513 else if (!petschavecuda && memtyperequested == CEED_MEM_DEVICE) 514 SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, 515 "PETSc was not built with CUDA. " 516 "Requested MemType CEED_MEM_DEVICE is not supported.", NULL); 517 518 // Determine size of process grid 519 ierr = MPI_Comm_size(comm, &size); CHKERRQ(ierr); 520 Split3(size, p, false); 521 522 // Find a nicely composite number of elements no less than localnodes 523 for (localelem = PetscMax(1, localnodes / (degree*degree*degree)); ; 524 localelem++) { 525 Split3(localelem, melem, true); 526 if (Max3(melem) / Min3(melem) <= 2) break; 527 } 528 529 // Find my location in the process grid 530 ierr = MPI_Comm_rank(comm, &rank); CHKERRQ(ierr); 531 for (int d=0,rankleft=rank; d<dim; d++) { 532 const int pstride[3] = {p[1] *p[2], p[2], 1}; 533 irank[d] = rankleft / pstride[d]; 534 rankleft -= irank[d] * pstride[d]; 535 } 536 537 GlobalNodes(p, irank, degree, melem, mnodes); 538 539 // Setup global vector 540 if (memtyperequested == CEED_MEM_DEVICE) { 541 ierr = VecSetType(X, VECCUDA); CHKERRQ(ierr); 542 } 543 ierr = VecCreate(comm, &X); CHKERRQ(ierr); 544 ierr = VecSetSizes(X, mnodes[0]*mnodes[1]*mnodes[2]*ncompu, PETSC_DECIDE); 545 CHKERRQ(ierr); 546 ierr = VecSetUp(X); CHKERRQ(ierr); 547 548 // Set up libCEED 549 CeedInit(ceedresource, &ceed); 550 551 // Print summary 552 CeedInt gsize; 553 ierr = VecGetSize(X, &gsize); CHKERRQ(ierr); 554 if (!test_mode) { 555 const char *usedresource; 556 CeedGetResource(ceed, &usedresource); 557 558 VecType vectype; 559 ierr = VecGetType(X, &vectype); CHKERRQ(ierr); 560 561 ierr = PetscPrintf(comm, 562 "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n" 563 " PETSc:\n" 564 " PETSc Vec Type : %s\n" 565 " libCEED:\n" 566 " libCEED Backend : %s\n" 567 " libCEED Backend MemType : %s\n" 568 " libCEED User Requested MemType : %s\n" 569 " Mesh:\n" 570 " Number of 1D Basis Nodes (p) : %d\n" 571 " Number of 1D Quadrature Points (q) : %d\n" 572 " Global nodes : %D\n" 573 " Process Decomposition : %D %D %D\n" 574 " Local Elements : %D = %D %D %D\n" 575 " Owned nodes : %D = %D %D %D\n" 576 " DoF per node : %D\n", 577 bpchoice+1, vectype, usedresource, 578 CeedMemTypes[memtypebackend], 579 (setmemtyperequest) ? 580 CeedMemTypes[memtyperequested] : "none", 581 P, Q, gsize/ncompu, p[0], p[1], p[2], localelem, 582 melem[0], melem[1], melem[2], 583 mnodes[0]*mnodes[1]*mnodes[2], mnodes[0], mnodes[1], 584 mnodes[2], ncompu); CHKERRQ(ierr); 585 } 586 587 { 588 lsize = 1; 589 for (int d=0; d<dim; d++) { 590 lnodes[d] = melem[d]*degree + 1; 591 lsize *= lnodes[d]; 592 } 593 if (memtyperequested == CEED_MEM_DEVICE) { 594 ierr = VecSetType(Xloc, VECCUDA); CHKERRQ(ierr); 595 } 596 ierr = VecCreate(PETSC_COMM_SELF, &Xloc); CHKERRQ(ierr); 597 ierr = VecSetSizes(Xloc, lsize*ncompu, PETSC_DECIDE); CHKERRQ(ierr); 598 ierr = VecSetUp(Xloc); CHKERRQ(ierr); 599 600 // Create local-to-global scatter 601 PetscInt *ltogind, *ltogind0, *locind, l0count; 602 IS ltogis, ltogis0, locis; 603 PetscInt gstart[2][2][2], gmnodes[2][2][2][dim]; 604 605 for (int i=0; i<2; i++) { 606 for (int j=0; j<2; j++) { 607 for (int k=0; k<2; k++) { 608 PetscInt ijkrank[3] = {irank[0]+i, irank[1]+j, irank[2]+k}; 609 gstart[i][j][k] = GlobalStart(p, ijkrank, degree, melem); 610 GlobalNodes(p, ijkrank, degree, melem, gmnodes[i][j][k]); 611 } 612 } 613 } 614 615 ierr = PetscMalloc1(lsize, <ogind); CHKERRQ(ierr); 616 ierr = PetscMalloc1(lsize, <ogind0); CHKERRQ(ierr); 617 ierr = PetscMalloc1(lsize, &locind); CHKERRQ(ierr); 618 l0count = 0; 619 for (PetscInt i=0,ir,ii; ir=i>=mnodes[0], ii=i-ir*mnodes[0], i<lnodes[0]; i++) 620 for (PetscInt j=0,jr,jj; jr=j>=mnodes[1], jj=j-jr*mnodes[1], j<lnodes[1]; j++) 621 for (PetscInt k=0,kr,kk; kr=k>=mnodes[2], kk=k-kr*mnodes[2], k<lnodes[2]; k++) { 622 PetscInt here = (i*lnodes[1]+j)*lnodes[2]+k; 623 ltogind[here] = 624 gstart[ir][jr][kr] + (ii*gmnodes[ir][jr][kr][1]+jj)*gmnodes[ir][jr][kr][2]+kk; 625 if ((irank[0] == 0 && i == 0) 626 || (irank[1] == 0 && j == 0) 627 || (irank[2] == 0 && k == 0) 628 || (irank[0]+1 == p[0] && i+1 == lnodes[0]) 629 || (irank[1]+1 == p[1] && j+1 == lnodes[1]) 630 || (irank[2]+1 == p[2] && k+1 == lnodes[2])) 631 continue; 632 ltogind0[l0count] = ltogind[here]; 633 locind[l0count++] = here; 634 } 635 ierr = ISCreateBlock(comm, ncompu, lsize, ltogind, PETSC_OWN_POINTER, 636 <ogis); CHKERRQ(ierr); 637 ierr = VecScatterCreate(Xloc, NULL, X, ltogis, <og); CHKERRQ(ierr); 638 CHKERRQ(ierr); 639 ierr = ISCreateBlock(comm, ncompu, l0count, ltogind0, PETSC_OWN_POINTER, 640 <ogis0); CHKERRQ(ierr); 641 ierr = ISCreateBlock(comm, ncompu, l0count, locind, PETSC_OWN_POINTER, 642 &locis); CHKERRQ(ierr); 643 ierr = VecScatterCreate(Xloc, locis, X, ltogis0, <og0); CHKERRQ(ierr); 644 { 645 // Create global-to-global scatter for Dirichlet values (everything not in 646 // ltogis0, which is the range of ltog0) 647 PetscInt xstart, xend, *indD, countD = 0; 648 IS isD; 649 const PetscScalar *x; 650 ierr = VecZeroEntries(Xloc); CHKERRQ(ierr); 651 ierr = VecSet(X, 1.0); CHKERRQ(ierr); 652 ierr = VecScatterBegin(ltog0, Xloc, X, INSERT_VALUES, SCATTER_FORWARD); 653 CHKERRQ(ierr); 654 ierr = VecScatterEnd(ltog0, Xloc, X, INSERT_VALUES, SCATTER_FORWARD); 655 CHKERRQ(ierr); 656 ierr = VecGetOwnershipRange(X, &xstart, &xend); CHKERRQ(ierr); 657 ierr = PetscMalloc1(xend-xstart, &indD); CHKERRQ(ierr); 658 ierr = VecGetArrayRead(X, &x); CHKERRQ(ierr); 659 for (PetscInt i=0; i<xend-xstart; i++) { 660 if (x[i] == 1.) indD[countD++] = xstart + i; 661 } 662 ierr = VecRestoreArrayRead(X, &x); CHKERRQ(ierr); 663 ierr = ISCreateGeneral(comm, countD, indD, PETSC_COPY_VALUES, &isD); 664 CHKERRQ(ierr); 665 ierr = PetscFree(indD); CHKERRQ(ierr); 666 ierr = VecScatterCreate(X, isD, X, isD, >ogD); CHKERRQ(ierr); 667 ierr = ISDestroy(&isD); CHKERRQ(ierr); 668 } 669 ierr = ISDestroy(<ogis); CHKERRQ(ierr); 670 ierr = ISDestroy(<ogis0); CHKERRQ(ierr); 671 ierr = ISDestroy(&locis); CHKERRQ(ierr); 672 } 673 674 // CEED bases 675 CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompu, P, Q, 676 bpOptions[bpchoice].qmode, &basisu); 677 CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, Q, 678 bpOptions[bpchoice].qmode, &basisx); 679 680 // CEED restrictions 681 CreateRestriction(ceed, CEED_INTERLACED, melem, P, ncompu, &Erestrictu); 682 CreateRestriction(ceed, CEED_NONINTERLACED, melem, 2, dim, &Erestrictx); 683 CeedInt nelem = melem[0]*melem[1]*melem[2]; 684 CeedElemRestrictionCreateStrided(ceed, nelem, Q*Q*Q, nelem*Q*Q*Q, ncompu, 685 CEED_STRIDES_BACKEND, &Erestrictui); 686 CeedElemRestrictionCreateStrided(ceed, nelem, Q*Q*Q, nelem*Q*Q*Q, 687 bpOptions[bpchoice].qdatasize, 688 CEED_STRIDES_BACKEND, &Erestrictqdi); 689 { 690 CeedScalar *xloc; 691 CeedInt shape[3] = {melem[0]+1, melem[1]+1, melem[2]+1}, len = 692 shape[0]*shape[1]*shape[2]; 693 xloc = malloc(len*ncompx*sizeof xloc[0]); 694 for (CeedInt i=0; i<shape[0]; i++) { 695 for (CeedInt j=0; j<shape[1]; j++) { 696 for (CeedInt k=0; k<shape[2]; k++) { 697 xloc[((i*shape[1]+j)*shape[2]+k) + 0*len] = 1.*(irank[0]*melem[0]+i) / 698 (p[0]*melem[0]); 699 xloc[((i*shape[1]+j)*shape[2]+k) + 1*len] = 1.*(irank[1]*melem[1]+j) / 700 (p[1]*melem[1]); 701 xloc[((i*shape[1]+j)*shape[2]+k) + 2*len] = 1.*(irank[2]*melem[2]+k) / 702 (p[2]*melem[2]); 703 } 704 } 705 } 706 CeedVectorCreate(ceed, len*ncompx, &xcoord); 707 CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_OWN_POINTER, xloc); 708 } 709 710 // Create the Qfunction that builds the operator quadrature data 711 CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].setupgeo, 712 bpOptions[bpchoice].setupgeofname, &qfsetupgeo); 713 CeedQFunctionAddInput(qfsetupgeo, "dx", ncompx*dim, CEED_EVAL_GRAD); 714 CeedQFunctionAddInput(qfsetupgeo, "weight", 1, CEED_EVAL_WEIGHT); 715 CeedQFunctionAddOutput(qfsetupgeo, "qdata", bpOptions[bpchoice].qdatasize, 716 CEED_EVAL_NONE); 717 718 // Create the Qfunction that sets up the RHS and true solution 719 CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].setuprhs, 720 bpOptions[bpchoice].setuprhsfname, &qfsetuprhs); 721 CeedQFunctionAddInput(qfsetuprhs, "x", ncompx, CEED_EVAL_INTERP); 722 CeedQFunctionAddInput(qfsetuprhs, "dx", ncompx*dim, CEED_EVAL_GRAD); 723 CeedQFunctionAddInput(qfsetuprhs, "weight", 1, CEED_EVAL_WEIGHT); 724 CeedQFunctionAddOutput(qfsetuprhs, "true_soln", ncompu, CEED_EVAL_NONE); 725 CeedQFunctionAddOutput(qfsetuprhs, "rhs", ncompu, CEED_EVAL_INTERP); 726 727 // Set up PDE operator 728 CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].apply, 729 bpOptions[bpchoice].applyfname, &qfapply); 730 // Add inputs and outputs 731 CeedInt inscale = bpOptions[bpchoice].inmode==CEED_EVAL_GRAD ? 3 : 1; 732 CeedInt outscale = bpOptions[bpchoice].outmode==CEED_EVAL_GRAD ? 3 : 1; 733 CeedQFunctionAddInput(qfapply, "u", ncompu*inscale, 734 bpOptions[bpchoice].inmode); 735 CeedQFunctionAddInput(qfapply, "qdata", bpOptions[bpchoice].qdatasize, 736 CEED_EVAL_NONE); 737 CeedQFunctionAddOutput(qfapply, "v", ncompu*outscale, 738 bpOptions[bpchoice].outmode); 739 740 // Create the error qfunction 741 CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].error, 742 bpOptions[bpchoice].errorfname, &qferror); 743 CeedQFunctionAddInput(qferror, "u", ncompu, CEED_EVAL_INTERP); 744 CeedQFunctionAddInput(qferror, "true_soln", ncompu, CEED_EVAL_NONE); 745 CeedQFunctionAddOutput(qferror, "error", ncompu, CEED_EVAL_NONE); 746 747 // Create the persistent vectors that will be needed in setup 748 CeedInt nqpts; 749 CeedBasisGetNumQuadraturePoints(basisu, &nqpts); 750 CeedVectorCreate(ceed, bpOptions[bpchoice].qdatasize*nelem*nqpts, &qdata); 751 CeedVectorCreate(ceed, nelem*nqpts*ncompu, &target); 752 CeedVectorCreate(ceed, lsize*ncompu, &rhsceed); 753 754 // Create the operator that builds the quadrature data for the ceed operator 755 CeedOperatorCreate(ceed, qfsetupgeo, CEED_QFUNCTION_NONE, 756 CEED_QFUNCTION_NONE, &opsetupgeo); 757 CeedOperatorSetField(opsetupgeo, "dx", Erestrictx, basisx, 758 CEED_VECTOR_ACTIVE); 759 CeedOperatorSetField(opsetupgeo, "weight", CEED_ELEMRESTRICTION_NONE, basisx, 760 CEED_VECTOR_NONE); 761 CeedOperatorSetField(opsetupgeo, "qdata", Erestrictqdi, 762 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 763 764 // Create the operator that builds the RHS and true solution 765 CeedOperatorCreate(ceed, qfsetuprhs, CEED_QFUNCTION_NONE, 766 CEED_QFUNCTION_NONE, &opsetuprhs); 767 CeedOperatorSetField(opsetuprhs, "x", Erestrictx, basisx, 768 CEED_VECTOR_ACTIVE); 769 CeedOperatorSetField(opsetuprhs, "dx", Erestrictx, basisx, 770 CEED_VECTOR_ACTIVE); 771 CeedOperatorSetField(opsetuprhs, "weight", CEED_ELEMRESTRICTION_NONE, basisx, 772 CEED_VECTOR_NONE); 773 CeedOperatorSetField(opsetuprhs, "true_soln", Erestrictui, 774 CEED_BASIS_COLLOCATED, target); 775 CeedOperatorSetField(opsetuprhs, "rhs", Erestrictu, basisu, 776 CEED_VECTOR_ACTIVE); 777 778 // Create the mass or diff operator 779 CeedOperatorCreate(ceed, qfapply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, 780 &opapply); 781 CeedOperatorSetField(opapply, "u", Erestrictu, basisu, CEED_VECTOR_ACTIVE); 782 CeedOperatorSetField(opapply, "qdata", Erestrictqdi, CEED_BASIS_COLLOCATED, 783 qdata); 784 CeedOperatorSetField(opapply, "v", Erestrictu, basisu, CEED_VECTOR_ACTIVE); 785 786 // Create the error operator 787 CeedOperatorCreate(ceed, qferror, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, 788 &operror); 789 CeedOperatorSetField(operror, "u", Erestrictu, basisu, CEED_VECTOR_ACTIVE); 790 CeedOperatorSetField(operror, "true_soln", Erestrictui, 791 CEED_BASIS_COLLOCATED, target); 792 CeedOperatorSetField(operror, "error", Erestrictui, CEED_BASIS_COLLOCATED, 793 CEED_VECTOR_ACTIVE); 794 795 // Set up Mat 796 ierr = PetscMalloc1(1, &user); CHKERRQ(ierr); 797 user->comm = comm; 798 user->ltog = ltog; 799 if (bpchoice != CEED_BP1 && bpchoice != CEED_BP2) { 800 user->ltog0 = ltog0; 801 user->gtogD = gtogD; 802 } 803 user->Xloc = Xloc; 804 ierr = VecDuplicate(Xloc, &user->Yloc); CHKERRQ(ierr); 805 CeedVectorCreate(ceed, lsize*ncompu, &user->xceed); 806 CeedVectorCreate(ceed, lsize*ncompu, &user->yceed); 807 user->op = opapply; 808 user->qdata = qdata; 809 user->ceed = ceed; 810 user->memtype = memtyperequested; 811 if (memtyperequested == CEED_MEM_HOST) { 812 user->VecGetArray = VecGetArray; 813 user->VecGetArrayRead = VecGetArrayRead; 814 user->VecRestoreArray = VecRestoreArray; 815 user->VecRestoreArrayRead = VecRestoreArrayRead; 816 } else { 817 user->VecGetArray = VecCUDAGetArray; 818 user->VecGetArrayRead = VecCUDAGetArrayRead; 819 user->VecRestoreArray = VecCUDARestoreArray; 820 user->VecRestoreArrayRead = VecCUDARestoreArrayRead; 821 } 822 823 ierr = MatCreateShell(comm, mnodes[0]*mnodes[1]*mnodes[2]*ncompu, 824 mnodes[0]*mnodes[1]*mnodes[2]*ncompu, 825 PETSC_DECIDE, PETSC_DECIDE, user, &mat); CHKERRQ(ierr); 826 if (bpchoice == CEED_BP1 || bpchoice == CEED_BP2) { 827 ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Mass); 828 CHKERRQ(ierr); 829 } else { 830 ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Diff); 831 CHKERRQ(ierr); 832 } 833 if (user->memtype == CEED_MEM_DEVICE) { 834 ierr = MatShellSetVecType(mat, VECCUDA); CHKERRQ(ierr); 835 } 836 837 // Get RHS vector 838 ierr = VecDuplicate(X, &rhs); CHKERRQ(ierr); 839 ierr = VecDuplicate(Xloc, &rhsloc); CHKERRQ(ierr); 840 ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr); 841 ierr = user->VecGetArray(rhsloc, &r); CHKERRQ(ierr); 842 CeedVectorSetArray(rhsceed, user->memtype, CEED_USE_POINTER, r); 843 844 // Setup qdata, rhs, and target 845 CeedOperatorApply(opsetupgeo, xcoord, qdata, CEED_REQUEST_IMMEDIATE); 846 CeedOperatorApply(opsetuprhs, xcoord, rhsceed, CEED_REQUEST_IMMEDIATE); 847 ierr = CeedVectorSyncArray(rhsceed, user->memtype); CHKERRQ(ierr); 848 CeedVectorDestroy(&xcoord); 849 850 // Gather RHS 851 ierr = user->VecRestoreArray(rhsloc, &r); CHKERRQ(ierr); 852 ierr = VecZeroEntries(rhs); CHKERRQ(ierr); 853 ierr = VecScatterBegin(ltog, rhsloc, rhs, ADD_VALUES, SCATTER_FORWARD); 854 CHKERRQ(ierr); 855 ierr = VecScatterEnd(ltog, rhsloc, rhs, ADD_VALUES, SCATTER_FORWARD); 856 CHKERRQ(ierr); 857 CeedVectorDestroy(&rhsceed); 858 859 ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr); 860 { 861 PC pc; 862 ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 863 if (bpchoice == CEED_BP1 || bpchoice == CEED_BP2) { 864 ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr); 865 ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr); 866 } else { 867 ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr); 868 } 869 ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 870 ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 871 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 872 PETSC_DEFAULT); CHKERRQ(ierr); 873 } 874 ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr); 875 ierr = KSPSetOperators(ksp, mat, mat); CHKERRQ(ierr); 876 // First run, if benchmarking 877 if (benchmark_mode) { 878 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1); 879 CHKERRQ(ierr); 880 my_rt_start = MPI_Wtime(); 881 ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 882 my_rt = MPI_Wtime() - my_rt_start; 883 ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm); 884 CHKERRQ(ierr); 885 // Set maxits based on first iteration timing 886 if (my_rt > 0.02) { 887 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5); 888 CHKERRQ(ierr); 889 } else { 890 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20); 891 CHKERRQ(ierr); 892 } 893 } 894 // Timed solve 895 ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr); 896 my_rt_start = MPI_Wtime(); 897 ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 898 my_rt = MPI_Wtime() - my_rt_start; 899 { 900 KSPType ksptype; 901 KSPConvergedReason reason; 902 PetscReal rnorm; 903 PetscInt its; 904 ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr); 905 ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr); 906 ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr); 907 ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr); 908 if (!test_mode || reason < 0 || rnorm > 1e-8) { 909 ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm); 910 CHKERRQ(ierr); 911 ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm); 912 CHKERRQ(ierr); 913 ierr = PetscPrintf(comm, 914 " KSP:\n" 915 " KSP Type : %s\n" 916 " KSP Convergence : %s\n" 917 " Total KSP Iterations : %D\n" 918 " Final rnorm : %e\n", 919 ksptype, KSPConvergedReasons[reason], its, 920 (double)rnorm); CHKERRQ(ierr); 921 ierr = PetscPrintf(comm, 922 " Performance:\n" 923 " CG Solve Time : %g (%g) sec\n" 924 " DoFs/Sec in CG : %g (%g) million\n", 925 rt_max, rt_min, 1e-6*gsize*its/rt_max, 926 1e-6*gsize*its/rt_min); CHKERRQ(ierr); 927 } 928 if (benchmark_mode && (!test_mode)) { 929 ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm); 930 CHKERRQ(ierr); 931 ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm); 932 CHKERRQ(ierr); 933 ierr = PetscPrintf(comm, 934 " Performance:\n" 935 " CG Solve Time : %g (%g) sec\n" 936 " DoFs/Sec in CG : %g (%g) million\n", 937 rt_max, rt_min, 1e-6*gsize*its/rt_max, 938 1e-6*gsize*its/rt_min); CHKERRQ(ierr); 939 } 940 } 941 942 { 943 PetscReal maxerror; 944 ierr = ComputeErrorMax(user, operror, X, target, &maxerror); CHKERRQ(ierr); 945 PetscReal tol = 5e-2; 946 if (!test_mode || maxerror > tol) { 947 ierr = PetscPrintf(comm, 948 " Pointwise Error (max) : %e\n", 949 (double)maxerror); CHKERRQ(ierr); 950 } 951 } 952 953 if (write_solution) { 954 PetscViewer vtkviewersoln; 955 956 ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr); 957 ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr); 958 ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr); 959 ierr = VecView(X, vtkviewersoln); CHKERRQ(ierr); 960 ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr); 961 } 962 963 ierr = VecDestroy(&rhs); CHKERRQ(ierr); 964 ierr = VecDestroy(&rhsloc); CHKERRQ(ierr); 965 ierr = VecDestroy(&X); CHKERRQ(ierr); 966 ierr = VecDestroy(&user->Xloc); CHKERRQ(ierr); 967 ierr = VecDestroy(&user->Yloc); CHKERRQ(ierr); 968 ierr = VecScatterDestroy(<og); CHKERRQ(ierr); 969 ierr = VecScatterDestroy(<og0); CHKERRQ(ierr); 970 ierr = VecScatterDestroy(>ogD); CHKERRQ(ierr); 971 ierr = MatDestroy(&mat); CHKERRQ(ierr); 972 ierr = KSPDestroy(&ksp); CHKERRQ(ierr); 973 974 CeedVectorDestroy(&user->xceed); 975 CeedVectorDestroy(&user->yceed); 976 CeedVectorDestroy(&user->qdata); 977 CeedVectorDestroy(&target); 978 CeedOperatorDestroy(&opsetupgeo); 979 CeedOperatorDestroy(&opsetuprhs); 980 CeedOperatorDestroy(&opapply); 981 CeedOperatorDestroy(&operror); 982 CeedElemRestrictionDestroy(&Erestrictu); 983 CeedElemRestrictionDestroy(&Erestrictx); 984 CeedElemRestrictionDestroy(&Erestrictui); 985 CeedElemRestrictionDestroy(&Erestrictqdi); 986 CeedQFunctionDestroy(&qfsetupgeo); 987 CeedQFunctionDestroy(&qfsetuprhs); 988 CeedQFunctionDestroy(&qfapply); 989 CeedQFunctionDestroy(&qferror); 990 CeedBasisDestroy(&basisu); 991 CeedBasisDestroy(&basisx); 992 CeedDestroy(&ceed); 993 ierr = PetscFree(user); CHKERRQ(ierr); 994 return PetscFinalize(); 995 } 996