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