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 // This function uses libCEED to compute the action of the mass matrix 262 static PetscErrorCode MatMult_Mass(Mat A, Vec X, Vec Y) { 263 PetscErrorCode ierr; 264 User user; 265 PetscScalar *x, *y; 266 267 PetscFunctionBeginUser; 268 ierr = MatShellGetContext(A, &user); CHKERRQ(ierr); 269 ierr = VecScatterBegin(user->ltog, X, user->Xloc, INSERT_VALUES, 270 SCATTER_REVERSE); CHKERRQ(ierr); 271 ierr = VecScatterEnd(user->ltog, X, user->Xloc, INSERT_VALUES, SCATTER_REVERSE); 272 CHKERRQ(ierr); 273 ierr = VecZeroEntries(user->Yloc); CHKERRQ(ierr); 274 275 ierr = VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr); 276 ierr = VecGetArray(user->Yloc, &y); CHKERRQ(ierr); 277 CeedVectorSetArray(user->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x); 278 CeedVectorSetArray(user->yceed, CEED_MEM_HOST, CEED_USE_POINTER, y); 279 280 CeedOperatorApply(user->op, user->xceed, user->yceed, 281 CEED_REQUEST_IMMEDIATE); 282 ierr = CeedVectorSyncArray(user->yceed, CEED_MEM_HOST); CHKERRQ(ierr); 283 284 ierr = VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr); 285 ierr = VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr); 286 287 if (Y) { 288 ierr = VecZeroEntries(Y); CHKERRQ(ierr); 289 ierr = VecScatterBegin(user->ltog, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD); 290 CHKERRQ(ierr); 291 ierr = VecScatterEnd(user->ltog, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD); 292 CHKERRQ(ierr); 293 } 294 PetscFunctionReturn(0); 295 } 296 297 // This function uses libCEED to compute the action of the Laplacian with 298 // Dirichlet boundary conditions 299 static PetscErrorCode MatMult_Diff(Mat A, Vec X, Vec Y) { 300 PetscErrorCode ierr; 301 User user; 302 PetscScalar *x, *y; 303 304 PetscFunctionBeginUser; 305 ierr = MatShellGetContext(A, &user); CHKERRQ(ierr); 306 307 // Global-to-local 308 ierr = VecScatterBegin(user->ltog0, X, user->Xloc, INSERT_VALUES, 309 SCATTER_REVERSE); CHKERRQ(ierr); 310 ierr = VecScatterEnd(user->ltog0, X, user->Xloc, INSERT_VALUES, 311 SCATTER_REVERSE); 312 CHKERRQ(ierr); 313 ierr = VecZeroEntries(user->Yloc); CHKERRQ(ierr); 314 315 // Setup CEED vectors 316 ierr = VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr); 317 ierr = VecGetArray(user->Yloc, &y); CHKERRQ(ierr); 318 CeedVectorSetArray(user->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x); 319 CeedVectorSetArray(user->yceed, CEED_MEM_HOST, CEED_USE_POINTER, y); 320 321 // Apply CEED operator 322 CeedOperatorApply(user->op, user->xceed, user->yceed, 323 CEED_REQUEST_IMMEDIATE); 324 ierr = CeedVectorSyncArray(user->yceed, CEED_MEM_HOST); CHKERRQ(ierr); 325 326 // Restore PETSc vectors 327 ierr = VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr); 328 ierr = VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr); 329 330 // Local-to-global 331 ierr = VecZeroEntries(Y); CHKERRQ(ierr); 332 ierr = VecScatterBegin(user->gtogD, X, Y, INSERT_VALUES, SCATTER_FORWARD); 333 CHKERRQ(ierr); 334 ierr = VecScatterEnd(user->gtogD, X, Y, INSERT_VALUES, SCATTER_FORWARD); 335 CHKERRQ(ierr); 336 ierr = VecScatterBegin(user->ltog0, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD); 337 CHKERRQ(ierr); 338 ierr = VecScatterEnd(user->ltog0, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD); 339 CHKERRQ(ierr); 340 341 PetscFunctionReturn(0); 342 } 343 344 // This function calculates the error in the final solution 345 static PetscErrorCode ComputeErrorMax(User user, CeedOperator op_error, Vec X, 346 CeedVector target, PetscReal *maxerror) { 347 PetscErrorCode ierr; 348 PetscScalar *x; 349 CeedVector collocated_error; 350 CeedInt length; 351 352 PetscFunctionBeginUser; 353 CeedVectorGetLength(target, &length); 354 CeedVectorCreate(user->ceed, length, &collocated_error); 355 356 // Global-to-local 357 ierr = VecScatterBegin(user->ltog, X, user->Xloc, INSERT_VALUES, 358 SCATTER_REVERSE); CHKERRQ(ierr); 359 ierr = VecScatterEnd(user->ltog, X, user->Xloc, INSERT_VALUES, SCATTER_REVERSE); 360 CHKERRQ(ierr); 361 362 // Setup CEED vector 363 ierr = VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr); 364 CeedVectorSetArray(user->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x); 365 366 // Apply CEED operator 367 CeedOperatorApply(op_error, user->xceed, collocated_error, 368 CEED_REQUEST_IMMEDIATE); 369 370 // Restore PETSc vector 371 VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr); 372 373 // Reduce max error 374 *maxerror = 0; 375 const CeedScalar *e; 376 CeedVectorGetArrayRead(collocated_error, CEED_MEM_HOST, &e); 377 for (CeedInt i=0; i<length; i++) { 378 *maxerror = PetscMax(*maxerror, PetscAbsScalar(e[i])); 379 } 380 CeedVectorRestoreArrayRead(collocated_error, &e); 381 ierr = MPI_Allreduce(MPI_IN_PLACE, maxerror, 382 1, MPIU_REAL, MPIU_MAX, user->comm); CHKERRQ(ierr); 383 384 // Cleanup 385 CeedVectorDestroy(&collocated_error); 386 387 PetscFunctionReturn(0); 388 } 389 390 int main(int argc, char **argv) { 391 PetscInt ierr; 392 MPI_Comm comm; 393 char ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self"; 394 double my_rt_start, my_rt, rt_min, rt_max; 395 PetscInt degree, qextra, localnodes, localelem, melem[3], mnodes[3], p[3], 396 irank[3], lnodes[3], lsize, ncompu = 1; 397 PetscScalar *r; 398 PetscBool test_mode, benchmark_mode, write_solution; 399 PetscMPIInt size, rank; 400 Vec X, Xloc, rhs, rhsloc; 401 Mat mat; 402 KSP ksp; 403 VecScatter ltog, ltog0, gtogD; 404 User user; 405 Ceed ceed; 406 CeedBasis basisx, basisu; 407 CeedElemRestriction Erestrictx, Erestrictu, Erestrictxi, Erestrictui, 408 Erestrictqdi; 409 CeedQFunction qf_setupgeo, qf_setuprhs, qf_apply, qf_error; 410 CeedOperator op_setupgeo, op_setuprhs, op_apply, op_error; 411 CeedVector xcoord, qdata, rhsceed, target; 412 CeedInt P, Q; 413 const CeedInt dim = 3, ncompx = 3; 414 bpType bpChoice; 415 416 ierr = PetscInitialize(&argc, &argv, NULL, help); 417 if (ierr) return ierr; 418 comm = PETSC_COMM_WORLD; 419 ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr); 420 bpChoice = CEED_BP1; 421 ierr = PetscOptionsEnum("-problem", 422 "CEED benchmark problem to solve", NULL, 423 bpTypes, (PetscEnum)bpChoice, (PetscEnum *)&bpChoice, 424 NULL); CHKERRQ(ierr); 425 ncompu = bpOptions[bpChoice].ncompu; 426 test_mode = PETSC_FALSE; 427 ierr = PetscOptionsBool("-test", 428 "Testing mode (do not print unless error is large)", 429 NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr); 430 benchmark_mode = PETSC_FALSE; 431 ierr = PetscOptionsBool("-benchmark", 432 "Benchmarking mode (prints benchmark statistics)", 433 NULL, benchmark_mode, &benchmark_mode, NULL); 434 CHKERRQ(ierr); 435 write_solution = PETSC_FALSE; 436 ierr = PetscOptionsBool("-write_solution", 437 "Write solution for visualization", 438 NULL, write_solution, &write_solution, NULL); 439 CHKERRQ(ierr); 440 degree = test_mode ? 3 : 1; 441 ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", 442 NULL, degree, °ree, NULL); CHKERRQ(ierr); 443 qextra = bpOptions[bpChoice].qextra; 444 ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", 445 NULL, qextra, &qextra, NULL); CHKERRQ(ierr); 446 ierr = PetscOptionsString("-ceed", "CEED resource specifier", 447 NULL, ceedresource, ceedresource, 448 sizeof(ceedresource), NULL); CHKERRQ(ierr); 449 localnodes = 1000; 450 ierr = PetscOptionsInt("-local", 451 "Target number of locally owned nodes per process", 452 NULL, localnodes, &localnodes, NULL); CHKERRQ(ierr); 453 ierr = PetscOptionsEnd(); CHKERRQ(ierr); 454 P = degree + 1; 455 Q = P + qextra; 456 457 // Determine size of process grid 458 ierr = MPI_Comm_size(comm, &size); CHKERRQ(ierr); 459 Split3(size, p, false); 460 461 // Find a nicely composite number of elements no less than localnodes 462 for (localelem = PetscMax(1, localnodes / (degree*degree*degree)); ; 463 localelem++) { 464 Split3(localelem, melem, true); 465 if (Max3(melem) / Min3(melem) <= 2) break; 466 } 467 468 // Find my location in the process grid 469 ierr = MPI_Comm_rank(comm, &rank); CHKERRQ(ierr); 470 for (int d=0,rankleft=rank; d<dim; d++) { 471 const int pstride[3] = {p[1] *p[2], p[2], 1}; 472 irank[d] = rankleft / pstride[d]; 473 rankleft -= irank[d] * pstride[d]; 474 } 475 476 GlobalNodes(p, irank, degree, melem, mnodes); 477 478 // Setup global vector 479 ierr = VecCreate(comm, &X); CHKERRQ(ierr); 480 ierr = VecSetSizes(X, mnodes[0]*mnodes[1]*mnodes[2]*ncompu, PETSC_DECIDE); 481 CHKERRQ(ierr); 482 ierr = VecSetUp(X); CHKERRQ(ierr); 483 484 // Print summary 485 if (!test_mode) { 486 CeedInt gsize; 487 ierr = VecGetSize(X, &gsize); CHKERRQ(ierr); 488 ierr = PetscPrintf(comm, 489 "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n" 490 " libCEED:\n" 491 " libCEED Backend : %s\n" 492 " Mesh:\n" 493 " Number of 1D Basis Nodes (p) : %d\n" 494 " Number of 1D Quadrature Points (q) : %d\n" 495 " Global nodes : %D\n" 496 " Process Decomposition : %D %D %D\n" 497 " Local Elements : %D = %D %D %D\n" 498 " Owned nodes : %D = %D %D %D\n", 499 bpChoice+1, ceedresource, P, Q, gsize/ncompu, p[0], 500 p[1], p[2], localelem, melem[0], melem[1], melem[2], 501 mnodes[0]*mnodes[1]*mnodes[2], mnodes[0], mnodes[1], mnodes[2]); 502 CHKERRQ(ierr); 503 } 504 505 { 506 lsize = 1; 507 for (int d=0; d<dim; d++) { 508 lnodes[d] = melem[d]*degree + 1; 509 lsize *= lnodes[d]; 510 } 511 ierr = VecCreate(PETSC_COMM_SELF, &Xloc); CHKERRQ(ierr); 512 ierr = VecSetSizes(Xloc, lsize*ncompu, PETSC_DECIDE); CHKERRQ(ierr); 513 ierr = VecSetUp(Xloc); CHKERRQ(ierr); 514 515 // Create local-to-global scatter 516 PetscInt *ltogind, *ltogind0, *locind, l0count; 517 IS ltogis, ltogis0, locis; 518 PetscInt gstart[2][2][2], gmnodes[2][2][2][dim]; 519 520 for (int i=0; i<2; i++) { 521 for (int j=0; j<2; j++) { 522 for (int k=0; k<2; k++) { 523 PetscInt ijkrank[3] = {irank[0]+i, irank[1]+j, irank[2]+k}; 524 gstart[i][j][k] = GlobalStart(p, ijkrank, degree, melem); 525 GlobalNodes(p, ijkrank, degree, melem, gmnodes[i][j][k]); 526 } 527 } 528 } 529 530 ierr = PetscMalloc1(lsize, <ogind); CHKERRQ(ierr); 531 ierr = PetscMalloc1(lsize, <ogind0); CHKERRQ(ierr); 532 ierr = PetscMalloc1(lsize, &locind); CHKERRQ(ierr); 533 l0count = 0; 534 for (PetscInt i=0,ir,ii; ir=i>=mnodes[0], ii=i-ir*mnodes[0], i<lnodes[0]; i++) 535 for (PetscInt j=0,jr,jj; jr=j>=mnodes[1], jj=j-jr*mnodes[1], j<lnodes[1]; j++) 536 for (PetscInt k=0,kr,kk; kr=k>=mnodes[2], kk=k-kr*mnodes[2], k<lnodes[2]; k++) { 537 PetscInt here = (i*lnodes[1]+j)*lnodes[2]+k; 538 ltogind[here] = 539 gstart[ir][jr][kr] + (ii*gmnodes[ir][jr][kr][1]+jj)*gmnodes[ir][jr][kr][2]+kk; 540 if ((irank[0] == 0 && i == 0) 541 || (irank[1] == 0 && j == 0) 542 || (irank[2] == 0 && k == 0) 543 || (irank[0]+1 == p[0] && i+1 == lnodes[0]) 544 || (irank[1]+1 == p[1] && j+1 == lnodes[1]) 545 || (irank[2]+1 == p[2] && k+1 == lnodes[2])) 546 continue; 547 ltogind0[l0count] = ltogind[here]; 548 locind[l0count++] = here; 549 } 550 ierr = ISCreateBlock(comm, ncompu, lsize, ltogind, PETSC_OWN_POINTER, 551 <ogis); CHKERRQ(ierr); 552 ierr = VecScatterCreate(Xloc, NULL, X, ltogis, <og); CHKERRQ(ierr); 553 CHKERRQ(ierr); 554 ierr = ISCreateBlock(comm, ncompu, l0count, ltogind0, PETSC_OWN_POINTER, 555 <ogis0); CHKERRQ(ierr); 556 ierr = ISCreateBlock(comm, ncompu, l0count, locind, PETSC_OWN_POINTER, 557 &locis); CHKERRQ(ierr); 558 ierr = VecScatterCreate(Xloc, locis, X, ltogis0, <og0); CHKERRQ(ierr); 559 { 560 // Create global-to-global scatter for Dirichlet values (everything not in 561 // ltogis0, which is the range of ltog0) 562 PetscInt xstart, xend, *indD, countD = 0; 563 IS isD; 564 const PetscScalar *x; 565 ierr = VecZeroEntries(Xloc); CHKERRQ(ierr); 566 ierr = VecSet(X, 1.0); CHKERRQ(ierr); 567 ierr = VecScatterBegin(ltog0, Xloc, X, INSERT_VALUES, SCATTER_FORWARD); 568 CHKERRQ(ierr); 569 ierr = VecScatterEnd(ltog0, Xloc, X, INSERT_VALUES, SCATTER_FORWARD); 570 CHKERRQ(ierr); 571 ierr = VecGetOwnershipRange(X, &xstart, &xend); CHKERRQ(ierr); 572 ierr = PetscMalloc1(xend-xstart, &indD); CHKERRQ(ierr); 573 ierr = VecGetArrayRead(X, &x); CHKERRQ(ierr); 574 for (PetscInt i=0; i<xend-xstart; i++) { 575 if (x[i] == 1.) indD[countD++] = xstart + i; 576 } 577 ierr = VecRestoreArrayRead(X, &x); CHKERRQ(ierr); 578 ierr = ISCreateGeneral(comm, countD, indD, PETSC_COPY_VALUES, &isD); 579 CHKERRQ(ierr); 580 ierr = PetscFree(indD); CHKERRQ(ierr); 581 ierr = VecScatterCreate(X, isD, X, isD, >ogD); CHKERRQ(ierr); 582 ierr = ISDestroy(&isD); CHKERRQ(ierr); 583 } 584 ierr = ISDestroy(<ogis); CHKERRQ(ierr); 585 ierr = ISDestroy(<ogis0); CHKERRQ(ierr); 586 ierr = ISDestroy(&locis); CHKERRQ(ierr); 587 } 588 589 // Set up libCEED 590 CeedInit(ceedresource, &ceed); 591 592 // CEED bases 593 CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompu, P, Q, 594 bpOptions[bpChoice].qmode, &basisu); 595 CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, Q, 596 bpOptions[bpChoice].qmode, &basisx); 597 598 // CEED restrictions 599 CreateRestriction(ceed, melem, P, ncompu, &Erestrictu); 600 CreateRestriction(ceed, melem, 2, dim, &Erestrictx); 601 CeedInt nelem = melem[0]*melem[1]*melem[2]; 602 CeedElemRestrictionCreateIdentity(ceed, nelem, Q*Q*Q, nelem*Q*Q*Q, ncompu, 603 &Erestrictui); 604 CeedElemRestrictionCreateIdentity(ceed, nelem, 605 Q*Q*Q, 606 nelem*Q*Q*Q, 607 bpOptions[bpChoice].qdatasize, &Erestrictqdi); 608 CeedElemRestrictionCreateIdentity(ceed, nelem, Q*Q*Q, nelem*Q*Q*Q, 1, 609 &Erestrictxi); 610 { 611 CeedScalar *xloc; 612 CeedInt shape[3] = {melem[0]+1, melem[1]+1, melem[2]+1}, len = 613 shape[0]*shape[1]*shape[2]; 614 xloc = malloc(len*ncompx*sizeof xloc[0]); 615 for (CeedInt i=0; i<shape[0]; i++) { 616 for (CeedInt j=0; j<shape[1]; j++) { 617 for (CeedInt k=0; k<shape[2]; k++) { 618 xloc[((i*shape[1]+j)*shape[2]+k) + 0*len] = 1.*(irank[0]*melem[0]+i) / 619 (p[0]*melem[0]); 620 xloc[((i*shape[1]+j)*shape[2]+k) + 1*len] = 1.*(irank[1]*melem[1]+j) / 621 (p[1]*melem[1]); 622 xloc[((i*shape[1]+j)*shape[2]+k) + 2*len] = 1.*(irank[2]*melem[2]+k) / 623 (p[2]*melem[2]); 624 } 625 } 626 } 627 CeedVectorCreate(ceed, len*ncompx, &xcoord); 628 CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_OWN_POINTER, xloc); 629 } 630 631 // Create the Qfunction that builds the operator quadrature data 632 CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].setupgeo, 633 bpOptions[bpChoice].setupgeofname, &qf_setupgeo); 634 CeedQFunctionAddInput(qf_setupgeo, "dx", ncompx*dim, CEED_EVAL_GRAD); 635 CeedQFunctionAddInput(qf_setupgeo, "weight", 1, CEED_EVAL_WEIGHT); 636 CeedQFunctionAddOutput(qf_setupgeo, "qdata", bpOptions[bpChoice].qdatasize, 637 CEED_EVAL_NONE); 638 639 // Create the Qfunction that sets up the RHS and true solution 640 CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].setuprhs, 641 bpOptions[bpChoice].setuprhsfname, &qf_setuprhs); 642 CeedQFunctionAddInput(qf_setuprhs, "x", ncompx, CEED_EVAL_INTERP); 643 CeedQFunctionAddInput(qf_setuprhs, "dx", ncompx*dim, CEED_EVAL_GRAD); 644 CeedQFunctionAddInput(qf_setuprhs, "weight", 1, CEED_EVAL_WEIGHT); 645 CeedQFunctionAddOutput(qf_setuprhs, "true_soln", ncompu, CEED_EVAL_NONE); 646 CeedQFunctionAddOutput(qf_setuprhs, "rhs", ncompu, CEED_EVAL_INTERP); 647 648 // Set up PDE operator 649 CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].apply, 650 bpOptions[bpChoice].applyfname, &qf_apply); 651 // Add inputs and outputs 652 CeedInt inscale = bpOptions[bpChoice].inmode==CEED_EVAL_GRAD ? 3 : 1; 653 CeedInt outscale = bpOptions[bpChoice].outmode==CEED_EVAL_GRAD ? 3 : 1; 654 CeedQFunctionAddInput(qf_apply, "u", ncompu*inscale, 655 bpOptions[bpChoice].inmode); 656 CeedQFunctionAddInput(qf_apply, "qdata", bpOptions[bpChoice].qdatasize, 657 CEED_EVAL_NONE); 658 CeedQFunctionAddOutput(qf_apply, "v", ncompu*outscale, 659 bpOptions[bpChoice].outmode); 660 661 // Create the error qfunction 662 CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].error, 663 bpOptions[bpChoice].errorfname, &qf_error); 664 CeedQFunctionAddInput(qf_error, "u", ncompu, CEED_EVAL_INTERP); 665 CeedQFunctionAddInput(qf_error, "true_soln", ncompu, CEED_EVAL_NONE); 666 CeedQFunctionAddOutput(qf_error, "error", ncompu, CEED_EVAL_NONE); 667 668 // Create the persistent vectors that will be needed in setup 669 CeedInt nqpts; 670 CeedBasisGetNumQuadraturePoints(basisu, &nqpts); 671 CeedVectorCreate(ceed, bpOptions[bpChoice].qdatasize*nelem*nqpts, &qdata); 672 CeedVectorCreate(ceed, nelem*nqpts*ncompu, &target); 673 CeedVectorCreate(ceed, lsize*ncompu, &rhsceed); 674 675 // Create the operator that builds the quadrature data for the ceed operator 676 CeedOperatorCreate(ceed, qf_setupgeo, NULL, NULL, &op_setupgeo); 677 CeedOperatorSetField(op_setupgeo, "dx", Erestrictx, CEED_NOTRANSPOSE, 678 basisx, CEED_VECTOR_ACTIVE); 679 CeedOperatorSetField(op_setupgeo, "weight", Erestrictxi, CEED_NOTRANSPOSE, 680 basisx, CEED_VECTOR_NONE); 681 CeedOperatorSetField(op_setupgeo, "qdata", Erestrictqdi, CEED_NOTRANSPOSE, 682 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 683 684 // Create the operator that builds the RHS and true solution 685 CeedOperatorCreate(ceed, qf_setuprhs, NULL, NULL, &op_setuprhs); 686 CeedOperatorSetField(op_setuprhs, "x", Erestrictx, CEED_NOTRANSPOSE, 687 basisx, CEED_VECTOR_ACTIVE); 688 CeedOperatorSetField(op_setuprhs, "dx", Erestrictx, CEED_NOTRANSPOSE, 689 basisx, CEED_VECTOR_ACTIVE); 690 CeedOperatorSetField(op_setuprhs, "weight", Erestrictxi, CEED_NOTRANSPOSE, 691 basisx, CEED_VECTOR_NONE); 692 CeedOperatorSetField(op_setuprhs, "true_soln", Erestrictui, CEED_NOTRANSPOSE, 693 CEED_BASIS_COLLOCATED, target); 694 CeedOperatorSetField(op_setuprhs, "rhs", Erestrictu, CEED_TRANSPOSE, 695 basisu, CEED_VECTOR_ACTIVE); 696 697 // Create the mass or diff operator 698 CeedOperatorCreate(ceed, qf_apply, NULL, NULL, &op_apply); 699 CeedOperatorSetField(op_apply, "u", Erestrictu, CEED_TRANSPOSE, 700 basisu, CEED_VECTOR_ACTIVE); 701 CeedOperatorSetField(op_apply, "qdata", Erestrictqdi, CEED_NOTRANSPOSE, 702 CEED_BASIS_COLLOCATED, qdata); 703 CeedOperatorSetField(op_apply, "v", Erestrictu, CEED_TRANSPOSE, 704 basisu, CEED_VECTOR_ACTIVE); 705 706 // Create the error operator 707 CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error); 708 CeedOperatorSetField(op_error, "u", Erestrictu, CEED_TRANSPOSE, 709 basisu, CEED_VECTOR_ACTIVE); 710 CeedOperatorSetField(op_error, "true_soln", Erestrictui, CEED_NOTRANSPOSE, 711 CEED_BASIS_COLLOCATED, target); 712 CeedOperatorSetField(op_error, "error", Erestrictui, CEED_NOTRANSPOSE, 713 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 714 715 // Set up Mat 716 ierr = PetscMalloc1(1, &user); CHKERRQ(ierr); 717 user->comm = comm; 718 user->ltog = ltog; 719 if (bpChoice != CEED_BP1 && bpChoice != CEED_BP2) { 720 user->ltog0 = ltog0; 721 user->gtogD = gtogD; 722 } 723 user->Xloc = Xloc; 724 ierr = VecDuplicate(Xloc, &user->Yloc); CHKERRQ(ierr); 725 CeedVectorCreate(ceed, lsize*ncompu, &user->xceed); 726 CeedVectorCreate(ceed, lsize*ncompu, &user->yceed); 727 user->op = op_apply; 728 user->qdata = qdata; 729 user->ceed = ceed; 730 731 ierr = MatCreateShell(comm, mnodes[0]*mnodes[1]*mnodes[2]*ncompu, 732 mnodes[0]*mnodes[1]*mnodes[2]*ncompu, 733 PETSC_DECIDE, PETSC_DECIDE, user, &mat); CHKERRQ(ierr); 734 if (bpChoice == CEED_BP1 || bpChoice == CEED_BP2) { 735 ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Mass); 736 CHKERRQ(ierr); 737 } else { 738 ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Diff); 739 CHKERRQ(ierr); 740 } 741 ierr = MatCreateVecs(mat, &rhs, NULL); CHKERRQ(ierr); 742 743 // Get RHS vector 744 ierr = VecDuplicate(Xloc, &rhsloc); CHKERRQ(ierr); 745 ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr); 746 ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr); 747 CeedVectorSetArray(rhsceed, CEED_MEM_HOST, CEED_USE_POINTER, r); 748 749 // Setup qdata, rhs, and target 750 CeedOperatorApply(op_setupgeo, xcoord, qdata, CEED_REQUEST_IMMEDIATE); 751 CeedOperatorApply(op_setuprhs, xcoord, rhsceed, CEED_REQUEST_IMMEDIATE); 752 ierr = CeedVectorSyncArray(rhsceed, CEED_MEM_HOST); CHKERRQ(ierr); 753 CeedVectorDestroy(&xcoord); 754 755 // Gather RHS 756 ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr); 757 ierr = VecZeroEntries(rhs); CHKERRQ(ierr); 758 ierr = VecScatterBegin(ltog, rhsloc, rhs, ADD_VALUES, SCATTER_FORWARD); 759 CHKERRQ(ierr); 760 ierr = VecScatterEnd(ltog, rhsloc, rhs, ADD_VALUES, SCATTER_FORWARD); 761 CHKERRQ(ierr); 762 CeedVectorDestroy(&rhsceed); 763 764 ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr); 765 { 766 PC pc; 767 ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 768 if (bpChoice == CEED_BP1 || bpChoice == CEED_BP2) { 769 ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr); 770 ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr); 771 } else { 772 ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr); 773 } 774 ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 775 ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 776 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 777 PETSC_DEFAULT); CHKERRQ(ierr); 778 } 779 ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr); 780 ierr = KSPSetOperators(ksp, mat, mat); CHKERRQ(ierr); 781 // First run, if benchmarking 782 if (benchmark_mode) { 783 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1); 784 CHKERRQ(ierr); 785 my_rt_start = MPI_Wtime(); 786 ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 787 my_rt = MPI_Wtime() - my_rt_start; 788 ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm); 789 CHKERRQ(ierr); 790 // Set maxits based on first iteration timing 791 if (my_rt > 0.02) { 792 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5); 793 CHKERRQ(ierr); 794 } else { 795 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20); 796 CHKERRQ(ierr); 797 } 798 } 799 // Timed solve 800 ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr); 801 my_rt_start = MPI_Wtime(); 802 ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 803 my_rt = MPI_Wtime() - my_rt_start; 804 { 805 KSPType ksptype; 806 KSPConvergedReason reason; 807 PetscReal rnorm; 808 PetscInt its; 809 ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr); 810 ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr); 811 ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr); 812 ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr); 813 if (!test_mode || reason < 0 || rnorm > 1e-8) { 814 ierr = PetscPrintf(comm, 815 " KSP:\n" 816 " KSP Type : %s\n" 817 " KSP Convergence : %s\n" 818 " Total KSP Iterations : %D\n" 819 " Final rnorm : %e\n", 820 ksptype, KSPConvergedReasons[reason], its, 821 (double)rnorm); CHKERRQ(ierr); 822 } 823 if (benchmark_mode && (!test_mode)) { 824 CeedInt gsize; 825 ierr = VecGetSize(X, &gsize); CHKERRQ(ierr); 826 ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm); 827 CHKERRQ(ierr); 828 ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm); 829 CHKERRQ(ierr); 830 ierr = PetscPrintf(comm, 831 " Performance:\n" 832 " CG Solve Time : %g (%g) sec\n" 833 " DoFs/Sec in CG : %g (%g) million\n", 834 rt_max, rt_min, 1e-6*gsize*its/rt_max, 835 1e-6*gsize*its/rt_min); CHKERRQ(ierr); 836 } 837 } 838 839 { 840 PetscReal maxerror; 841 ierr = ComputeErrorMax(user, op_error, X, target, &maxerror); CHKERRQ(ierr); 842 PetscReal tol = (bpChoice == CEED_BP1 || bpChoice == CEED_BP2) ? 5e-3 : 5e-2; 843 if (!test_mode || maxerror > tol) { 844 ierr = PetscPrintf(comm, 845 " Pointwise Error (max) : %e\n", 846 (double)maxerror); CHKERRQ(ierr); 847 } 848 } 849 850 if (write_solution) { 851 PetscViewer vtkviewersoln; 852 853 ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr); 854 ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr); 855 ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr); 856 ierr = VecView(X, vtkviewersoln); CHKERRQ(ierr); 857 ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr); 858 } 859 860 ierr = VecDestroy(&rhs); CHKERRQ(ierr); 861 ierr = VecDestroy(&rhsloc); CHKERRQ(ierr); 862 ierr = VecDestroy(&X); CHKERRQ(ierr); 863 ierr = VecDestroy(&user->Xloc); CHKERRQ(ierr); 864 ierr = VecDestroy(&user->Yloc); CHKERRQ(ierr); 865 ierr = VecScatterDestroy(<og); CHKERRQ(ierr); 866 ierr = VecScatterDestroy(<og0); CHKERRQ(ierr); 867 ierr = VecScatterDestroy(>ogD); CHKERRQ(ierr); 868 ierr = MatDestroy(&mat); CHKERRQ(ierr); 869 ierr = KSPDestroy(&ksp); CHKERRQ(ierr); 870 871 CeedVectorDestroy(&user->xceed); 872 CeedVectorDestroy(&user->yceed); 873 CeedVectorDestroy(&user->qdata); 874 CeedVectorDestroy(&target); 875 CeedOperatorDestroy(&op_setupgeo); 876 CeedOperatorDestroy(&op_setuprhs); 877 CeedOperatorDestroy(&op_apply); 878 CeedOperatorDestroy(&op_error); 879 CeedElemRestrictionDestroy(&Erestrictu); 880 CeedElemRestrictionDestroy(&Erestrictx); 881 CeedElemRestrictionDestroy(&Erestrictui); 882 CeedElemRestrictionDestroy(&Erestrictxi); 883 CeedElemRestrictionDestroy(&Erestrictqdi); 884 CeedQFunctionDestroy(&qf_setupgeo); 885 CeedQFunctionDestroy(&qf_setuprhs); 886 CeedQFunctionDestroy(&qf_apply); 887 CeedQFunctionDestroy(&qf_error); 888 CeedBasisDestroy(&basisu); 889 CeedBasisDestroy(&basisx); 890 CeedDestroy(&ceed); 891 ierr = PetscFree(user); CHKERRQ(ierr); 892 return PetscFinalize(); 893 } 894