xref: /petsc/src/dm/dt/fe/tests/ex1.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown static const char help[] = "Performance Tests for FE Integration";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown #include <petscfe.h>
5c4762a1bSJed Brown #include <petscds.h>
6c4762a1bSJed Brown 
7c4762a1bSJed Brown typedef struct {
8c4762a1bSJed Brown   PetscInt  dim;     /* The topological dimension */
9c4762a1bSJed Brown   PetscBool simplex; /* True for simplices, false for hexes */
10c4762a1bSJed Brown   PetscInt  its;     /* Number of replications for timing */
11c4762a1bSJed Brown   PetscInt  cbs;     /* Number of cells in an integration block */
12c4762a1bSJed Brown } AppCtx;
13c4762a1bSJed Brown 
14c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
15c4762a1bSJed Brown {
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   PetscFunctionBeginUser;
18c4762a1bSJed Brown   options->dim     = 2;
19c4762a1bSJed Brown   options->simplex = PETSC_TRUE;
20c4762a1bSJed Brown   options->its     = 1;
21c4762a1bSJed Brown   options->cbs     = 8;
22c4762a1bSJed Brown 
23d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "FE Integration Performance Options", "PETSCFE");
249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-dim", "The topological dimension", "ex1.c", options->dim, &options->dim, NULL));
259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-simplex", "Simplex or hex cells", "ex1.c", options->simplex, &options->simplex, NULL));
269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-its", "The number of replications for timing", "ex1.c", options->its, &options->its, NULL));
279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-cbs", "The number of cells in an integration block", "ex1.c", options->cbs, &options->cbs, NULL));
28d0609cedSBarry Smith   PetscOptionsEnd();
29c4762a1bSJed Brown   PetscFunctionReturn(0);
30c4762a1bSJed Brown }
31c4762a1bSJed Brown 
32c4762a1bSJed Brown static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
33c4762a1bSJed Brown {
34c4762a1bSJed Brown   PetscInt d;
35c4762a1bSJed Brown   *u = 0.0;
36c4762a1bSJed Brown   for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0*PETSC_PI*x[d]);
37c4762a1bSJed Brown   return 0;
38c4762a1bSJed Brown }
39c4762a1bSJed Brown 
40c4762a1bSJed Brown static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
41c4762a1bSJed Brown                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
42c4762a1bSJed Brown                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
43c4762a1bSJed Brown                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
44c4762a1bSJed Brown {
45c4762a1bSJed Brown   PetscInt d;
46c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f0[0] += -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]);
47c4762a1bSJed Brown }
48c4762a1bSJed Brown 
49c4762a1bSJed Brown static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
50c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
51c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
52c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
53c4762a1bSJed Brown {
54c4762a1bSJed Brown   PetscInt d;
55c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
56c4762a1bSJed Brown }
57c4762a1bSJed Brown 
58c4762a1bSJed Brown static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
59c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
60c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
61c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
62c4762a1bSJed Brown {
63c4762a1bSJed Brown   PetscInt d;
64c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
65c4762a1bSJed Brown }
66c4762a1bSJed Brown 
67c4762a1bSJed Brown static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
68c4762a1bSJed Brown {
69c4762a1bSJed Brown   PetscDS        prob;
7045480ffeSMatthew G. Knepley   DMLabel        label;
71c4762a1bSJed Brown   const PetscInt id = 1;
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   PetscFunctionBeginUser;
749566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &prob));
759566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(prob, 0, f0_trig_u, f1_u));
769566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu));
779566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(prob, 0, trig_u, user));
789566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
799566063dSJacob Faibussowitsch   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) trig_u, NULL, user, NULL));
80c4762a1bSJed Brown   PetscFunctionReturn(0);
81c4762a1bSJed Brown }
82c4762a1bSJed Brown 
83c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
84c4762a1bSJed Brown {
85c4762a1bSJed Brown   DM             cdm = dm;
86c4762a1bSJed Brown   PetscFE        fe;
87c4762a1bSJed Brown   char           prefix[PETSC_MAX_PATH_LEN];
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   PetscFunctionBeginUser;
90c4762a1bSJed Brown   /* Create finite element */
919566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
929566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject) dm), user->dim, 1, user->simplex, name ? prefix : NULL, -1, &fe));
939566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) fe, name));
94c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
959566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe));
969566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
979566063dSJacob Faibussowitsch   PetscCall((*setup)(dm, user));
98c4762a1bSJed Brown   while (cdm) {
999566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm,cdm));
100c4762a1bSJed Brown     /* TODO: Check whether the boundary of coarse meshes is marked */
1019566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
102c4762a1bSJed Brown   }
1039566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
104c4762a1bSJed Brown   PetscFunctionReturn(0);
105c4762a1bSJed Brown }
106c4762a1bSJed Brown 
107c4762a1bSJed Brown static PetscErrorCode PetscContainerUserDestroy_PetscFEGeom(void *ctx)
108c4762a1bSJed Brown {
109c4762a1bSJed Brown   PetscFEGeom   *geom = (PetscFEGeom *) ctx;
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   PetscFunctionBegin;
1129566063dSJacob Faibussowitsch   PetscCall(PetscFEGeomDestroy(&geom));
113c4762a1bSJed Brown   PetscFunctionReturn(0);
114c4762a1bSJed Brown }
115c4762a1bSJed Brown 
116c4762a1bSJed Brown PetscErrorCode CellRangeGetFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
117c4762a1bSJed Brown {
118c4762a1bSJed Brown   char            composeStr[33] = {0};
119c4762a1bSJed Brown   PetscObjectId   id;
120c4762a1bSJed Brown   PetscContainer  container;
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   PetscFunctionBegin;
1239566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetId((PetscObject) quad, &id));
12463a3b9bcSJacob Faibussowitsch   PetscCall(PetscSNPrintf(composeStr, 32, "CellRangeGetFEGeom_%" PetscInt64_FMT "\n", id));
1259566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject) cellIS, composeStr, (PetscObject *) &container));
126c4762a1bSJed Brown   if (container) {
1279566063dSJacob Faibussowitsch     PetscCall(PetscContainerGetPointer(container, (void **) geom));
128c4762a1bSJed Brown   } else {
1299566063dSJacob Faibussowitsch     PetscCall(DMFieldCreateFEGeom(coordField, cellIS, quad, faceData, geom));
1309566063dSJacob Faibussowitsch     PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container));
1319566063dSJacob Faibussowitsch     PetscCall(PetscContainerSetPointer(container, (void *) *geom));
1329566063dSJacob Faibussowitsch     PetscCall(PetscContainerSetUserDestroy(container, PetscContainerUserDestroy_PetscFEGeom));
1339566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject) cellIS, composeStr, (PetscObject) container));
1349566063dSJacob Faibussowitsch     PetscCall(PetscContainerDestroy(&container));
135c4762a1bSJed Brown   }
136c4762a1bSJed Brown   PetscFunctionReturn(0);
137c4762a1bSJed Brown }
138c4762a1bSJed Brown 
139c4762a1bSJed Brown PetscErrorCode CellRangeRestoreFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
140c4762a1bSJed Brown {
141c4762a1bSJed Brown   PetscFunctionBegin;
142c4762a1bSJed Brown   *geom = NULL;
143c4762a1bSJed Brown   PetscFunctionReturn(0);
144c4762a1bSJed Brown }
145c4762a1bSJed Brown 
146c4762a1bSJed Brown static PetscErrorCode CreateFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms)
147c4762a1bSJed Brown {
148c4762a1bSJed Brown   DMField        coordField;
149c4762a1bSJed Brown   PetscInt       Nf, f, maxDegree;
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   PetscFunctionBeginUser;
152c4762a1bSJed Brown   *affineQuad = NULL;
153c4762a1bSJed Brown   *affineGeom = NULL;
154c4762a1bSJed Brown   *quads      = NULL;
155c4762a1bSJed Brown   *geoms      = NULL;
1569566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
1579566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateField(dm, &coordField));
1589566063dSJacob Faibussowitsch   PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
159c4762a1bSJed Brown   if (maxDegree <= 1) {
1609566063dSJacob Faibussowitsch     PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, affineQuad));
1619566063dSJacob Faibussowitsch     if (*affineQuad) PetscCall(CellRangeGetFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom));
162c4762a1bSJed Brown   } else {
1639566063dSJacob Faibussowitsch     PetscCall(PetscCalloc2(Nf, quads, Nf, geoms));
164c4762a1bSJed Brown     for (f = 0; f < Nf; ++f) {
165c4762a1bSJed Brown       PetscFE fe;
166c4762a1bSJed Brown 
1679566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *) &fe));
1689566063dSJacob Faibussowitsch       PetscCall(PetscFEGetQuadrature(fe, &(*quads)[f]));
1699566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject) (*quads)[f]));
1709566063dSJacob Faibussowitsch       PetscCall(CellRangeGetFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f]));
171c4762a1bSJed Brown     }
172c4762a1bSJed Brown   }
173c4762a1bSJed Brown   PetscFunctionReturn(0);
174c4762a1bSJed Brown }
175c4762a1bSJed Brown 
176c4762a1bSJed Brown static PetscErrorCode DestroyFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms)
177c4762a1bSJed Brown {
178c4762a1bSJed Brown   DMField        coordField;
179c4762a1bSJed Brown   PetscInt       Nf, f;
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   PetscFunctionBeginUser;
1829566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
1839566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateField(dm, &coordField));
184c4762a1bSJed Brown   if (*affineQuad) {
1859566063dSJacob Faibussowitsch     PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom));
1869566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureDestroy(affineQuad));
187c4762a1bSJed Brown   } else {
188c4762a1bSJed Brown     for (f = 0; f < Nf; ++f) {
1899566063dSJacob Faibussowitsch       PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f]));
1909566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureDestroy(&(*quads)[f]));
191c4762a1bSJed Brown     }
1929566063dSJacob Faibussowitsch     PetscCall(PetscFree2(*quads, *geoms));
193c4762a1bSJed Brown   }
194c4762a1bSJed Brown   PetscFunctionReturn(0);
195c4762a1bSJed Brown }
196c4762a1bSJed Brown 
197c4762a1bSJed Brown static PetscErrorCode TestIntegration(DM dm, PetscInt cbs, PetscInt its)
198c4762a1bSJed Brown {
199c4762a1bSJed Brown   PetscDS         ds;
200c4762a1bSJed Brown   PetscFEGeom    *chunkGeom = NULL;
201c4762a1bSJed Brown   PetscQuadrature affineQuad,  *quads = NULL;
202c4762a1bSJed Brown   PetscFEGeom    *affineGeom, **geoms = NULL;
203c4762a1bSJed Brown   PetscScalar    *u, *elemVec;
204c4762a1bSJed Brown   IS              cellIS;
205c4762a1bSJed Brown   PetscInt        depth, cStart, cEnd, cell, chunkSize = cbs, Nch = 0, Nf, f, totDim, i, k;
206956f8c0dSBarry Smith #if defined(PETSC_USE_LOG)
207c4762a1bSJed Brown   PetscLogStage   stage;
208c4762a1bSJed Brown   PetscLogEvent   event;
209956f8c0dSBarry Smith #endif
210c4762a1bSJed Brown 
211c4762a1bSJed Brown   PetscFunctionBeginUser;
2129566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("PetscFE Residual Integration Test", &stage));
2139566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("FEIntegRes", PETSCFE_CLASSID, &event));
2149566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stage));
2159566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
2169566063dSJacob Faibussowitsch   PetscCall(DMGetStratumIS(dm, "depth", depth, &cellIS));
2179566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
2189566063dSJacob Faibussowitsch   PetscCall(DMGetCellDS(dm, cStart, &ds));
2199566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
2209566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
2219566063dSJacob Faibussowitsch   PetscCall(CreateFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms));
2229566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(chunkSize*totDim, &u, chunkSize*totDim, &elemVec));
223c4762a1bSJed Brown   /* Assumptions:
224c4762a1bSJed Brown     - Single field
225c4762a1bSJed Brown     - No input data
226c4762a1bSJed Brown     - No auxiliary data
227c4762a1bSJed Brown     - No time-dependence
228c4762a1bSJed Brown   */
229c4762a1bSJed Brown   for (i = 0; i < its; ++i) {
230c4762a1bSJed Brown     for (cell = cStart; cell < cEnd; cell += chunkSize, ++Nch) {
231c4762a1bSJed Brown       const PetscInt cS = cell, cE = PetscMin(cS + chunkSize, cEnd), Ne = cE - cS;
232c4762a1bSJed Brown 
2339566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(elemVec, chunkSize*totDim));
234c4762a1bSJed Brown       /* TODO Replace with DMPlexGetCellFields() */
235c4762a1bSJed Brown       for (k = 0; k < chunkSize*totDim; ++k) u[k] = 1.0;
236c4762a1bSJed Brown       for (f = 0; f < Nf; ++f) {
23706ad1575SMatthew G. Knepley         PetscFormKey key;
238c4762a1bSJed Brown         PetscFEGeom     *geom = affineGeom ? affineGeom : geoms[f];
239c4762a1bSJed Brown         /* PetscQuadrature quad = affineQuad ? affineQuad : quads[f]; */
240c4762a1bSJed Brown 
2410816aa78SBarry Smith         key.label = NULL; key.value = 0; key.field = f; key.part = 0;
2429566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomGetChunk(geom, cS, cE, &chunkGeom));
2439566063dSJacob Faibussowitsch         PetscCall(PetscLogEventBegin(event,0,0,0,0));
2449566063dSJacob Faibussowitsch         PetscCall(PetscFEIntegrateResidual(ds, key, Ne, chunkGeom, u, NULL, NULL, NULL, 0.0, elemVec));
2459566063dSJacob Faibussowitsch         PetscCall(PetscLogEventEnd(event,0,0,0,0));
246c4762a1bSJed Brown       }
247c4762a1bSJed Brown     }
248c4762a1bSJed Brown   }
2499566063dSJacob Faibussowitsch   PetscCall(PetscFEGeomRestoreChunk(affineGeom, cStart, cEnd, &chunkGeom));
2509566063dSJacob Faibussowitsch   PetscCall(DestroyFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms));
2519566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cellIS));
2529566063dSJacob Faibussowitsch   PetscCall(PetscFree2(u, elemVec));
2539566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
254956f8c0dSBarry Smith #if defined(PETSC_USE_LOG)
255c4762a1bSJed Brown   {
256c4762a1bSJed Brown     const char        *title = "Petsc FE Residual Integration";
257c4762a1bSJed Brown     PetscEventPerfInfo eventInfo;
258c4762a1bSJed Brown     PetscInt           N = (cEnd - cStart)*Nf*its;
259c4762a1bSJed Brown     PetscReal          flopRate, cellRate;
260c4762a1bSJed Brown 
2619566063dSJacob Faibussowitsch     PetscCall(PetscLogEventGetPerfInfo(stage, event, &eventInfo));
262c4762a1bSJed Brown     flopRate = eventInfo.time != 0.0 ? eventInfo.flops/eventInfo.time : 0.0;
263c4762a1bSJed Brown     cellRate = eventInfo.time != 0.0 ? N/eventInfo.time : 0.0;
26463a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject) dm), "%s: %" PetscInt_FMT " integrals %" PetscInt_FMT " chunks %" PetscInt_FMT " reps\n  Cell rate: %.2f/s flop rate: %.2f MF/s\n", title, N, Nch, its, (double)cellRate, (double)(flopRate/1.e6)));
265c4762a1bSJed Brown   }
266956f8c0dSBarry Smith #endif
267c4762a1bSJed Brown   PetscFunctionReturn(0);
268c4762a1bSJed Brown }
269c4762a1bSJed Brown 
270c4762a1bSJed Brown static PetscErrorCode TestIntegration2(DM dm, PetscInt cbs, PetscInt its)
271c4762a1bSJed Brown {
272c4762a1bSJed Brown   Vec             X, F;
273956f8c0dSBarry Smith #if defined(PETSC_USE_LOG)
274c4762a1bSJed Brown   PetscLogStage   stage;
275956f8c0dSBarry Smith #endif
276c4762a1bSJed Brown   PetscInt        i;
277c4762a1bSJed Brown 
278c4762a1bSJed Brown   PetscFunctionBeginUser;
2799566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("DMPlex Residual Integration Test", &stage));
2809566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stage));
2819566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &X));
2829566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &F));
283c4762a1bSJed Brown   for (i = 0; i < its; ++i) {
2849566063dSJacob Faibussowitsch     PetscCall(DMPlexSNESComputeResidualFEM(dm, X, F, NULL));
285c4762a1bSJed Brown   }
2869566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &X));
2879566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &F));
2889566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
289956f8c0dSBarry Smith #if defined(PETSC_USE_LOG)
290c4762a1bSJed Brown   {
291c4762a1bSJed Brown     const char         *title = "DMPlex Residual Integration";
292c4762a1bSJed Brown     PetscEventPerfInfo eventInfo;
293c4762a1bSJed Brown     PetscReal          flopRate, cellRate;
294c4762a1bSJed Brown     PetscInt           cStart, cEnd, Nf, N;
295956f8c0dSBarry Smith     PetscLogEvent      event;
296c4762a1bSJed Brown 
2979566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
2989566063dSJacob Faibussowitsch     PetscCall(DMGetNumFields(dm, &Nf));
2999566063dSJacob Faibussowitsch     PetscCall(PetscLogEventGetId("DMPlexResidualFE", &event));
3009566063dSJacob Faibussowitsch     PetscCall(PetscLogEventGetPerfInfo(stage, event, &eventInfo));
301c4762a1bSJed Brown     N        = (cEnd - cStart)*Nf*eventInfo.count;
302c4762a1bSJed Brown     flopRate = eventInfo.time != 0.0 ? eventInfo.flops/eventInfo.time : 0.0;
303c4762a1bSJed Brown     cellRate = eventInfo.time != 0.0 ? N/eventInfo.time : 0.0;
30463a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject) dm), "%s: %" PetscInt_FMT " integrals %d reps\n  Cell rate: %.2f/s flop rate: %.2f MF/s\n", title, N, eventInfo.count, (double)cellRate, (double)(flopRate/1.e6)));
305c4762a1bSJed Brown   }
306956f8c0dSBarry Smith #endif
307c4762a1bSJed Brown   PetscFunctionReturn(0);
308c4762a1bSJed Brown }
309c4762a1bSJed Brown 
310c4762a1bSJed Brown int main(int argc, char **argv)
311c4762a1bSJed Brown {
312c4762a1bSJed Brown   DM             dm;
313c4762a1bSJed Brown   AppCtx         ctx;
314c4762a1bSJed Brown   PetscMPIInt    size;
315c4762a1bSJed Brown 
316*327415f7SBarry Smith   PetscFunctionBeginUser;
3179566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
31908401ef6SPierre Jolivet   PetscCheck(size <= 1,PETSC_COMM_WORLD, PETSC_ERR_SUP, "This is a uniprocessor example only.");
3209566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
3219566063dSJacob Faibussowitsch   PetscCall(PetscLogDefaultBegin());
3229566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
3239566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
3249566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
3259566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) dm, "Mesh"));
3269566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject) dm, NULL, "-dm_view"));
3279566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &ctx));
3289566063dSJacob Faibussowitsch   PetscCall(TestIntegration(dm, ctx.cbs, ctx.its));
3299566063dSJacob Faibussowitsch   PetscCall(TestIntegration2(dm, ctx.cbs, ctx.its));
3309566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
3319566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
332b122ec5aSJacob Faibussowitsch   return 0;
333c4762a1bSJed Brown }
334c4762a1bSJed Brown 
335c4762a1bSJed Brown /*TEST
336c4762a1bSJed Brown   test:
337c4762a1bSJed Brown     suffix: 0
338c4762a1bSJed Brown     requires: triangle
339c4762a1bSJed Brown     args: -dm_view
340c4762a1bSJed Brown 
341c4762a1bSJed Brown   test:
342c4762a1bSJed Brown     suffix: 1
343c4762a1bSJed Brown     requires: triangle
344c4762a1bSJed Brown     args: -dm_view -potential_petscspace_degree 1
345c4762a1bSJed Brown 
346c4762a1bSJed Brown   test:
347c4762a1bSJed Brown     suffix: 2
348c4762a1bSJed Brown     requires: triangle
349c4762a1bSJed Brown     args: -dm_view -potential_petscspace_degree 2
350c4762a1bSJed Brown 
351c4762a1bSJed Brown   test:
352c4762a1bSJed Brown     suffix: 3
353c4762a1bSJed Brown     requires: triangle
354c4762a1bSJed Brown     args: -dm_view -potential_petscspace_degree 3
355c4762a1bSJed Brown TEST*/
356