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