xref: /libCEED/examples/petsc/bps.c (revision e07e9ddc8cd0eebb660553e99198f283c75572d0)
10c59ef15SJeremy L Thompson //                        libCEED + PETSc Example: CEED BPs
20c59ef15SJeremy L Thompson //
30c59ef15SJeremy L Thompson // This example demonstrates a simple usage of libCEED with PETSc to solve the
40c59ef15SJeremy L Thompson // CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps.
50c59ef15SJeremy L Thompson //
60c59ef15SJeremy L Thompson // The code is intentionally "raw", using only low-level communication
70c59ef15SJeremy L Thompson // primitives.
80c59ef15SJeremy L Thompson //
90c59ef15SJeremy L Thompson // Build with:
100c59ef15SJeremy L Thompson //
110c59ef15SJeremy L Thompson //     make bps [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
120c59ef15SJeremy L Thompson //
130c59ef15SJeremy L Thompson // Sample runs:
140c59ef15SJeremy L Thompson //
150c59ef15SJeremy L Thompson //     bps -problem bp1
160c59ef15SJeremy L Thompson //     bps -problem bp2 -ceed /cpu/self
170c59ef15SJeremy L Thompson //     bps -problem bp3 -ceed /gpu/occa
180c59ef15SJeremy L Thompson //     bps -problem bp4 -ceed /cpu/occa
190c59ef15SJeremy L Thompson //     bps -problem bp5 -ceed /omp/occa
200c59ef15SJeremy L Thompson //     bps -problem bp6 -ceed /ocl/occa
210c59ef15SJeremy L Thompson //
220c59ef15SJeremy L Thompson //TESTARGS -ceed {ceed_resource} -test -problem bp3
230c59ef15SJeremy L Thompson 
240c59ef15SJeremy L Thompson /// @file
250c59ef15SJeremy L Thompson /// CEED BPs example using PETSc
260c59ef15SJeremy L Thompson const char help[] = "Solve CEED BPs using PETSc\n";
270c59ef15SJeremy L Thompson 
280c59ef15SJeremy L Thompson #include <stdbool.h>
290c59ef15SJeremy L Thompson #include <string.h>
300c59ef15SJeremy L Thompson #include "common.h"
310c59ef15SJeremy L Thompson #include "bp1.h"
320c59ef15SJeremy L Thompson #include "bp2.h"
330c59ef15SJeremy L Thompson #include "bp3.h"
340c59ef15SJeremy L Thompson #include "bp4.h"
350c59ef15SJeremy L Thompson 
360c59ef15SJeremy L Thompson #define PATH(BASE) __DIR__ #BASE
370c59ef15SJeremy L Thompson 
380c59ef15SJeremy L Thompson static void Split3(PetscInt size, PetscInt m[3], bool reverse) {
390c59ef15SJeremy L Thompson   for (PetscInt d=0,sizeleft=size; d<3; d++) {
400c59ef15SJeremy L Thompson     PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(sizeleft, 1./(3 - d)));
410c59ef15SJeremy L Thompson     while (try * (sizeleft / try) != sizeleft) try++;
420c59ef15SJeremy L Thompson     m[reverse ? 2-d : d] = try;
430c59ef15SJeremy L Thompson     sizeleft /= try;
440c59ef15SJeremy L Thompson   }
450c59ef15SJeremy L Thompson }
460c59ef15SJeremy L Thompson 
470c59ef15SJeremy L Thompson static PetscInt Max3(const PetscInt a[3]) {
480c59ef15SJeremy L Thompson   return PetscMax(a[0], PetscMax(a[1], a[2]));
490c59ef15SJeremy L Thompson }
500c59ef15SJeremy L Thompson static PetscInt Min3(const PetscInt a[3]) {
510c59ef15SJeremy L Thompson   return PetscMin(a[0], PetscMin(a[1], a[2]));
520c59ef15SJeremy L Thompson }
530c59ef15SJeremy L Thompson static void GlobalDof(const PetscInt p[3], const PetscInt irank[3],
540c59ef15SJeremy L Thompson                       PetscInt degree, const PetscInt melem[3],
550c59ef15SJeremy L Thompson                       PetscInt mdof[3]) {
560c59ef15SJeremy L Thompson   for (int d=0; d<3; d++)
570c59ef15SJeremy L Thompson     mdof[d] = degree*melem[d] + (irank[d] == p[d]-1);
580c59ef15SJeremy L Thompson }
590c59ef15SJeremy L Thompson static PetscInt GlobalStart(const PetscInt p[3], const PetscInt irank[3],
600c59ef15SJeremy L Thompson                             PetscInt degree, const PetscInt melem[3]) {
610c59ef15SJeremy L Thompson   PetscInt start = 0;
620c59ef15SJeremy L Thompson   // Dumb brute-force is easier to read
630c59ef15SJeremy L Thompson   for (PetscInt i=0; i<p[0]; i++) {
640c59ef15SJeremy L Thompson     for (PetscInt j=0; j<p[1]; j++) {
650c59ef15SJeremy L Thompson       for (PetscInt k=0; k<p[2]; k++) {
660c59ef15SJeremy L Thompson         PetscInt mdof[3], ijkrank[] = {i,j,k};
670c59ef15SJeremy L Thompson         if (i == irank[0] && j == irank[1] && k == irank[2]) return start;
680c59ef15SJeremy L Thompson         GlobalDof(p, ijkrank, degree, melem, mdof);
690c59ef15SJeremy L Thompson         start += mdof[0] * mdof[1] * mdof[2];
700c59ef15SJeremy L Thompson       }
710c59ef15SJeremy L Thompson     }
720c59ef15SJeremy L Thompson   }
730c59ef15SJeremy L Thompson   return -1;
740c59ef15SJeremy L Thompson }
750c59ef15SJeremy L Thompson static int CreateRestriction(Ceed ceed, const CeedInt melem[3],
760c59ef15SJeremy L Thompson                              CeedInt P, CeedInt ncomp,
770c59ef15SJeremy L Thompson                              CeedElemRestriction *Erestrict) {
780c59ef15SJeremy L Thompson   const PetscInt Nelem = melem[0]*melem[1]*melem[2];
790c59ef15SJeremy L Thompson   PetscInt mdof[3], *idx, *idxp;
800c59ef15SJeremy L Thompson 
810c59ef15SJeremy L Thompson   for (int d=0; d<3; d++) mdof[d] = melem[d]*(P-1) + 1;
820c59ef15SJeremy L Thompson   idxp = idx = malloc(Nelem*P*P*P*sizeof idx[0]);
830c59ef15SJeremy L Thompson   for (CeedInt i=0; i<melem[0]; i++) {
840c59ef15SJeremy L Thompson     for (CeedInt j=0; j<melem[1]; j++) {
850c59ef15SJeremy L Thompson       for (CeedInt k=0; k<melem[2]; k++,idxp += P*P*P) {
860c59ef15SJeremy L Thompson         for (CeedInt ii=0; ii<P; ii++) {
870c59ef15SJeremy L Thompson           for (CeedInt jj=0; jj<P; jj++) {
880c59ef15SJeremy L Thompson             for (CeedInt kk=0; kk<P; kk++) {
890c59ef15SJeremy L Thompson               if (0) { // This is the C-style (i,j,k) ordering that I prefer
900c59ef15SJeremy L Thompson                 idxp[(ii*P+jj)*P+kk] = (((i*(P-1)+ii)*mdof[1]
910c59ef15SJeremy L Thompson                                          + (j*(P-1)+jj))*mdof[2]
920c59ef15SJeremy L Thompson                                         + (k*(P-1)+kk));
930c59ef15SJeremy L Thompson               } else { // (k,j,i) ordering for consistency with MFEM example
940c59ef15SJeremy L Thompson                 idxp[ii+P*(jj+P*kk)] = (((i*(P-1)+ii)*mdof[1]
950c59ef15SJeremy L Thompson                                          + (j*(P-1)+jj))*mdof[2]
960c59ef15SJeremy L Thompson                                         + (k*(P-1)+kk));
970c59ef15SJeremy L Thompson               }
980c59ef15SJeremy L Thompson             }
990c59ef15SJeremy L Thompson           }
1000c59ef15SJeremy L Thompson         }
1010c59ef15SJeremy L Thompson       }
1020c59ef15SJeremy L Thompson     }
1030c59ef15SJeremy L Thompson   }
1040c59ef15SJeremy L Thompson   CeedElemRestrictionCreate(ceed, Nelem, P*P*P, mdof[0]*mdof[1]*mdof[2], ncomp,
1050c59ef15SJeremy L Thompson                             CEED_MEM_HOST, CEED_OWN_POINTER, idx, Erestrict);
1060c59ef15SJeremy L Thompson   PetscFunctionReturn(0);
1070c59ef15SJeremy L Thompson }
1080c59ef15SJeremy L Thompson 
1090c59ef15SJeremy L Thompson // Data for PETSc
1100c59ef15SJeremy L Thompson typedef struct User_ *User;
1110c59ef15SJeremy L Thompson struct User_ {
1120c59ef15SJeremy L Thompson   MPI_Comm comm;
1130c59ef15SJeremy L Thompson   VecScatter ltog;              // Scatter for all entries
1140c59ef15SJeremy L Thompson   VecScatter ltog0;             // Skip Dirichlet values
1150c59ef15SJeremy L Thompson   VecScatter gtogD;             // global-to-global; only Dirichlet values
1160c59ef15SJeremy L Thompson   Vec Xloc, Yloc;
1170c59ef15SJeremy L Thompson   CeedVector xceed, yceed;
1180c59ef15SJeremy L Thompson   CeedOperator op;
1190c59ef15SJeremy L Thompson   CeedVector rho;
1200c59ef15SJeremy L Thompson   Ceed ceed;
1210c59ef15SJeremy L Thompson };
1220c59ef15SJeremy L Thompson 
1230c59ef15SJeremy L Thompson // BP Options
1240c59ef15SJeremy L Thompson typedef enum {
1250c59ef15SJeremy L Thompson   CEED_BP1 = 0, CEED_BP2 = 1, CEED_BP3 = 2,
1260c59ef15SJeremy L Thompson   CEED_BP4 = 3, CEED_BP5 = 4, CEED_BP6 = 5
1270c59ef15SJeremy L Thompson } bpType;
1280c59ef15SJeremy L Thompson static const char *const bpTypes[] = {"bp1","bp2","bp3","bp4","bp5","bp6",
1290c59ef15SJeremy L Thompson                                       "bpType","CEED_BP",0};
1300c59ef15SJeremy L Thompson 
1310c59ef15SJeremy L Thompson // BP specific data
1320c59ef15SJeremy L Thompson typedef struct {
1330c59ef15SJeremy L Thompson   CeedInt vscale, qdatasize, qextra;
1340c59ef15SJeremy L Thompson   CeedQFunctionUser setup, apply, error;
1350c59ef15SJeremy L Thompson   const char setupfname[PETSC_MAX_PATH_LEN], applyfname[PETSC_MAX_PATH_LEN],
1360c59ef15SJeremy L Thompson              errorfname[PETSC_MAX_PATH_LEN];
1370c59ef15SJeremy L Thompson   CeedEvalMode inmode, outmode;
138*e07e9ddcSjeremylt   CeedQuadMode qmode;
1390c59ef15SJeremy L Thompson } bpData;
1400c59ef15SJeremy L Thompson 
1410c59ef15SJeremy L Thompson bpData bpOptions[6] = {
1420c59ef15SJeremy L Thompson   [CEED_BP1] = {
1430c59ef15SJeremy L Thompson     .vscale = 1,
1440c59ef15SJeremy L Thompson     .qdatasize = 1,
145194c25f7Sjeremylt     .qextra = 1,
1460c59ef15SJeremy L Thompson     .setup = SetupMass,
1470c59ef15SJeremy L Thompson     .apply = Mass,
1480c59ef15SJeremy L Thompson     .error = Error,
1490c59ef15SJeremy L Thompson     .setupfname = PATH(bp1.h:SetupMass),
1500c59ef15SJeremy L Thompson     .applyfname = PATH(bp1.h:Mass),
1510c59ef15SJeremy L Thompson     .errorfname = PATH(common.h:Error),
1520c59ef15SJeremy L Thompson     .inmode = CEED_EVAL_INTERP,
153*e07e9ddcSjeremylt     .outmode = CEED_EVAL_INTERP,
154*e07e9ddcSjeremylt     .qmode = CEED_GAUSS},
1550c59ef15SJeremy L Thompson   [CEED_BP2] = {
1560c59ef15SJeremy L Thompson     .vscale = 3,
1570c59ef15SJeremy L Thompson     .qdatasize = 1,
158194c25f7Sjeremylt     .qextra = 1,
1590c59ef15SJeremy L Thompson     .setup = SetupMass3,
1600c59ef15SJeremy L Thompson     .apply = Mass3,
1610c59ef15SJeremy L Thompson     .error = Error3,
1620c59ef15SJeremy L Thompson     .setupfname = PATH(bp2.h:SetupMass3),
1630c59ef15SJeremy L Thompson     .applyfname = PATH(bp2.h:Mass3),
1640c59ef15SJeremy L Thompson     .errorfname = PATH(common.h:Error3),
1650c59ef15SJeremy L Thompson     .inmode = CEED_EVAL_INTERP,
166*e07e9ddcSjeremylt     .outmode = CEED_EVAL_INTERP,
167*e07e9ddcSjeremylt     .qmode = CEED_GAUSS},
1680c59ef15SJeremy L Thompson   [CEED_BP3] = {
1690c59ef15SJeremy L Thompson     .vscale = 1,
1700c59ef15SJeremy L Thompson     .qdatasize = 6,
171194c25f7Sjeremylt     .qextra = 1,
1720c59ef15SJeremy L Thompson     .setup = SetupDiff,
1730c59ef15SJeremy L Thompson     .apply = Diff,
1740c59ef15SJeremy L Thompson     .error = Error,
1750c59ef15SJeremy L Thompson     .setupfname = PATH(bp3.h:SetupDiff),
1760c59ef15SJeremy L Thompson     .applyfname = PATH(bp3.h:Diff),
1770c59ef15SJeremy L Thompson     .errorfname = PATH(common.h:Error),
1780c59ef15SJeremy L Thompson     .inmode = CEED_EVAL_GRAD,
179*e07e9ddcSjeremylt     .outmode = CEED_EVAL_GRAD,
180*e07e9ddcSjeremylt     .qmode = CEED_GAUSS},
1810c59ef15SJeremy L Thompson   [CEED_BP4] = {
1820c59ef15SJeremy L Thompson     .vscale = 3,
1830c59ef15SJeremy L Thompson     .qdatasize = 6,
184194c25f7Sjeremylt     .qextra = 1,
1850c59ef15SJeremy L Thompson     .setup = SetupDiff3,
1860c59ef15SJeremy L Thompson     .apply = Diff3,
1874b5b4ec1Sjeremylt     .error = Error3,
1880c59ef15SJeremy L Thompson     .setupfname = PATH(bp4.h:SetupDiff3),
1890c59ef15SJeremy L Thompson     .applyfname = PATH(bp4.h:Diff),
1900c59ef15SJeremy L Thompson     .errorfname = PATH(common.h:Error3),
1910c59ef15SJeremy L Thompson     .inmode = CEED_EVAL_GRAD,
192*e07e9ddcSjeremylt     .outmode = CEED_EVAL_GRAD,
193*e07e9ddcSjeremylt     .qmode = CEED_GAUSS},
1940c59ef15SJeremy L Thompson   [CEED_BP5] = {
1950c59ef15SJeremy L Thompson     .vscale = 1,
1960c59ef15SJeremy L Thompson     .qdatasize = 6,
197194c25f7Sjeremylt     .qextra = 0,
1980c59ef15SJeremy L Thompson     .setup = SetupDiff,
1990c59ef15SJeremy L Thompson     .apply = Diff,
2000c59ef15SJeremy L Thompson     .error = Error,
2010c59ef15SJeremy L Thompson     .setupfname = PATH(bp3.h:SetupDiff),
2020c59ef15SJeremy L Thompson     .applyfname = PATH(bp3.h:Diff),
2030c59ef15SJeremy L Thompson     .errorfname = PATH(common.h:Error),
2040c59ef15SJeremy L Thompson     .inmode = CEED_EVAL_GRAD,
205*e07e9ddcSjeremylt     .outmode = CEED_EVAL_GRAD,
206*e07e9ddcSjeremylt     .qmode = CEED_GAUSS_LOBATTO},
2070c59ef15SJeremy L Thompson   [CEED_BP6] = {
2080c59ef15SJeremy L Thompson     .vscale = 3,
2090c59ef15SJeremy L Thompson     .qdatasize = 6,
210194c25f7Sjeremylt     .qextra = 0,
2110c59ef15SJeremy L Thompson     .setup = SetupDiff3,
2120c59ef15SJeremy L Thompson     .apply = Diff3,
2134b5b4ec1Sjeremylt     .error = Error3,
2140c59ef15SJeremy L Thompson     .setupfname = PATH(bp4.h:SetupDiff3),
2150c59ef15SJeremy L Thompson     .applyfname = PATH(bp4.h:Diff),
2160c59ef15SJeremy L Thompson     .errorfname = PATH(common.h:Error3),
2170c59ef15SJeremy L Thompson     .inmode = CEED_EVAL_GRAD,
218*e07e9ddcSjeremylt     .outmode = CEED_EVAL_GRAD,
219*e07e9ddcSjeremylt     .qmode = CEED_GAUSS_LOBATTO}
2200c59ef15SJeremy L Thompson };
2210c59ef15SJeremy L Thompson 
2220c59ef15SJeremy L Thompson // This function uses libCEED to compute the action of the mass matrix
2230c59ef15SJeremy L Thompson static PetscErrorCode MatMult_Mass(Mat A, Vec X, Vec Y) {
2240c59ef15SJeremy L Thompson   PetscErrorCode ierr;
2250c59ef15SJeremy L Thompson   User user;
2260c59ef15SJeremy L Thompson   PetscScalar *x, *y;
2270c59ef15SJeremy L Thompson 
2280c59ef15SJeremy L Thompson   PetscFunctionBeginUser;
2290c59ef15SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
2300c59ef15SJeremy L Thompson   ierr = VecScatterBegin(user->ltog, X, user->Xloc, INSERT_VALUES,
2310c59ef15SJeremy L Thompson                          SCATTER_REVERSE); CHKERRQ(ierr);
2320c59ef15SJeremy L Thompson   ierr = VecScatterEnd(user->ltog, X, user->Xloc, INSERT_VALUES, SCATTER_REVERSE);
2330c59ef15SJeremy L Thompson   CHKERRQ(ierr);
2340c59ef15SJeremy L Thompson   ierr = VecZeroEntries(user->Yloc); CHKERRQ(ierr);
2350c59ef15SJeremy L Thompson 
2360c59ef15SJeremy L Thompson   ierr = VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr);
2370c59ef15SJeremy L Thompson   ierr = VecGetArray(user->Yloc, &y); CHKERRQ(ierr);
2380c59ef15SJeremy L Thompson   CeedVectorSetArray(user->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
2390c59ef15SJeremy L Thompson   CeedVectorSetArray(user->yceed, CEED_MEM_HOST, CEED_USE_POINTER, y);
2400c59ef15SJeremy L Thompson 
2410c59ef15SJeremy L Thompson   CeedOperatorApply(user->op, user->xceed, user->yceed,
2420c59ef15SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
2430c59ef15SJeremy L Thompson   ierr = CeedVectorSyncArray(user->yceed, CEED_MEM_HOST); CHKERRQ(ierr);
2440c59ef15SJeremy L Thompson 
2450c59ef15SJeremy L Thompson   ierr = VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr);
2460c59ef15SJeremy L Thompson   ierr = VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr);
2470c59ef15SJeremy L Thompson 
2480c59ef15SJeremy L Thompson   if (Y) {
2490c59ef15SJeremy L Thompson     ierr = VecZeroEntries(Y); CHKERRQ(ierr);
2500c59ef15SJeremy L Thompson     ierr = VecScatterBegin(user->ltog, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD);
2510c59ef15SJeremy L Thompson     CHKERRQ(ierr);
2520c59ef15SJeremy L Thompson     ierr = VecScatterEnd(user->ltog, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD);
2530c59ef15SJeremy L Thompson     CHKERRQ(ierr);
2540c59ef15SJeremy L Thompson   }
2550c59ef15SJeremy L Thompson   PetscFunctionReturn(0);
2560c59ef15SJeremy L Thompson }
2570c59ef15SJeremy L Thompson 
2580c59ef15SJeremy L Thompson // This function uses libCEED to compute the action of the Laplacian with
2590c59ef15SJeremy L Thompson // Dirichlet boundary conditions
2600c59ef15SJeremy L Thompson static PetscErrorCode MatMult_Diff(Mat A, Vec X, Vec Y) {
2610c59ef15SJeremy L Thompson   PetscErrorCode ierr;
2620c59ef15SJeremy L Thompson   User user;
2630c59ef15SJeremy L Thompson   PetscScalar *x, *y;
2640c59ef15SJeremy L Thompson 
2650c59ef15SJeremy L Thompson   PetscFunctionBeginUser;
2660c59ef15SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
2670c59ef15SJeremy L Thompson   ierr = VecScatterBegin(user->ltog0, X, user->Xloc, INSERT_VALUES,
2680c59ef15SJeremy L Thompson                          SCATTER_REVERSE); CHKERRQ(ierr);
2690c59ef15SJeremy L Thompson   ierr = VecScatterEnd(user->ltog0, X, user->Xloc, INSERT_VALUES,
2700c59ef15SJeremy L Thompson                        SCATTER_REVERSE);
2710c59ef15SJeremy L Thompson   CHKERRQ(ierr);
2720c59ef15SJeremy L Thompson   ierr = VecZeroEntries(user->Yloc); CHKERRQ(ierr);
2730c59ef15SJeremy L Thompson 
2740c59ef15SJeremy L Thompson   ierr = VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr);
2750c59ef15SJeremy L Thompson   ierr = VecGetArray(user->Yloc, &y); CHKERRQ(ierr);
2760c59ef15SJeremy L Thompson   CeedVectorSetArray(user->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
2770c59ef15SJeremy L Thompson   CeedVectorSetArray(user->yceed, CEED_MEM_HOST, CEED_USE_POINTER, y);
2780c59ef15SJeremy L Thompson 
2790c59ef15SJeremy L Thompson   CeedOperatorApply(user->op, user->xceed, user->yceed,
2800c59ef15SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
2810c59ef15SJeremy L Thompson   ierr = CeedVectorSyncArray(user->yceed, CEED_MEM_HOST); CHKERRQ(ierr);
2820c59ef15SJeremy L Thompson 
2830c59ef15SJeremy L Thompson   ierr = VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr);
2840c59ef15SJeremy L Thompson   ierr = VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr);
2850c59ef15SJeremy L Thompson 
2860c59ef15SJeremy L Thompson   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
2870c59ef15SJeremy L Thompson   ierr = VecScatterBegin(user->gtogD, X, Y, INSERT_VALUES, SCATTER_FORWARD);
2880c59ef15SJeremy L Thompson   CHKERRQ(ierr);
2890c59ef15SJeremy L Thompson   ierr = VecScatterEnd(user->gtogD, X, Y, INSERT_VALUES, SCATTER_FORWARD);
2900c59ef15SJeremy L Thompson   CHKERRQ(ierr);
2910c59ef15SJeremy L Thompson   ierr = VecScatterBegin(user->ltog0, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD);
2920c59ef15SJeremy L Thompson   CHKERRQ(ierr);
2930c59ef15SJeremy L Thompson   ierr = VecScatterEnd(user->ltog0, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD);
2940c59ef15SJeremy L Thompson   CHKERRQ(ierr);
2950c59ef15SJeremy L Thompson   PetscFunctionReturn(0);
2960c59ef15SJeremy L Thompson }
2970c59ef15SJeremy L Thompson 
2980c59ef15SJeremy L Thompson // This function calculates the error in the final solution
2990c59ef15SJeremy L Thompson static PetscErrorCode ComputeErrorMax(User user, CeedOperator op_error, Vec X,
3000c59ef15SJeremy L Thompson                                       CeedVector target, PetscReal *maxerror) {
3010c59ef15SJeremy L Thompson   PetscErrorCode ierr;
3020c59ef15SJeremy L Thompson   PetscScalar *x;
3030c59ef15SJeremy L Thompson   CeedVector collocated_error;
3040c59ef15SJeremy L Thompson   CeedInt length;
3050c59ef15SJeremy L Thompson 
3060c59ef15SJeremy L Thompson   PetscFunctionBeginUser;
3070c59ef15SJeremy L Thompson   CeedVectorGetLength(target, &length);
3080c59ef15SJeremy L Thompson   CeedVectorCreate(user->ceed, length, &collocated_error);
3090c59ef15SJeremy L Thompson   ierr = VecScatterBegin(user->ltog, X, user->Xloc, INSERT_VALUES,
3100c59ef15SJeremy L Thompson                          SCATTER_REVERSE); CHKERRQ(ierr);
3110c59ef15SJeremy L Thompson   ierr = VecScatterEnd(user->ltog, X, user->Xloc, INSERT_VALUES, SCATTER_REVERSE);
3120c59ef15SJeremy L Thompson   CHKERRQ(ierr);
3130c59ef15SJeremy L Thompson   ierr = VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr);
3140c59ef15SJeremy L Thompson   CeedVectorSetArray(user->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
3150c59ef15SJeremy L Thompson   CeedOperatorApply(op_error, user->xceed, collocated_error,
3160c59ef15SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
3170c59ef15SJeremy L Thompson   VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr);
3180c59ef15SJeremy L Thompson 
3190c59ef15SJeremy L Thompson   *maxerror = 0;
3200c59ef15SJeremy L Thompson   const CeedScalar *e;
3210c59ef15SJeremy L Thompson   CeedVectorGetArrayRead(collocated_error, CEED_MEM_HOST, &e);
3220c59ef15SJeremy L Thompson   for (CeedInt i=0; i<length; i++) {
3230c59ef15SJeremy L Thompson     *maxerror = PetscMax(*maxerror, PetscAbsScalar(e[i]));
3240c59ef15SJeremy L Thompson   }
3250c59ef15SJeremy L Thompson   CeedVectorRestoreArrayRead(collocated_error, &e);
3260c59ef15SJeremy L Thompson   ierr = MPI_Allreduce(MPI_IN_PLACE, &maxerror,
3270c59ef15SJeremy L Thompson                        1, MPIU_SCALAR, MPIU_MAX, user->comm); CHKERRQ(ierr);
3280c59ef15SJeremy L Thompson   CeedVectorDestroy(&collocated_error);
3290c59ef15SJeremy L Thompson   PetscFunctionReturn(0);
3300c59ef15SJeremy L Thompson }
3310c59ef15SJeremy L Thompson 
3320c59ef15SJeremy L Thompson int main(int argc, char **argv) {
3330c59ef15SJeremy L Thompson   PetscInt ierr;
3340c59ef15SJeremy L Thompson   MPI_Comm comm;
3350c59ef15SJeremy L Thompson   char ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self";
3360c59ef15SJeremy L Thompson   PetscInt degree, qextra, localdof, localelem, melem[3], mdof[3], p[3],
3370c59ef15SJeremy L Thompson            irank[3], ldof[3], lsize, vscale = 1;
3380c59ef15SJeremy L Thompson   PetscScalar *r;
3390c59ef15SJeremy L Thompson   PetscBool test_mode, benchmark_mode;
3400c59ef15SJeremy L Thompson   PetscMPIInt size, rank;
3414b5b4ec1Sjeremylt   VecScatter ltog, ltog0, gtogD;
3420c59ef15SJeremy L Thompson   Ceed ceed;
3430c59ef15SJeremy L Thompson   CeedBasis basisx, basisu;
3440c59ef15SJeremy L Thompson   CeedElemRestriction Erestrictx, Erestrictu, Erestrictxi, Erestrictui,
3450c59ef15SJeremy L Thompson                       Erestrictqdi;
3460c59ef15SJeremy L Thompson   CeedQFunction qf_setup, qf_apply, qf_error;
3470c59ef15SJeremy L Thompson   CeedOperator op_setup, op_apply, op_error;
3480c59ef15SJeremy L Thompson   CeedVector xcoord, rho, rhsceed, target;
3490c59ef15SJeremy L Thompson   CeedInt P, Q;
3500c59ef15SJeremy L Thompson   Vec X, Xloc, rhs, rhsloc;
3510c59ef15SJeremy L Thompson   Mat mat;
3520c59ef15SJeremy L Thompson   KSP ksp;
3530c59ef15SJeremy L Thompson   User user;
3540c59ef15SJeremy L Thompson   double my_rt_start, my_rt, rt_min, rt_max;
3550c59ef15SJeremy L Thompson   bpType bpChoice;
3560c59ef15SJeremy L Thompson 
3570c59ef15SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
3580c59ef15SJeremy L Thompson   if (ierr) return ierr;
3590c59ef15SJeremy L Thompson   comm = PETSC_COMM_WORLD;
3600c59ef15SJeremy L Thompson   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr);
3610c59ef15SJeremy L Thompson   bpChoice = CEED_BP1;
3620c59ef15SJeremy L Thompson   ierr = PetscOptionsEnum("-problem",
3630c59ef15SJeremy L Thompson                           "CEED benchmark problem to solve", NULL,
3640c59ef15SJeremy L Thompson                           bpTypes, (PetscEnum)bpChoice, (PetscEnum*)&bpChoice,
3650c59ef15SJeremy L Thompson                           NULL); CHKERRQ(ierr);
3660c59ef15SJeremy L Thompson   vscale = bpOptions[bpChoice].vscale;
3670c59ef15SJeremy L Thompson   test_mode = PETSC_FALSE;
3680c59ef15SJeremy L Thompson   ierr = PetscOptionsBool("-test",
3690c59ef15SJeremy L Thompson                           "Testing mode (do not print unless error is large)",
3700c59ef15SJeremy L Thompson                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
3710c59ef15SJeremy L Thompson   benchmark_mode = PETSC_FALSE;
3720c59ef15SJeremy L Thompson   ierr = PetscOptionsBool("-benchmark",
3730c59ef15SJeremy L Thompson                           "Benchmarking mode (prints benchmark statistics)",
3740c59ef15SJeremy L Thompson                           NULL, benchmark_mode, &benchmark_mode, NULL);
3750c59ef15SJeremy L Thompson   CHKERRQ(ierr);
3760c59ef15SJeremy L Thompson   degree = test_mode ? 3 : 1;
3770c59ef15SJeremy L Thompson   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
3780c59ef15SJeremy L Thompson                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
3790c59ef15SJeremy L Thompson   qextra = bpOptions[bpChoice].qextra;
3800c59ef15SJeremy L Thompson   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
3810c59ef15SJeremy L Thompson                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
3820c59ef15SJeremy L Thompson   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
3830c59ef15SJeremy L Thompson                             NULL, ceedresource, ceedresource,
3840c59ef15SJeremy L Thompson                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
3850c59ef15SJeremy L Thompson   localdof = 1000;
3860c59ef15SJeremy L Thompson   ierr = PetscOptionsInt("-local",
3870c59ef15SJeremy L Thompson                          "Target number of locally owned degrees of freedom per process",
3880c59ef15SJeremy L Thompson                          NULL, localdof, &localdof, NULL); CHKERRQ(ierr);
3890c59ef15SJeremy L Thompson   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
3903ce86546SJeremy L Thompson   P = degree + 1;
3913ce86546SJeremy L Thompson   Q = P + qextra;
3920c59ef15SJeremy L Thompson 
3930c59ef15SJeremy L Thompson   // Determine size of process grid
3940c59ef15SJeremy L Thompson   ierr = MPI_Comm_size(comm, &size); CHKERRQ(ierr);
3950c59ef15SJeremy L Thompson   Split3(size, p, false);
3960c59ef15SJeremy L Thompson 
3970c59ef15SJeremy L Thompson   // Find a nicely composite number of elements no less than localdof
3980c59ef15SJeremy L Thompson   for (localelem = PetscMax(1, localdof / (degree*degree*degree)); ;
3990c59ef15SJeremy L Thompson        localelem++) {
4000c59ef15SJeremy L Thompson     Split3(localelem, melem, true);
4010c59ef15SJeremy L Thompson     if (Max3(melem) / Min3(melem) <= 2) break;
4020c59ef15SJeremy L Thompson   }
4030c59ef15SJeremy L Thompson 
4040c59ef15SJeremy L Thompson   // Find my location in the process grid
4050c59ef15SJeremy L Thompson   ierr = MPI_Comm_rank(comm, &rank); CHKERRQ(ierr);
4060c59ef15SJeremy L Thompson   for (int d=0,rankleft=rank; d<3; d++) {
4070c59ef15SJeremy L Thompson     const int pstride[3] = {p[1] *p[2], p[2], 1};
4080c59ef15SJeremy L Thompson     irank[d] = rankleft / pstride[d];
4090c59ef15SJeremy L Thompson     rankleft -= irank[d] * pstride[d];
4100c59ef15SJeremy L Thompson   }
4110c59ef15SJeremy L Thompson 
4120c59ef15SJeremy L Thompson   GlobalDof(p, irank, degree, melem, mdof);
4130c59ef15SJeremy L Thompson 
4140c59ef15SJeremy L Thompson   ierr = VecCreate(comm, &X); CHKERRQ(ierr);
4150c59ef15SJeremy L Thompson   ierr = VecSetSizes(X, mdof[0]*mdof[1]*mdof[2]*vscale, PETSC_DECIDE);
4160c59ef15SJeremy L Thompson   CHKERRQ(ierr);
4170c59ef15SJeremy L Thompson   ierr = VecSetUp(X); CHKERRQ(ierr);
4180c59ef15SJeremy L Thompson 
4190c59ef15SJeremy L Thompson   if (!test_mode) {
4200c59ef15SJeremy L Thompson     CeedInt gsize;
4210c59ef15SJeremy L Thompson     ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
4223ce86546SJeremy L Thompson     ierr = PetscPrintf(comm,
4233ce86546SJeremy L Thompson                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n"
4243ce86546SJeremy L Thompson                        "  libCEED:\n"
4253ce86546SJeremy L Thompson                        "    libCEED Backend                    : %s\n"
4263ce86546SJeremy L Thompson                        "  Mesh:\n"
4273ce86546SJeremy L Thompson                        "    Number of 1D Basis Nodes (p)       : %d\n"
4283ce86546SJeremy L Thompson                        "    Number of 1D Quadrature Points (q) : %d\n"
4293ce86546SJeremy L Thompson                        "    Global DOFs                        : %D\n"
4303ce86546SJeremy L Thompson                        "    Process Decomposition              : %D %D %D\n"
4313ce86546SJeremy L Thompson                        "    Local Elements                     : %D = %D %D %D\n"
4323ce86546SJeremy L Thompson                        "    Owned DOFs                         : %D = %D %D %D\n",
4333ce86546SJeremy L Thompson                        bpChoice+1, ceedresource, P, Q,  gsize/vscale, p[0],
4343ce86546SJeremy L Thompson                        p[1], p[2], localelem, melem[0], melem[1], melem[2],
4353ce86546SJeremy L Thompson                        mdof[0]*mdof[1]*mdof[2], mdof[0], mdof[1], mdof[2]);
4363ce86546SJeremy L Thompson     CHKERRQ(ierr);
4370c59ef15SJeremy L Thompson   }
4380c59ef15SJeremy L Thompson 
4390c59ef15SJeremy L Thompson   {
4400c59ef15SJeremy L Thompson     lsize = 1;
4410c59ef15SJeremy L Thompson     for (int d=0; d<3; d++) {
4420c59ef15SJeremy L Thompson       ldof[d] = melem[d]*degree + 1;
4430c59ef15SJeremy L Thompson       lsize *= ldof[d];
4440c59ef15SJeremy L Thompson     }
4450c59ef15SJeremy L Thompson     ierr = VecCreate(PETSC_COMM_SELF, &Xloc); CHKERRQ(ierr);
4460c59ef15SJeremy L Thompson     ierr = VecSetSizes(Xloc, lsize*vscale, PETSC_DECIDE); CHKERRQ(ierr);
4470c59ef15SJeremy L Thompson     ierr = VecSetUp(Xloc); CHKERRQ(ierr);
4480c59ef15SJeremy L Thompson 
4490c59ef15SJeremy L Thompson     // Create local-to-global scatter
4500c59ef15SJeremy L Thompson     PetscInt *ltogind, *ltogind0, *locind, l0count;
4510c59ef15SJeremy L Thompson     IS ltogis, ltogis0, locis;
4520c59ef15SJeremy L Thompson     PetscInt gstart[2][2][2], gmdof[2][2][2][3];
4530c59ef15SJeremy L Thompson 
4540c59ef15SJeremy L Thompson     for (int i=0; i<2; i++) {
4550c59ef15SJeremy L Thompson       for (int j=0; j<2; j++) {
4560c59ef15SJeremy L Thompson         for (int k=0; k<2; k++) {
4570c59ef15SJeremy L Thompson           PetscInt ijkrank[3] = {irank[0]+i, irank[1]+j, irank[2]+k};
4580c59ef15SJeremy L Thompson           gstart[i][j][k] = GlobalStart(p, ijkrank, degree, melem);
4590c59ef15SJeremy L Thompson           GlobalDof(p, ijkrank, degree, melem, gmdof[i][j][k]);
4600c59ef15SJeremy L Thompson         }
4610c59ef15SJeremy L Thompson       }
4620c59ef15SJeremy L Thompson     }
4630c59ef15SJeremy L Thompson 
4640c59ef15SJeremy L Thompson     ierr = PetscMalloc1(lsize, &ltogind); CHKERRQ(ierr);
4650c59ef15SJeremy L Thompson     ierr = PetscMalloc1(lsize, &ltogind0); CHKERRQ(ierr);
4660c59ef15SJeremy L Thompson     ierr = PetscMalloc1(lsize, &locind); CHKERRQ(ierr);
4670c59ef15SJeremy L Thompson     l0count = 0;
4680c59ef15SJeremy L Thompson     for (PetscInt i=0,ir,ii; ir=i>=mdof[0], ii=i-ir*mdof[0], i<ldof[0]; i++) {
4690c59ef15SJeremy L Thompson       for (PetscInt j=0,jr,jj; jr=j>=mdof[1], jj=j-jr*mdof[1], j<ldof[1]; j++) {
4700c59ef15SJeremy L Thompson         for (PetscInt k=0,kr,kk; kr=k>=mdof[2], kk=k-kr*mdof[2], k<ldof[2]; k++) {
4710c59ef15SJeremy L Thompson           PetscInt here = (i*ldof[1]+j)*ldof[2]+k;
4720c59ef15SJeremy L Thompson           ltogind[here] =
4730c59ef15SJeremy L Thompson             gstart[ir][jr][kr] + (ii*gmdof[ir][jr][kr][1]+jj)*gmdof[ir][jr][kr][2]+kk;
4740c59ef15SJeremy L Thompson           if ((irank[0] == 0 && i == 0)
4750c59ef15SJeremy L Thompson               || (irank[1] == 0 && j == 0)
4760c59ef15SJeremy L Thompson               || (irank[2] == 0 && k == 0)
4770c59ef15SJeremy L Thompson               || (irank[0]+1 == p[0] && i+1 == ldof[0])
4780c59ef15SJeremy L Thompson               || (irank[1]+1 == p[1] && j+1 == ldof[1])
4790c59ef15SJeremy L Thompson               || (irank[2]+1 == p[2] && k+1 == ldof[2]))
4800c59ef15SJeremy L Thompson             continue;
4810c59ef15SJeremy L Thompson           ltogind0[l0count] = ltogind[here];
4820c59ef15SJeremy L Thompson           locind[l0count++] = here;
4830c59ef15SJeremy L Thompson         }
4840c59ef15SJeremy L Thompson       }
4850c59ef15SJeremy L Thompson     }
4860c59ef15SJeremy L Thompson     ierr = ISCreateBlock(comm, vscale, lsize, ltogind, PETSC_OWN_POINTER,
4870c59ef15SJeremy L Thompson                          &ltogis); CHKERRQ(ierr);
4880c59ef15SJeremy L Thompson     ierr = VecScatterCreate(Xloc, NULL, X, ltogis, &ltog); CHKERRQ(ierr);
4890c59ef15SJeremy L Thompson     CHKERRQ(ierr);
4900c59ef15SJeremy L Thompson     ierr = ISCreateBlock(comm, vscale, l0count, ltogind0, PETSC_OWN_POINTER,
4910c59ef15SJeremy L Thompson                          &ltogis0); CHKERRQ(ierr);
4920c59ef15SJeremy L Thompson     ierr = ISCreateBlock(comm, vscale, l0count, locind, PETSC_OWN_POINTER,
4930c59ef15SJeremy L Thompson                          &locis); CHKERRQ(ierr);
4940c59ef15SJeremy L Thompson     ierr = VecScatterCreate(Xloc, locis, X, ltogis0, &ltog0); CHKERRQ(ierr);
4950c59ef15SJeremy L Thompson     {
4960c59ef15SJeremy L Thompson       // Create global-to-global scatter for Dirichlet values (everything not in
4970c59ef15SJeremy L Thompson       // ltogis0, which is the range of ltog0)
4980c59ef15SJeremy L Thompson       PetscInt xstart, xend, *indD, countD = 0;
4990c59ef15SJeremy L Thompson       IS isD;
5000c59ef15SJeremy L Thompson       const PetscScalar *x;
5010c59ef15SJeremy L Thompson       ierr = VecZeroEntries(Xloc); CHKERRQ(ierr);
5020c59ef15SJeremy L Thompson       ierr = VecSet(X, 1.0); CHKERRQ(ierr);
5030c59ef15SJeremy L Thompson       ierr = VecScatterBegin(ltog0, Xloc, X, INSERT_VALUES, SCATTER_FORWARD);
5040c59ef15SJeremy L Thompson       CHKERRQ(ierr);
5050c59ef15SJeremy L Thompson       ierr = VecScatterEnd(ltog0, Xloc, X, INSERT_VALUES, SCATTER_FORWARD);
5060c59ef15SJeremy L Thompson       CHKERRQ(ierr);
5070c59ef15SJeremy L Thompson       ierr = VecGetOwnershipRange(X, &xstart, &xend); CHKERRQ(ierr);
5080c59ef15SJeremy L Thompson       ierr = PetscMalloc1(xend-xstart, &indD); CHKERRQ(ierr);
5090c59ef15SJeremy L Thompson       ierr = VecGetArrayRead(X, &x); CHKERRQ(ierr);
5100c59ef15SJeremy L Thompson       for (PetscInt i=0; i<xend-xstart; i++) {
5110c59ef15SJeremy L Thompson         if (x[i] == 1.) indD[countD++] = xstart + i;
5120c59ef15SJeremy L Thompson       }
5130c59ef15SJeremy L Thompson       ierr = VecRestoreArrayRead(X, &x); CHKERRQ(ierr);
5140c59ef15SJeremy L Thompson       ierr = ISCreateGeneral(comm, countD, indD, PETSC_COPY_VALUES, &isD);
5150c59ef15SJeremy L Thompson       CHKERRQ(ierr);
5160c59ef15SJeremy L Thompson       ierr = PetscFree(indD); CHKERRQ(ierr);
5170c59ef15SJeremy L Thompson       ierr = VecScatterCreate(X, isD, X, isD, &gtogD); CHKERRQ(ierr);
5180c59ef15SJeremy L Thompson       ierr = ISDestroy(&isD); CHKERRQ(ierr);
5190c59ef15SJeremy L Thompson     }
5200c59ef15SJeremy L Thompson     ierr = ISDestroy(&ltogis); CHKERRQ(ierr);
5210c59ef15SJeremy L Thompson     ierr = ISDestroy(&ltogis0); CHKERRQ(ierr);
5220c59ef15SJeremy L Thompson     ierr = ISDestroy(&locis); CHKERRQ(ierr);
5230c59ef15SJeremy L Thompson   }
5240c59ef15SJeremy L Thompson 
5250c59ef15SJeremy L Thompson   // Set up libCEED
5260c59ef15SJeremy L Thompson   CeedInit(ceedresource, &ceed);
527*e07e9ddcSjeremylt   CeedBasisCreateTensorH1Lagrange(ceed, 3, vscale, P, Q,
528*e07e9ddcSjeremylt                                   bpOptions[bpChoice].qmode, &basisu);
529*e07e9ddcSjeremylt   CeedBasisCreateTensorH1Lagrange(ceed, 3, 3, 2, Q,
530*e07e9ddcSjeremylt                                   bpOptions[bpChoice].qmode, &basisx);
5310c59ef15SJeremy L Thompson 
5320c59ef15SJeremy L Thompson   CreateRestriction(ceed, melem, P, vscale, &Erestrictu);
5330c59ef15SJeremy L Thompson   CreateRestriction(ceed, melem, 2, 3, &Erestrictx);
5340c59ef15SJeremy L Thompson   CeedInt nelem = melem[0]*melem[1]*melem[2];
5350c59ef15SJeremy L Thompson   CeedElemRestrictionCreateIdentity(ceed, nelem, Q*Q*Q, nelem*Q*Q*Q, vscale,
5360c59ef15SJeremy L Thompson                                     &Erestrictui);
5370c59ef15SJeremy L Thompson   CeedElemRestrictionCreateIdentity(ceed, nelem,
5380c59ef15SJeremy L Thompson                                     bpOptions[bpChoice].qdatasize*Q*Q*Q,
5390c59ef15SJeremy L Thompson                                     bpOptions[bpChoice].qdatasize*nelem*Q*Q*Q,
5400c59ef15SJeremy L Thompson                                     1, &Erestrictqdi);
5410c59ef15SJeremy L Thompson   CeedElemRestrictionCreateIdentity(ceed, nelem, Q*Q*Q, nelem*Q*Q*Q, 1,
5420c59ef15SJeremy L Thompson                                     &Erestrictxi);
5430c59ef15SJeremy L Thompson   {
5440c59ef15SJeremy L Thompson     CeedScalar *xloc;
5450c59ef15SJeremy L Thompson     CeedInt shape[3] = {melem[0]+1, melem[1]+1, melem[2]+1}, len =
5460c59ef15SJeremy L Thompson                          shape[0]*shape[1]*shape[2];
5470c59ef15SJeremy L Thompson     xloc = malloc(len*3*sizeof xloc[0]);
5480c59ef15SJeremy L Thompson     for (CeedInt i=0; i<shape[0]; i++) {
5490c59ef15SJeremy L Thompson       for (CeedInt j=0; j<shape[1]; j++) {
5500c59ef15SJeremy L Thompson         for (CeedInt k=0; k<shape[2]; k++) {
5510c59ef15SJeremy L Thompson           xloc[((i*shape[1]+j)*shape[2]+k) + 0*len] = 1.*(irank[0]*melem[0]+i) /
5520c59ef15SJeremy L Thompson               (p[0]*melem[0]);
5530c59ef15SJeremy L Thompson           xloc[((i*shape[1]+j)*shape[2]+k) + 1*len] = 1.*(irank[1]*melem[1]+j) /
5540c59ef15SJeremy L Thompson               (p[1]*melem[1]);
5550c59ef15SJeremy L Thompson           xloc[((i*shape[1]+j)*shape[2]+k) + 2*len] = 1.*(irank[2]*melem[2]+k) /
5560c59ef15SJeremy L Thompson               (p[2]*melem[2]);
5570c59ef15SJeremy L Thompson         }
5580c59ef15SJeremy L Thompson       }
5590c59ef15SJeremy L Thompson     }
5600c59ef15SJeremy L Thompson     CeedVectorCreate(ceed, len*3, &xcoord);
5610c59ef15SJeremy L Thompson     CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_OWN_POINTER, xloc);
5620c59ef15SJeremy L Thompson   }
5630c59ef15SJeremy L Thompson 
5640c59ef15SJeremy L Thompson   // Create the Q-function that builds the operator (i.e. computes its
5654b5b4ec1Sjeremylt   // quadrature data) and set its context data
5660c59ef15SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].setup,
5670c59ef15SJeremy L Thompson                               bpOptions[bpChoice].setupfname, &qf_setup);
5680c59ef15SJeremy L Thompson   CeedQFunctionAddInput(qf_setup, "x", 3, CEED_EVAL_INTERP);
5690c59ef15SJeremy L Thompson   CeedQFunctionAddInput(qf_setup, "dx", 3, CEED_EVAL_GRAD);
5700c59ef15SJeremy L Thompson   CeedQFunctionAddInput(qf_setup, "weight", 1, CEED_EVAL_WEIGHT);
5710c59ef15SJeremy L Thompson   CeedQFunctionAddOutput(qf_setup, "rho", bpOptions[bpChoice].qdatasize,
5720c59ef15SJeremy L Thompson                          CEED_EVAL_NONE);
5730c59ef15SJeremy L Thompson   CeedQFunctionAddOutput(qf_setup, "true_soln", vscale, CEED_EVAL_NONE);
5740c59ef15SJeremy L Thompson   CeedQFunctionAddOutput(qf_setup, "rhs", vscale, CEED_EVAL_INTERP);
5750c59ef15SJeremy L Thompson 
5760c59ef15SJeremy L Thompson   // Set up PDE operator
5770c59ef15SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].apply,
5780c59ef15SJeremy L Thompson                               bpOptions[bpChoice].applyfname, &qf_apply);
5790c59ef15SJeremy L Thompson   // Add inputs and outputs
5800c59ef15SJeremy L Thompson   CeedQFunctionAddInput(qf_apply, "u", vscale, bpOptions[bpChoice].inmode);
5810c59ef15SJeremy L Thompson   CeedQFunctionAddInput(qf_apply, "rho", bpOptions[bpChoice].qdatasize,
5820c59ef15SJeremy L Thompson                         CEED_EVAL_NONE);
5830c59ef15SJeremy L Thompson   CeedQFunctionAddOutput(qf_apply, "v", vscale, bpOptions[bpChoice].outmode);
5840c59ef15SJeremy L Thompson 
5850c59ef15SJeremy L Thompson   // Create the error qfunction
5860c59ef15SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].error,
5870c59ef15SJeremy L Thompson                               bpOptions[bpChoice].errorfname, &qf_error);
5880c59ef15SJeremy L Thompson   CeedQFunctionAddInput(qf_error, "u", vscale, CEED_EVAL_INTERP);
5890c59ef15SJeremy L Thompson   CeedQFunctionAddInput(qf_error, "true_soln", vscale, CEED_EVAL_NONE);
5900c59ef15SJeremy L Thompson   CeedQFunctionAddOutput(qf_error, "error", vscale, CEED_EVAL_NONE);
5910c59ef15SJeremy L Thompson 
5920c59ef15SJeremy L Thompson   // Create the persistent vectors that will be needed in setup
5930c59ef15SJeremy L Thompson   CeedInt Nqpts, Nelem = melem[0]*melem[1]*melem[2];
5940c59ef15SJeremy L Thompson   CeedBasisGetNumQuadraturePoints(basisu, &Nqpts);
5950c59ef15SJeremy L Thompson   CeedVectorCreate(ceed, bpOptions[bpChoice].qdatasize*Nelem*Nqpts, &rho);
5960c59ef15SJeremy L Thompson   CeedVectorCreate(ceed, Nelem*Nqpts*vscale, &target);
5970c59ef15SJeremy L Thompson   CeedVectorCreate(ceed, lsize*vscale, &rhsceed);
5980c59ef15SJeremy L Thompson 
5994b5b4ec1Sjeremylt   // Create the operator that builds the quadrature data for the ceed operator
6000c59ef15SJeremy L Thompson   CeedOperatorCreate(ceed, qf_setup, NULL, NULL, &op_setup);
6010c59ef15SJeremy L Thompson   CeedOperatorSetField(op_setup, "x", Erestrictx, CEED_NOTRANSPOSE,
6020c59ef15SJeremy L Thompson                        basisx, CEED_VECTOR_ACTIVE);
6030c59ef15SJeremy L Thompson   CeedOperatorSetField(op_setup, "dx", Erestrictx, CEED_NOTRANSPOSE,
6040c59ef15SJeremy L Thompson                        basisx, CEED_VECTOR_ACTIVE);
6050c59ef15SJeremy L Thompson   CeedOperatorSetField(op_setup, "weight", Erestrictxi, CEED_NOTRANSPOSE,
6060c59ef15SJeremy L Thompson                        basisx, CEED_VECTOR_NONE);
6070c59ef15SJeremy L Thompson   CeedOperatorSetField(op_setup, "rho", Erestrictqdi, CEED_NOTRANSPOSE,
6080c59ef15SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
6090c59ef15SJeremy L Thompson   CeedOperatorSetField(op_setup, "true_soln", Erestrictui, CEED_NOTRANSPOSE,
6100c59ef15SJeremy L Thompson                        CEED_BASIS_COLLOCATED, target);
6110c59ef15SJeremy L Thompson   CeedOperatorSetField(op_setup, "rhs", Erestrictu, CEED_TRANSPOSE,
6120c59ef15SJeremy L Thompson                        basisu, rhsceed);
6130c59ef15SJeremy L Thompson 
6144b5b4ec1Sjeremylt   // Create the mass or diff operator
6150c59ef15SJeremy L Thompson   CeedOperatorCreate(ceed, qf_apply, NULL, NULL, &op_apply);
6160c59ef15SJeremy L Thompson   CeedOperatorSetField(op_apply, "u", Erestrictu, CEED_TRANSPOSE,
6170c59ef15SJeremy L Thompson                        basisu, CEED_VECTOR_ACTIVE);
6180c59ef15SJeremy L Thompson   CeedOperatorSetField(op_apply, "rho", Erestrictqdi, CEED_NOTRANSPOSE,
6190c59ef15SJeremy L Thompson                        CEED_BASIS_COLLOCATED, rho);
6200c59ef15SJeremy L Thompson   CeedOperatorSetField(op_apply, "v", Erestrictu, CEED_TRANSPOSE,
6210c59ef15SJeremy L Thompson                        basisu, CEED_VECTOR_ACTIVE);
6220c59ef15SJeremy L Thompson 
6230c59ef15SJeremy L Thompson   // Create the error operator
6240c59ef15SJeremy L Thompson   CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error);
6250c59ef15SJeremy L Thompson   CeedOperatorSetField(op_error, "u", Erestrictu, CEED_TRANSPOSE,
6260c59ef15SJeremy L Thompson                        basisu, CEED_VECTOR_ACTIVE);
6270c59ef15SJeremy L Thompson   CeedOperatorSetField(op_error, "true_soln", Erestrictui, CEED_NOTRANSPOSE,
6280c59ef15SJeremy L Thompson                        CEED_BASIS_COLLOCATED, target);
6290c59ef15SJeremy L Thompson   CeedOperatorSetField(op_error, "error", Erestrictui, CEED_NOTRANSPOSE,
6300c59ef15SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
6310c59ef15SJeremy L Thompson 
6320c59ef15SJeremy L Thompson 
6330c59ef15SJeremy L Thompson   // Set up Mat
6340c59ef15SJeremy L Thompson   ierr = PetscMalloc1(1, &user); CHKERRQ(ierr);
6350c59ef15SJeremy L Thompson   user->comm = comm;
6360c59ef15SJeremy L Thompson   user->ltog = ltog;
6370c59ef15SJeremy L Thompson   if (bpChoice != CEED_BP1 && bpChoice != CEED_BP2) {
6380c59ef15SJeremy L Thompson     user->ltog0 = ltog0;
6390c59ef15SJeremy L Thompson     user->gtogD = gtogD;
6400c59ef15SJeremy L Thompson   }
6410c59ef15SJeremy L Thompson   user->Xloc = Xloc;
6420c59ef15SJeremy L Thompson   ierr = VecDuplicate(Xloc, &user->Yloc); CHKERRQ(ierr);
6430c59ef15SJeremy L Thompson   CeedVectorCreate(ceed, lsize*vscale, &user->xceed);
6440c59ef15SJeremy L Thompson   CeedVectorCreate(ceed, lsize*vscale, &user->yceed);
6450c59ef15SJeremy L Thompson   user->op = op_apply;
6460c59ef15SJeremy L Thompson   user->rho = rho;
6470c59ef15SJeremy L Thompson   user->ceed = ceed;
6480c59ef15SJeremy L Thompson 
6490c59ef15SJeremy L Thompson   ierr = MatCreateShell(comm, mdof[0]*mdof[1]*mdof[2]*vscale,
6500c59ef15SJeremy L Thompson                         mdof[0]*mdof[1]*mdof[2]*vscale,
6510c59ef15SJeremy L Thompson                         PETSC_DECIDE, PETSC_DECIDE, user, &mat); CHKERRQ(ierr);
6520c59ef15SJeremy L Thompson   if (bpChoice == CEED_BP1 || bpChoice == CEED_BP2) {
6530c59ef15SJeremy L Thompson     ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Mass);
6540c59ef15SJeremy L Thompson     CHKERRQ(ierr);
6550c59ef15SJeremy L Thompson   } else {
6560c59ef15SJeremy L Thompson     ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Diff);
6570c59ef15SJeremy L Thompson     CHKERRQ(ierr);
6580c59ef15SJeremy L Thompson   }
6590c59ef15SJeremy L Thompson   ierr = MatCreateVecs(mat, &rhs, NULL); CHKERRQ(ierr);
6600c59ef15SJeremy L Thompson 
6610c59ef15SJeremy L Thompson   // Get RHS vector
6620c59ef15SJeremy L Thompson   ierr = VecDuplicate(Xloc, &rhsloc); CHKERRQ(ierr);
6630c59ef15SJeremy L Thompson   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
6640c59ef15SJeremy L Thompson   ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr);
6650c59ef15SJeremy L Thompson   CeedVectorSetArray(rhsceed, CEED_MEM_HOST, CEED_USE_POINTER, r);
6660c59ef15SJeremy L Thompson 
6670c59ef15SJeremy L Thompson   // Setup rho, rhs, and target
6680c59ef15SJeremy L Thompson   CeedOperatorApply(op_setup, xcoord, rho, CEED_REQUEST_IMMEDIATE);
6690c59ef15SJeremy L Thompson   ierr = CeedVectorSyncArray(rhsceed, CEED_MEM_HOST); CHKERRQ(ierr);
6700c59ef15SJeremy L Thompson   CeedVectorDestroy(&xcoord);
6710c59ef15SJeremy L Thompson 
6720c59ef15SJeremy L Thompson   // Gather RHS
6730c59ef15SJeremy L Thompson   ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr);
6740c59ef15SJeremy L Thompson   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
6750c59ef15SJeremy L Thompson   ierr = VecScatterBegin(ltog, rhsloc, rhs, ADD_VALUES, SCATTER_FORWARD);
6760c59ef15SJeremy L Thompson   CHKERRQ(ierr);
6770c59ef15SJeremy L Thompson   ierr = VecScatterEnd(ltog, rhsloc, rhs, ADD_VALUES, SCATTER_FORWARD);
6780c59ef15SJeremy L Thompson   CHKERRQ(ierr);
6790c59ef15SJeremy L Thompson   CeedVectorDestroy(&rhsceed);
6800c59ef15SJeremy L Thompson 
6810c59ef15SJeremy L Thompson   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
6820c59ef15SJeremy L Thompson   {
6830c59ef15SJeremy L Thompson     PC pc;
6840c59ef15SJeremy L Thompson     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
6850c59ef15SJeremy L Thompson     if (bpChoice == CEED_BP1 || bpChoice == CEED_BP2) {
6860c59ef15SJeremy L Thompson       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
6870c59ef15SJeremy L Thompson       ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr);
6880c59ef15SJeremy L Thompson     } else {
6890c59ef15SJeremy L Thompson       ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr);
6900c59ef15SJeremy L Thompson     }
6910c59ef15SJeremy L Thompson     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
6920c59ef15SJeremy L Thompson     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
6930c59ef15SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
6940c59ef15SJeremy L Thompson                             PETSC_DEFAULT); CHKERRQ(ierr);
6950c59ef15SJeremy L Thompson   }
6960c59ef15SJeremy L Thompson   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
6970c59ef15SJeremy L Thompson   ierr = KSPSetOperators(ksp, mat, mat); CHKERRQ(ierr);
6980c59ef15SJeremy L Thompson   // First run, if benchmarking
6990c59ef15SJeremy L Thompson   if (benchmark_mode) {
7000c59ef15SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
7010c59ef15SJeremy L Thompson     CHKERRQ(ierr);
7020c59ef15SJeremy L Thompson     my_rt_start = MPI_Wtime();
7030c59ef15SJeremy L Thompson     ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
7040c59ef15SJeremy L Thompson     my_rt = MPI_Wtime() - my_rt_start;
7050c59ef15SJeremy L Thompson     // Set maxits based on first iteration timing
7060c59ef15SJeremy L Thompson     if (my_rt > 0.02) {
7070c59ef15SJeremy L Thompson       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
7080c59ef15SJeremy L Thompson       CHKERRQ(ierr);
7090c59ef15SJeremy L Thompson     } else {
7100c59ef15SJeremy L Thompson       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
7110c59ef15SJeremy L Thompson       CHKERRQ(ierr);
7120c59ef15SJeremy L Thompson     }
7130c59ef15SJeremy L Thompson   }
7140c59ef15SJeremy L Thompson   // Timed solve
7150c59ef15SJeremy L Thompson   my_rt_start = MPI_Wtime();
7160c59ef15SJeremy L Thompson   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
7170c59ef15SJeremy L Thompson   my_rt = MPI_Wtime() - my_rt_start;
7180c59ef15SJeremy L Thompson   {
7190c59ef15SJeremy L Thompson     KSPType ksptype;
7200c59ef15SJeremy L Thompson     KSPConvergedReason reason;
7210c59ef15SJeremy L Thompson     PetscReal rnorm;
7220c59ef15SJeremy L Thompson     PetscInt its;
7230c59ef15SJeremy L Thompson     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
7240c59ef15SJeremy L Thompson     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
7250c59ef15SJeremy L Thompson     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
7260c59ef15SJeremy L Thompson     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
7270c59ef15SJeremy L Thompson     if (!test_mode || reason < 0 || rnorm > 1e-8) {
7283ce86546SJeremy L Thompson       ierr = PetscPrintf(comm,
7293ce86546SJeremy L Thompson                          "  KSP:\n"
7303ce86546SJeremy L Thompson                          "    KSP Type                           : %s\n"
7313ce86546SJeremy L Thompson                          "    KSP Convergence                    : %s\n"
7323ce86546SJeremy L Thompson                          "    Total KSP Iterations               : %D\n"
7333ce86546SJeremy L Thompson                          "    Final rnorm                        : %e\n",
7343ce86546SJeremy L Thompson                          ksptype, KSPConvergedReasons[reason], its,
7353ce86546SJeremy L Thompson                          (double)rnorm); CHKERRQ(ierr);
7360c59ef15SJeremy L Thompson     }
7370c59ef15SJeremy L Thompson     if (benchmark_mode && (!test_mode)) {
7380c59ef15SJeremy L Thompson       CeedInt gsize;
7390c59ef15SJeremy L Thompson       ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
7400c59ef15SJeremy L Thompson       MPI_Reduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, 0, comm);
7410c59ef15SJeremy L Thompson       MPI_Reduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, 0, comm);
7420c59ef15SJeremy L Thompson       ierr = PetscPrintf(comm,
7433ce86546SJeremy L Thompson                          "  Performance:\n"
7443ce86546SJeremy L Thompson                          "    CG Solve Time                      : %g (%g) sec\n"
7453ce86546SJeremy L Thompson                          "    DOFs/Sec in CG                     : %g (%g) million\n",
7463ce86546SJeremy L Thompson                          rt_max, rt_min, 1e-6*gsize*its/rt_max,
7473ce86546SJeremy L Thompson                          1e-6*gsize*its/rt_min); CHKERRQ(ierr);
7480c59ef15SJeremy L Thompson     }
7490c59ef15SJeremy L Thompson   }
7500c59ef15SJeremy L Thompson 
7510c59ef15SJeremy L Thompson   {
7520c59ef15SJeremy L Thompson     PetscReal maxerror;
7530c59ef15SJeremy L Thompson     ierr = ComputeErrorMax(user, op_error, X, target, &maxerror); CHKERRQ(ierr);
7540c59ef15SJeremy L Thompson     PetscReal tol = (bpChoice == CEED_BP1 || bpChoice == CEED_BP2) ? 5e-3 : 5e-2;
7550c59ef15SJeremy L Thompson     if (!test_mode || maxerror > tol) {
7563ce86546SJeremy L Thompson       ierr = PetscPrintf(comm,
7573ce86546SJeremy L Thompson                          "    Pointwise Error (max)              : %e\n",
7583ce86546SJeremy L Thompson                          (double)maxerror); CHKERRQ(ierr);
7590c59ef15SJeremy L Thompson     }
7600c59ef15SJeremy L Thompson   }
7610c59ef15SJeremy L Thompson 
7620c59ef15SJeremy L Thompson   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
7630c59ef15SJeremy L Thompson   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
7640c59ef15SJeremy L Thompson   ierr = VecDestroy(&X); CHKERRQ(ierr);
7650c59ef15SJeremy L Thompson   ierr = VecDestroy(&user->Xloc); CHKERRQ(ierr);
7660c59ef15SJeremy L Thompson   ierr = VecDestroy(&user->Yloc); CHKERRQ(ierr);
7670c59ef15SJeremy L Thompson   ierr = VecScatterDestroy(&ltog); CHKERRQ(ierr);
7680c59ef15SJeremy L Thompson   ierr = VecScatterDestroy(&ltog0); CHKERRQ(ierr);
7690c59ef15SJeremy L Thompson   ierr = VecScatterDestroy(&gtogD); CHKERRQ(ierr);
7700c59ef15SJeremy L Thompson   ierr = MatDestroy(&mat); CHKERRQ(ierr);
7710c59ef15SJeremy L Thompson   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
7720c59ef15SJeremy L Thompson 
7730c59ef15SJeremy L Thompson   CeedVectorDestroy(&user->xceed);
7740c59ef15SJeremy L Thompson   CeedVectorDestroy(&user->yceed);
7750c59ef15SJeremy L Thompson   CeedVectorDestroy(&user->rho);
7760c59ef15SJeremy L Thompson   CeedVectorDestroy(&target);
7770c59ef15SJeremy L Thompson   CeedOperatorDestroy(&op_setup);
7780c59ef15SJeremy L Thompson   CeedOperatorDestroy(&op_apply);
7790c59ef15SJeremy L Thompson   CeedOperatorDestroy(&op_error);
7800c59ef15SJeremy L Thompson   CeedElemRestrictionDestroy(&Erestrictu);
7810c59ef15SJeremy L Thompson   CeedElemRestrictionDestroy(&Erestrictx);
7820c59ef15SJeremy L Thompson   CeedElemRestrictionDestroy(&Erestrictui);
7830c59ef15SJeremy L Thompson   CeedElemRestrictionDestroy(&Erestrictxi);
7840c59ef15SJeremy L Thompson   CeedElemRestrictionDestroy(&Erestrictqdi);
7850c59ef15SJeremy L Thompson   CeedQFunctionDestroy(&qf_setup);
7860c59ef15SJeremy L Thompson   CeedQFunctionDestroy(&qf_apply);
7870c59ef15SJeremy L Thompson   CeedQFunctionDestroy(&qf_error);
7880c59ef15SJeremy L Thompson   CeedBasisDestroy(&basisu);
7890c59ef15SJeremy L Thompson   CeedBasisDestroy(&basisx);
7900c59ef15SJeremy L Thompson   CeedDestroy(&ceed);
7910c59ef15SJeremy L Thompson   ierr = PetscFree(user); CHKERRQ(ierr);
7920c59ef15SJeremy L Thompson   return PetscFinalize();
7930c59ef15SJeremy L Thompson }
794