xref: /petsc/src/dm/dt/fe/tests/ex1.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static const char help[] = "Performance Tests for FE Integration";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown #include <petscdmplex.h>
4*c4762a1bSJed Brown #include <petscfe.h>
5*c4762a1bSJed Brown #include <petscds.h>
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown typedef struct {
8*c4762a1bSJed Brown   PetscInt  dim;     /* The topological dimension */
9*c4762a1bSJed Brown   PetscBool simplex; /* True for simplices, false for hexes */
10*c4762a1bSJed Brown   PetscInt  its;     /* Number of replications for timing */
11*c4762a1bSJed Brown   PetscInt  cbs;     /* Number of cells in an integration block */
12*c4762a1bSJed Brown } AppCtx;
13*c4762a1bSJed Brown 
14*c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
15*c4762a1bSJed Brown {
16*c4762a1bSJed Brown   PetscErrorCode ierr;
17*c4762a1bSJed Brown 
18*c4762a1bSJed Brown   PetscFunctionBeginUser;
19*c4762a1bSJed Brown   options->dim     = 2;
20*c4762a1bSJed Brown   options->simplex = PETSC_TRUE;
21*c4762a1bSJed Brown   options->its     = 1;
22*c4762a1bSJed Brown   options->cbs     = 8;
23*c4762a1bSJed Brown 
24*c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "FE Integration Performance Options", "PETSCFE");CHKERRQ(ierr);
25*c4762a1bSJed Brown   ierr = PetscOptionsInt("-dim", "The topological dimension", "ex1.c", options->dim, &options->dim, NULL);CHKERRQ(ierr);
26*c4762a1bSJed Brown   ierr = PetscOptionsBool("-simplex", "Simplex or hex cells", "ex1.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
27*c4762a1bSJed Brown   ierr = PetscOptionsInt("-its", "The number of replications for timing", "ex1.c", options->its, &options->its, NULL);CHKERRQ(ierr);
28*c4762a1bSJed Brown   ierr = PetscOptionsInt("-cbs", "The number of cells in an integration block", "ex1.c", options->cbs, &options->cbs, NULL);CHKERRQ(ierr);
29*c4762a1bSJed Brown   ierr = PetscOptionsEnd();
30*c4762a1bSJed Brown   PetscFunctionReturn(0);
31*c4762a1bSJed Brown }
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
34*c4762a1bSJed Brown {
35*c4762a1bSJed Brown   PetscInt d;
36*c4762a1bSJed Brown   *u = 0.0;
37*c4762a1bSJed Brown   for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0*PETSC_PI*x[d]);
38*c4762a1bSJed Brown   return 0;
39*c4762a1bSJed Brown }
40*c4762a1bSJed Brown 
41*c4762a1bSJed Brown static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
42*c4762a1bSJed Brown                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
43*c4762a1bSJed Brown                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
44*c4762a1bSJed Brown                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
45*c4762a1bSJed Brown {
46*c4762a1bSJed Brown   PetscInt d;
47*c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f0[0] += -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]);
48*c4762a1bSJed Brown }
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
51*c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
52*c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
53*c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
54*c4762a1bSJed Brown {
55*c4762a1bSJed Brown   PetscInt d;
56*c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
57*c4762a1bSJed Brown }
58*c4762a1bSJed Brown 
59*c4762a1bSJed Brown static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
60*c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
61*c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
62*c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
63*c4762a1bSJed Brown {
64*c4762a1bSJed Brown   PetscInt d;
65*c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
66*c4762a1bSJed Brown }
67*c4762a1bSJed Brown 
68*c4762a1bSJed Brown static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
69*c4762a1bSJed Brown {
70*c4762a1bSJed Brown   PetscDS        prob;
71*c4762a1bSJed Brown   const PetscInt id = 1;
72*c4762a1bSJed Brown   PetscErrorCode ierr;
73*c4762a1bSJed Brown 
74*c4762a1bSJed Brown   PetscFunctionBeginUser;
75*c4762a1bSJed Brown   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
76*c4762a1bSJed Brown   ierr = PetscDSSetResidual(prob, 0, f0_trig_u, f1_u);CHKERRQ(ierr);
77*c4762a1bSJed Brown   ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr);
78*c4762a1bSJed Brown   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) trig_u, 1, &id, user);CHKERRQ(ierr);
79*c4762a1bSJed Brown   ierr = PetscDSSetExactSolution(prob, 0, trig_u, user);CHKERRQ(ierr);
80*c4762a1bSJed Brown   PetscFunctionReturn(0);
81*c4762a1bSJed Brown }
82*c4762a1bSJed Brown 
83*c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
84*c4762a1bSJed Brown {
85*c4762a1bSJed Brown   DM             cdm = dm;
86*c4762a1bSJed Brown   PetscFE        fe;
87*c4762a1bSJed Brown   char           prefix[PETSC_MAX_PATH_LEN];
88*c4762a1bSJed Brown   PetscErrorCode ierr;
89*c4762a1bSJed Brown 
90*c4762a1bSJed Brown   PetscFunctionBeginUser;
91*c4762a1bSJed Brown   /* Create finite element */
92*c4762a1bSJed Brown   ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name);CHKERRQ(ierr);
93*c4762a1bSJed Brown   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), user->dim, 1, user->simplex, name ? prefix : NULL, -1, &fe);CHKERRQ(ierr);
94*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr);
95*c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
96*c4762a1bSJed Brown   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
97*c4762a1bSJed Brown   ierr = DMCreateDS(dm);CHKERRQ(ierr);
98*c4762a1bSJed Brown   ierr = (*setup)(dm, user);CHKERRQ(ierr);
99*c4762a1bSJed Brown   while (cdm) {
100*c4762a1bSJed Brown     ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr);
101*c4762a1bSJed Brown     /* TODO: Check whether the boundary of coarse meshes is marked */
102*c4762a1bSJed Brown     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
103*c4762a1bSJed Brown   }
104*c4762a1bSJed Brown   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
105*c4762a1bSJed Brown   PetscFunctionReturn(0);
106*c4762a1bSJed Brown }
107*c4762a1bSJed Brown 
108*c4762a1bSJed Brown static PetscErrorCode PetscContainerUserDestroy_PetscFEGeom(void *ctx)
109*c4762a1bSJed Brown {
110*c4762a1bSJed Brown   PetscFEGeom   *geom = (PetscFEGeom *) ctx;
111*c4762a1bSJed Brown   PetscErrorCode ierr;
112*c4762a1bSJed Brown 
113*c4762a1bSJed Brown   PetscFunctionBegin;
114*c4762a1bSJed Brown   ierr = PetscFEGeomDestroy(&geom);CHKERRQ(ierr);
115*c4762a1bSJed Brown   PetscFunctionReturn(0);
116*c4762a1bSJed Brown }
117*c4762a1bSJed Brown 
118*c4762a1bSJed Brown PetscErrorCode CellRangeGetFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
119*c4762a1bSJed Brown {
120*c4762a1bSJed Brown   char            composeStr[33] = {0};
121*c4762a1bSJed Brown   PetscObjectId   id;
122*c4762a1bSJed Brown   PetscContainer  container;
123*c4762a1bSJed Brown   PetscErrorCode  ierr;
124*c4762a1bSJed Brown 
125*c4762a1bSJed Brown   PetscFunctionBegin;
126*c4762a1bSJed Brown   ierr = PetscObjectGetId((PetscObject) quad, &id);CHKERRQ(ierr);
127*c4762a1bSJed Brown   ierr = PetscSNPrintf(composeStr, 32, "CellRangeGetFEGeom_%x\n", id);CHKERRQ(ierr);
128*c4762a1bSJed Brown   ierr = PetscObjectQuery((PetscObject) cellIS, composeStr, (PetscObject *) &container);CHKERRQ(ierr);
129*c4762a1bSJed Brown   if (container) {
130*c4762a1bSJed Brown     ierr = PetscContainerGetPointer(container, (void **) geom);CHKERRQ(ierr);
131*c4762a1bSJed Brown   } else {
132*c4762a1bSJed Brown     ierr = DMFieldCreateFEGeom(coordField, cellIS, quad, faceData, geom);CHKERRQ(ierr);
133*c4762a1bSJed Brown     ierr = PetscContainerCreate(PETSC_COMM_SELF, &container);CHKERRQ(ierr);
134*c4762a1bSJed Brown     ierr = PetscContainerSetPointer(container, (void *) *geom);CHKERRQ(ierr);
135*c4762a1bSJed Brown     ierr = PetscContainerSetUserDestroy(container, PetscContainerUserDestroy_PetscFEGeom);CHKERRQ(ierr);
136*c4762a1bSJed Brown     ierr = PetscObjectCompose((PetscObject) cellIS, composeStr, (PetscObject) container);CHKERRQ(ierr);
137*c4762a1bSJed Brown     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
138*c4762a1bSJed Brown   }
139*c4762a1bSJed Brown   PetscFunctionReturn(0);
140*c4762a1bSJed Brown }
141*c4762a1bSJed Brown 
142*c4762a1bSJed Brown PetscErrorCode CellRangeRestoreFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
143*c4762a1bSJed Brown {
144*c4762a1bSJed Brown   PetscFunctionBegin;
145*c4762a1bSJed Brown   *geom = NULL;
146*c4762a1bSJed Brown   PetscFunctionReturn(0);
147*c4762a1bSJed Brown }
148*c4762a1bSJed Brown 
149*c4762a1bSJed Brown static PetscErrorCode CreateFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms)
150*c4762a1bSJed Brown {
151*c4762a1bSJed Brown   DMField        coordField;
152*c4762a1bSJed Brown   PetscInt       Nf, f, maxDegree;
153*c4762a1bSJed Brown   PetscErrorCode ierr;
154*c4762a1bSJed Brown 
155*c4762a1bSJed Brown   PetscFunctionBeginUser;
156*c4762a1bSJed Brown   *affineQuad = NULL;
157*c4762a1bSJed Brown   *affineGeom = NULL;
158*c4762a1bSJed Brown   *quads      = NULL;
159*c4762a1bSJed Brown   *geoms      = NULL;
160*c4762a1bSJed Brown   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
161*c4762a1bSJed Brown   ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr);
162*c4762a1bSJed Brown   ierr = DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree);CHKERRQ(ierr);
163*c4762a1bSJed Brown   if (maxDegree <= 1) {
164*c4762a1bSJed Brown     ierr = DMFieldCreateDefaultQuadrature(coordField, cellIS, affineQuad);CHKERRQ(ierr);
165*c4762a1bSJed Brown     if (*affineQuad) {ierr = CellRangeGetFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom);CHKERRQ(ierr);}
166*c4762a1bSJed Brown   } else {
167*c4762a1bSJed Brown     ierr = PetscCalloc2(Nf, quads, Nf, geoms);CHKERRQ(ierr);
168*c4762a1bSJed Brown     for (f = 0; f < Nf; ++f) {
169*c4762a1bSJed Brown       PetscFE fe;
170*c4762a1bSJed Brown 
171*c4762a1bSJed Brown       ierr = PetscDSGetDiscretization(ds, f, (PetscObject *) &fe);CHKERRQ(ierr);
172*c4762a1bSJed Brown       ierr = PetscFEGetQuadrature(fe, &(*quads)[f]);CHKERRQ(ierr);
173*c4762a1bSJed Brown       ierr = PetscObjectReference((PetscObject) (*quads)[f]);CHKERRQ(ierr);
174*c4762a1bSJed Brown       ierr = CellRangeGetFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f]);CHKERRQ(ierr);
175*c4762a1bSJed Brown     }
176*c4762a1bSJed Brown   }
177*c4762a1bSJed Brown   PetscFunctionReturn(0);
178*c4762a1bSJed Brown }
179*c4762a1bSJed Brown 
180*c4762a1bSJed Brown static PetscErrorCode DestroyFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms)
181*c4762a1bSJed Brown {
182*c4762a1bSJed Brown   DMField        coordField;
183*c4762a1bSJed Brown   PetscInt       Nf, f;
184*c4762a1bSJed Brown   PetscErrorCode ierr;
185*c4762a1bSJed Brown 
186*c4762a1bSJed Brown   PetscFunctionBeginUser;
187*c4762a1bSJed Brown   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
188*c4762a1bSJed Brown   ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr);
189*c4762a1bSJed Brown   if (*affineQuad) {
190*c4762a1bSJed Brown     ierr = CellRangeRestoreFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom);CHKERRQ(ierr);
191*c4762a1bSJed Brown     ierr = PetscQuadratureDestroy(affineQuad);CHKERRQ(ierr);
192*c4762a1bSJed Brown   } else {
193*c4762a1bSJed Brown     for (f = 0; f < Nf; ++f) {
194*c4762a1bSJed Brown       ierr = CellRangeRestoreFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f]);CHKERRQ(ierr);
195*c4762a1bSJed Brown       ierr = PetscQuadratureDestroy(&(*quads)[f]);CHKERRQ(ierr);
196*c4762a1bSJed Brown     }
197*c4762a1bSJed Brown     ierr = PetscFree2(*quads, *geoms);CHKERRQ(ierr);
198*c4762a1bSJed Brown   }
199*c4762a1bSJed Brown   PetscFunctionReturn(0);
200*c4762a1bSJed Brown }
201*c4762a1bSJed Brown 
202*c4762a1bSJed Brown static PetscErrorCode TestIntegration(DM dm, PetscInt cbs, PetscInt its)
203*c4762a1bSJed Brown {
204*c4762a1bSJed Brown   PetscDS         ds;
205*c4762a1bSJed Brown   PetscFEGeom    *chunkGeom = NULL;
206*c4762a1bSJed Brown   PetscQuadrature affineQuad,  *quads = NULL;
207*c4762a1bSJed Brown   PetscFEGeom    *affineGeom, **geoms = NULL;
208*c4762a1bSJed Brown   PetscScalar    *u, *elemVec;
209*c4762a1bSJed Brown   IS              cellIS;
210*c4762a1bSJed Brown   PetscInt        depth, cStart, cEnd, cell, chunkSize = cbs, Nch = 0, Nf, f, totDim, i, k;
211*c4762a1bSJed Brown   PetscLogStage   stage;
212*c4762a1bSJed Brown   PetscLogEvent   event;
213*c4762a1bSJed Brown   PetscErrorCode  ierr;
214*c4762a1bSJed Brown 
215*c4762a1bSJed Brown   PetscFunctionBeginUser;
216*c4762a1bSJed Brown   ierr = PetscLogStageRegister("PetscFE Residual Integration Test", &stage);CHKERRQ(ierr);
217*c4762a1bSJed Brown   ierr = PetscLogEventRegister("FEIntegRes", PETSCFE_CLASSID, &event);CHKERRQ(ierr);
218*c4762a1bSJed Brown   ierr = PetscLogStagePush(stage);CHKERRQ(ierr);
219*c4762a1bSJed Brown   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
220*c4762a1bSJed Brown   ierr = DMGetStratumIS(dm, "depth", depth, &cellIS);CHKERRQ(ierr);
221*c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
222*c4762a1bSJed Brown   ierr = DMGetCellDS(dm, cStart, &ds);CHKERRQ(ierr);
223*c4762a1bSJed Brown   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
224*c4762a1bSJed Brown   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
225*c4762a1bSJed Brown   ierr = CreateFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms);CHKERRQ(ierr);;
226*c4762a1bSJed Brown   ierr = PetscMalloc2(chunkSize*totDim, &u, chunkSize*totDim, &elemVec);CHKERRQ(ierr);
227*c4762a1bSJed Brown   /* Assumptions:
228*c4762a1bSJed Brown     - Single field
229*c4762a1bSJed Brown     - No input data
230*c4762a1bSJed Brown     - No auxiliary data
231*c4762a1bSJed Brown     - No time-dependence
232*c4762a1bSJed Brown   */
233*c4762a1bSJed Brown   for (i = 0; i < its; ++i) {
234*c4762a1bSJed Brown     for (cell = cStart; cell < cEnd; cell += chunkSize, ++Nch) {
235*c4762a1bSJed Brown       const PetscInt cS = cell, cE = PetscMin(cS + chunkSize, cEnd), Ne = cE - cS;
236*c4762a1bSJed Brown 
237*c4762a1bSJed Brown       ierr = PetscArrayzero(elemVec, chunkSize*totDim);CHKERRQ(ierr);
238*c4762a1bSJed Brown       /* TODO Replace with DMPlexGetCellFields() */
239*c4762a1bSJed Brown       for (k = 0; k < chunkSize*totDim; ++k) u[k] = 1.0;
240*c4762a1bSJed Brown       for (f = 0; f < Nf; ++f) {
241*c4762a1bSJed Brown         PetscFEGeom    *geom = affineGeom ? affineGeom : geoms[f];
242*c4762a1bSJed Brown         /* PetscQuadrature quad = affineQuad ? affineQuad : quads[f]; */
243*c4762a1bSJed Brown 
244*c4762a1bSJed Brown         ierr = PetscFEGeomGetChunk(geom, cS, cE, &chunkGeom);CHKERRQ(ierr);
245*c4762a1bSJed Brown         ierr = PetscLogEventBegin(event,0,0,0,0);CHKERRQ(ierr);
246*c4762a1bSJed Brown         ierr = PetscFEIntegrateResidual(ds, f, Ne, chunkGeom, u, NULL, NULL, NULL, 0.0, elemVec);CHKERRQ(ierr);
247*c4762a1bSJed Brown         ierr = PetscLogEventEnd(event,0,0,0,0);CHKERRQ(ierr);
248*c4762a1bSJed Brown       }
249*c4762a1bSJed Brown     }
250*c4762a1bSJed Brown   }
251*c4762a1bSJed Brown   ierr = PetscFEGeomRestoreChunk(affineGeom, cStart, cEnd, &chunkGeom);CHKERRQ(ierr);
252*c4762a1bSJed Brown   ierr = DestroyFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms);CHKERRQ(ierr);;
253*c4762a1bSJed Brown   ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
254*c4762a1bSJed Brown   ierr = PetscFree2(u, elemVec);CHKERRQ(ierr);
255*c4762a1bSJed Brown   ierr = PetscLogStagePop();CHKERRQ(ierr);
256*c4762a1bSJed Brown   {
257*c4762a1bSJed Brown     const char        *title = "Petsc FE Residual Integration";
258*c4762a1bSJed Brown     PetscEventPerfInfo eventInfo;
259*c4762a1bSJed Brown     PetscInt           N = (cEnd - cStart)*Nf*its;
260*c4762a1bSJed Brown     PetscReal          flopRate, cellRate;
261*c4762a1bSJed Brown 
262*c4762a1bSJed Brown     ierr = PetscLogEventGetPerfInfo(stage, event, &eventInfo);CHKERRQ(ierr);
263*c4762a1bSJed Brown     flopRate = eventInfo.time != 0.0 ? eventInfo.flops/eventInfo.time : 0.0;
264*c4762a1bSJed Brown     cellRate = eventInfo.time != 0.0 ? N/eventInfo.time : 0.0;
265*c4762a1bSJed Brown     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "%s: %D integrals %D chunks %D reps\n  Cell rate: %.2f/s flop rate: %.2f MF/s\n", title, N, Nch, its, (double)cellRate, (double)(flopRate/1.e6));CHKERRQ(ierr);
266*c4762a1bSJed Brown   }
267*c4762a1bSJed Brown   PetscFunctionReturn(0);
268*c4762a1bSJed Brown }
269*c4762a1bSJed Brown 
270*c4762a1bSJed Brown static PetscErrorCode TestIntegration2(DM dm, PetscInt cbs, PetscInt its)
271*c4762a1bSJed Brown {
272*c4762a1bSJed Brown   Vec             X, F;
273*c4762a1bSJed Brown   PetscLogStage   stage;
274*c4762a1bSJed Brown   PetscLogEvent   event;
275*c4762a1bSJed Brown   PetscInt        i;
276*c4762a1bSJed Brown   PetscErrorCode  ierr;
277*c4762a1bSJed Brown 
278*c4762a1bSJed Brown   PetscFunctionBeginUser;
279*c4762a1bSJed Brown   ierr = PetscLogStageRegister("DMPlex Residual Integration Test", &stage);CHKERRQ(ierr);
280*c4762a1bSJed Brown   ierr = PetscLogEventGetId("DMPlexResidualFE", &event);CHKERRQ(ierr);
281*c4762a1bSJed Brown   ierr = PetscLogStagePush(stage);CHKERRQ(ierr);
282*c4762a1bSJed Brown   ierr = DMGetLocalVector(dm, &X);CHKERRQ(ierr);
283*c4762a1bSJed Brown   ierr = DMGetLocalVector(dm, &F);CHKERRQ(ierr);
284*c4762a1bSJed Brown   for (i = 0; i < its; ++i) {
285*c4762a1bSJed Brown     ierr = DMPlexSNESComputeResidualFEM(dm, X, F, NULL);CHKERRQ(ierr);
286*c4762a1bSJed Brown   }
287*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(dm, &X);CHKERRQ(ierr);
288*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(dm, &F);CHKERRQ(ierr);
289*c4762a1bSJed Brown   ierr = PetscLogStagePop();CHKERRQ(ierr);
290*c4762a1bSJed Brown   {
291*c4762a1bSJed Brown     const char        *title = "DMPlex Residual Integration";
292*c4762a1bSJed Brown     PetscEventPerfInfo eventInfo;
293*c4762a1bSJed Brown     PetscReal          flopRate, cellRate;
294*c4762a1bSJed Brown     PetscInt           cStart, cEnd, Nf, N;
295*c4762a1bSJed Brown 
296*c4762a1bSJed Brown     ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
297*c4762a1bSJed Brown     ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
298*c4762a1bSJed Brown     ierr = PetscLogEventGetPerfInfo(stage, event, &eventInfo);CHKERRQ(ierr);
299*c4762a1bSJed Brown     N        = (cEnd - cStart)*Nf*eventInfo.count;
300*c4762a1bSJed Brown     flopRate = eventInfo.time != 0.0 ? eventInfo.flops/eventInfo.time : 0.0;
301*c4762a1bSJed Brown     cellRate = eventInfo.time != 0.0 ? N/eventInfo.time : 0.0;
302*c4762a1bSJed Brown     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "%s: %D integrals %D reps\n  Cell rate: %.2f/s flop rate: %.2f MF/s\n", title, N, eventInfo.count, (double)cellRate, (double)(flopRate/1.e6));CHKERRQ(ierr);
303*c4762a1bSJed Brown   }
304*c4762a1bSJed Brown   PetscFunctionReturn(0);
305*c4762a1bSJed Brown }
306*c4762a1bSJed Brown 
307*c4762a1bSJed Brown int main(int argc, char **argv)
308*c4762a1bSJed Brown {
309*c4762a1bSJed Brown   DM             dm;
310*c4762a1bSJed Brown   AppCtx         ctx;
311*c4762a1bSJed Brown   PetscMPIInt    size;
312*c4762a1bSJed Brown   PetscErrorCode ierr;
313*c4762a1bSJed Brown 
314*c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help); if (ierr) return ierr;
315*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRQ(ierr);
316*c4762a1bSJed Brown   if (size > 1) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "This is a uniprocessor example only.");
317*c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &ctx);CHKERRQ(ierr);
318*c4762a1bSJed Brown   ierr = PetscLogDefaultBegin();CHKERRQ(ierr);
319*c4762a1bSJed Brown   ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, ctx.dim, ctx.simplex, NULL, NULL, NULL, NULL, PETSC_TRUE, &dm);CHKERRQ(ierr);
320*c4762a1bSJed Brown   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
321*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) dm, "Mesh");CHKERRQ(ierr);
322*c4762a1bSJed Brown   ierr = PetscObjectViewFromOptions((PetscObject) dm, NULL, "-dm_view");CHKERRQ(ierr);
323*c4762a1bSJed Brown   ierr = SetupDiscretization(dm, "potential", SetupPrimalProblem, &ctx);CHKERRQ(ierr);
324*c4762a1bSJed Brown   ierr = TestIntegration(dm, ctx.cbs, ctx.its);CHKERRQ(ierr);
325*c4762a1bSJed Brown   ierr = TestIntegration2(dm, ctx.cbs, ctx.its);CHKERRQ(ierr);
326*c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
327*c4762a1bSJed Brown   ierr = PetscFinalize();
328*c4762a1bSJed Brown   return ierr;
329*c4762a1bSJed Brown }
330*c4762a1bSJed Brown 
331*c4762a1bSJed Brown /*TEST
332*c4762a1bSJed Brown   test:
333*c4762a1bSJed Brown     suffix: 0
334*c4762a1bSJed Brown     requires: triangle
335*c4762a1bSJed Brown     args: -dm_view
336*c4762a1bSJed Brown 
337*c4762a1bSJed Brown   test:
338*c4762a1bSJed Brown     suffix: 1
339*c4762a1bSJed Brown     requires: triangle
340*c4762a1bSJed Brown     args: -dm_view -potential_petscspace_degree 1
341*c4762a1bSJed Brown 
342*c4762a1bSJed Brown   test:
343*c4762a1bSJed Brown     suffix: 2
344*c4762a1bSJed Brown     requires: triangle
345*c4762a1bSJed Brown     args: -dm_view -potential_petscspace_degree 2
346*c4762a1bSJed Brown 
347*c4762a1bSJed Brown   test:
348*c4762a1bSJed Brown     suffix: 3
349*c4762a1bSJed Brown     requires: triangle
350*c4762a1bSJed Brown     args: -dm_view -potential_petscspace_degree 3
351*c4762a1bSJed Brown TEST*/
352