xref: /petsc/src/dm/dt/fe/tests/ex1.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
149371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
15c4762a1bSJed Brown   PetscFunctionBeginUser;
16c4762a1bSJed Brown   options->dim     = 2;
17c4762a1bSJed Brown   options->simplex = PETSC_TRUE;
18c4762a1bSJed Brown   options->its     = 1;
19c4762a1bSJed Brown   options->cbs     = 8;
20c4762a1bSJed Brown 
21d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "FE Integration Performance Options", "PETSCFE");
229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-dim", "The topological dimension", "ex1.c", options->dim, &options->dim, NULL));
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-simplex", "Simplex or hex cells", "ex1.c", options->simplex, &options->simplex, NULL));
249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-its", "The number of replications for timing", "ex1.c", options->its, &options->its, NULL));
259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-cbs", "The number of cells in an integration block", "ex1.c", options->cbs, &options->cbs, NULL));
26d0609cedSBarry Smith   PetscOptionsEnd();
27c4762a1bSJed Brown   PetscFunctionReturn(0);
28c4762a1bSJed Brown }
29c4762a1bSJed Brown 
309371c9d4SSatish Balay static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
31c4762a1bSJed Brown   PetscInt d;
32c4762a1bSJed Brown   *u = 0.0;
33c4762a1bSJed Brown   for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0 * PETSC_PI * x[d]);
34c4762a1bSJed Brown   return 0;
35c4762a1bSJed Brown }
36c4762a1bSJed Brown 
379371c9d4SSatish Balay static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
38c4762a1bSJed Brown   PetscInt d;
39c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f0[0] += -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]);
40c4762a1bSJed Brown }
41c4762a1bSJed Brown 
429371c9d4SSatish Balay static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) {
43c4762a1bSJed Brown   PetscInt d;
44c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
45c4762a1bSJed Brown }
46c4762a1bSJed Brown 
479371c9d4SSatish Balay static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) {
48c4762a1bSJed Brown   PetscInt d;
49c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
50c4762a1bSJed Brown }
51c4762a1bSJed Brown 
529371c9d4SSatish Balay static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) {
53c4762a1bSJed Brown   PetscDS        prob;
5445480ffeSMatthew G. Knepley   DMLabel        label;
55c4762a1bSJed Brown   const PetscInt id = 1;
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   PetscFunctionBeginUser;
589566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &prob));
599566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(prob, 0, f0_trig_u, f1_u));
609566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu));
619566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(prob, 0, trig_u, user));
629566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
639566063dSJacob Faibussowitsch   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))trig_u, NULL, user, NULL));
64c4762a1bSJed Brown   PetscFunctionReturn(0);
65c4762a1bSJed Brown }
66c4762a1bSJed Brown 
679371c9d4SSatish Balay static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) {
68c4762a1bSJed Brown   DM      cdm = dm;
69c4762a1bSJed Brown   PetscFE fe;
70c4762a1bSJed Brown   char    prefix[PETSC_MAX_PATH_LEN];
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   PetscFunctionBeginUser;
73c4762a1bSJed Brown   /* Create finite element */
749566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
759566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), user->dim, 1, user->simplex, name ? prefix : NULL, -1, &fe));
769566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe, name));
77c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
789566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
799566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
809566063dSJacob Faibussowitsch   PetscCall((*setup)(dm, user));
81c4762a1bSJed Brown   while (cdm) {
829566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm, cdm));
83c4762a1bSJed Brown     /* TODO: Check whether the boundary of coarse meshes is marked */
849566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
85c4762a1bSJed Brown   }
869566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
87c4762a1bSJed Brown   PetscFunctionReturn(0);
88c4762a1bSJed Brown }
89c4762a1bSJed Brown 
909371c9d4SSatish Balay static PetscErrorCode PetscContainerUserDestroy_PetscFEGeom(void *ctx) {
91c4762a1bSJed Brown   PetscFEGeom *geom = (PetscFEGeom *)ctx;
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   PetscFunctionBegin;
949566063dSJacob Faibussowitsch   PetscCall(PetscFEGeomDestroy(&geom));
95c4762a1bSJed Brown   PetscFunctionReturn(0);
96c4762a1bSJed Brown }
97c4762a1bSJed Brown 
989371c9d4SSatish Balay PetscErrorCode CellRangeGetFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) {
99c4762a1bSJed Brown   char           composeStr[33] = {0};
100c4762a1bSJed Brown   PetscObjectId  id;
101c4762a1bSJed Brown   PetscContainer container;
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   PetscFunctionBegin;
1049566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetId((PetscObject)quad, &id));
10563a3b9bcSJacob Faibussowitsch   PetscCall(PetscSNPrintf(composeStr, 32, "CellRangeGetFEGeom_%" PetscInt64_FMT "\n", id));
1069566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)cellIS, composeStr, (PetscObject *)&container));
107c4762a1bSJed Brown   if (container) {
1089566063dSJacob Faibussowitsch     PetscCall(PetscContainerGetPointer(container, (void **)geom));
109c4762a1bSJed Brown   } else {
1109566063dSJacob Faibussowitsch     PetscCall(DMFieldCreateFEGeom(coordField, cellIS, quad, faceData, geom));
1119566063dSJacob Faibussowitsch     PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container));
1129566063dSJacob Faibussowitsch     PetscCall(PetscContainerSetPointer(container, (void *)*geom));
1139566063dSJacob Faibussowitsch     PetscCall(PetscContainerSetUserDestroy(container, PetscContainerUserDestroy_PetscFEGeom));
1149566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)cellIS, composeStr, (PetscObject)container));
1159566063dSJacob Faibussowitsch     PetscCall(PetscContainerDestroy(&container));
116c4762a1bSJed Brown   }
117c4762a1bSJed Brown   PetscFunctionReturn(0);
118c4762a1bSJed Brown }
119c4762a1bSJed Brown 
1209371c9d4SSatish Balay PetscErrorCode CellRangeRestoreFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) {
121c4762a1bSJed Brown   PetscFunctionBegin;
122c4762a1bSJed Brown   *geom = NULL;
123c4762a1bSJed Brown   PetscFunctionReturn(0);
124c4762a1bSJed Brown }
125c4762a1bSJed Brown 
1269371c9d4SSatish Balay static PetscErrorCode CreateFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms) {
127c4762a1bSJed Brown   DMField  coordField;
128c4762a1bSJed Brown   PetscInt Nf, f, maxDegree;
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   PetscFunctionBeginUser;
131c4762a1bSJed Brown   *affineQuad = NULL;
132c4762a1bSJed Brown   *affineGeom = NULL;
133c4762a1bSJed Brown   *quads      = NULL;
134c4762a1bSJed Brown   *geoms      = NULL;
1359566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
1369566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateField(dm, &coordField));
1379566063dSJacob Faibussowitsch   PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
138c4762a1bSJed Brown   if (maxDegree <= 1) {
1399566063dSJacob Faibussowitsch     PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, affineQuad));
1409566063dSJacob Faibussowitsch     if (*affineQuad) PetscCall(CellRangeGetFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom));
141c4762a1bSJed Brown   } else {
1429566063dSJacob Faibussowitsch     PetscCall(PetscCalloc2(Nf, quads, Nf, geoms));
143c4762a1bSJed Brown     for (f = 0; f < Nf; ++f) {
144c4762a1bSJed Brown       PetscFE fe;
145c4762a1bSJed Brown 
1469566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
1479566063dSJacob Faibussowitsch       PetscCall(PetscFEGetQuadrature(fe, &(*quads)[f]));
1489566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)(*quads)[f]));
1499566063dSJacob Faibussowitsch       PetscCall(CellRangeGetFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f]));
150c4762a1bSJed Brown     }
151c4762a1bSJed Brown   }
152c4762a1bSJed Brown   PetscFunctionReturn(0);
153c4762a1bSJed Brown }
154c4762a1bSJed Brown 
1559371c9d4SSatish Balay static PetscErrorCode DestroyFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms) {
156c4762a1bSJed Brown   DMField  coordField;
157c4762a1bSJed Brown   PetscInt Nf, f;
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   PetscFunctionBeginUser;
1609566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
1619566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateField(dm, &coordField));
162c4762a1bSJed Brown   if (*affineQuad) {
1639566063dSJacob Faibussowitsch     PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom));
1649566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureDestroy(affineQuad));
165c4762a1bSJed Brown   } else {
166c4762a1bSJed Brown     for (f = 0; f < Nf; ++f) {
1679566063dSJacob Faibussowitsch       PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f]));
1689566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureDestroy(&(*quads)[f]));
169c4762a1bSJed Brown     }
1709566063dSJacob Faibussowitsch     PetscCall(PetscFree2(*quads, *geoms));
171c4762a1bSJed Brown   }
172c4762a1bSJed Brown   PetscFunctionReturn(0);
173c4762a1bSJed Brown }
174c4762a1bSJed Brown 
1759371c9d4SSatish Balay static PetscErrorCode TestIntegration(DM dm, PetscInt cbs, PetscInt its) {
176c4762a1bSJed Brown   PetscDS         ds;
177c4762a1bSJed Brown   PetscFEGeom    *chunkGeom = NULL;
178c4762a1bSJed Brown   PetscQuadrature affineQuad, *quads  = NULL;
179c4762a1bSJed Brown   PetscFEGeom    *affineGeom, **geoms = NULL;
180c4762a1bSJed Brown   PetscScalar    *u, *elemVec;
181c4762a1bSJed Brown   IS              cellIS;
182c4762a1bSJed Brown   PetscInt        depth, cStart, cEnd, cell, chunkSize = cbs, Nch = 0, Nf, f, totDim, i, k;
183956f8c0dSBarry Smith #if defined(PETSC_USE_LOG)
184c4762a1bSJed Brown   PetscLogStage stage;
185c4762a1bSJed Brown   PetscLogEvent event;
186956f8c0dSBarry Smith #endif
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   PetscFunctionBeginUser;
1899566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("PetscFE Residual Integration Test", &stage));
1909566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("FEIntegRes", PETSCFE_CLASSID, &event));
1919566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stage));
1929566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
1939566063dSJacob Faibussowitsch   PetscCall(DMGetStratumIS(dm, "depth", depth, &cellIS));
1949566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1959566063dSJacob Faibussowitsch   PetscCall(DMGetCellDS(dm, cStart, &ds));
1969566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
1979566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
1989566063dSJacob Faibussowitsch   PetscCall(CreateFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms));
1999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(chunkSize * totDim, &u, chunkSize * totDim, &elemVec));
200c4762a1bSJed Brown   /* Assumptions:
201c4762a1bSJed Brown     - Single field
202c4762a1bSJed Brown     - No input data
203c4762a1bSJed Brown     - No auxiliary data
204c4762a1bSJed Brown     - No time-dependence
205c4762a1bSJed Brown   */
206c4762a1bSJed Brown   for (i = 0; i < its; ++i) {
207c4762a1bSJed Brown     for (cell = cStart; cell < cEnd; cell += chunkSize, ++Nch) {
208c4762a1bSJed Brown       const PetscInt cS = cell, cE = PetscMin(cS + chunkSize, cEnd), Ne = cE - cS;
209c4762a1bSJed Brown 
2109566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(elemVec, chunkSize * totDim));
211c4762a1bSJed Brown       /* TODO Replace with DMPlexGetCellFields() */
212c4762a1bSJed Brown       for (k = 0; k < chunkSize * totDim; ++k) u[k] = 1.0;
213c4762a1bSJed Brown       for (f = 0; f < Nf; ++f) {
21406ad1575SMatthew G. Knepley         PetscFormKey key;
215c4762a1bSJed Brown         PetscFEGeom *geom = affineGeom ? affineGeom : geoms[f];
216c4762a1bSJed Brown         /* PetscQuadrature quad = affineQuad ? affineQuad : quads[f]; */
217c4762a1bSJed Brown 
2189371c9d4SSatish Balay         key.label = NULL;
2199371c9d4SSatish Balay         key.value = 0;
2209371c9d4SSatish Balay         key.field = f;
2219371c9d4SSatish Balay         key.part  = 0;
2229566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomGetChunk(geom, cS, cE, &chunkGeom));
2239566063dSJacob Faibussowitsch         PetscCall(PetscLogEventBegin(event, 0, 0, 0, 0));
2249566063dSJacob Faibussowitsch         PetscCall(PetscFEIntegrateResidual(ds, key, Ne, chunkGeom, u, NULL, NULL, NULL, 0.0, elemVec));
2259566063dSJacob Faibussowitsch         PetscCall(PetscLogEventEnd(event, 0, 0, 0, 0));
226c4762a1bSJed Brown       }
227c4762a1bSJed Brown     }
228c4762a1bSJed Brown   }
2299566063dSJacob Faibussowitsch   PetscCall(PetscFEGeomRestoreChunk(affineGeom, cStart, cEnd, &chunkGeom));
2309566063dSJacob Faibussowitsch   PetscCall(DestroyFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms));
2319566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cellIS));
2329566063dSJacob Faibussowitsch   PetscCall(PetscFree2(u, elemVec));
2339566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
234956f8c0dSBarry Smith #if defined(PETSC_USE_LOG)
235c4762a1bSJed Brown   {
236c4762a1bSJed Brown     const char        *title = "Petsc FE Residual Integration";
237c4762a1bSJed Brown     PetscEventPerfInfo eventInfo;
238c4762a1bSJed Brown     PetscInt           N = (cEnd - cStart) * Nf * its;
239c4762a1bSJed Brown     PetscReal          flopRate, cellRate;
240c4762a1bSJed Brown 
2419566063dSJacob Faibussowitsch     PetscCall(PetscLogEventGetPerfInfo(stage, event, &eventInfo));
242c4762a1bSJed Brown     flopRate = eventInfo.time != 0.0 ? eventInfo.flops / eventInfo.time : 0.0;
243c4762a1bSJed Brown     cellRate = eventInfo.time != 0.0 ? N / eventInfo.time : 0.0;
24463a3b9bcSJacob 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)));
245c4762a1bSJed Brown   }
246956f8c0dSBarry Smith #endif
247c4762a1bSJed Brown   PetscFunctionReturn(0);
248c4762a1bSJed Brown }
249c4762a1bSJed Brown 
2509371c9d4SSatish Balay static PetscErrorCode TestIntegration2(DM dm, PetscInt cbs, PetscInt its) {
251c4762a1bSJed Brown   Vec X, F;
252956f8c0dSBarry Smith #if defined(PETSC_USE_LOG)
253c4762a1bSJed Brown   PetscLogStage stage;
254956f8c0dSBarry Smith #endif
255c4762a1bSJed Brown   PetscInt i;
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   PetscFunctionBeginUser;
2589566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("DMPlex Residual Integration Test", &stage));
2599566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stage));
2609566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &X));
2619566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &F));
262*48a46eb9SPierre Jolivet   for (i = 0; i < its; ++i) PetscCall(DMPlexSNESComputeResidualFEM(dm, X, F, NULL));
2639566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &X));
2649566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &F));
2659566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
266956f8c0dSBarry Smith #if defined(PETSC_USE_LOG)
267c4762a1bSJed Brown   {
268c4762a1bSJed Brown     const char        *title = "DMPlex Residual Integration";
269c4762a1bSJed Brown     PetscEventPerfInfo eventInfo;
270c4762a1bSJed Brown     PetscReal          flopRate, cellRate;
271c4762a1bSJed Brown     PetscInt           cStart, cEnd, Nf, N;
272956f8c0dSBarry Smith     PetscLogEvent      event;
273c4762a1bSJed Brown 
2749566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
2759566063dSJacob Faibussowitsch     PetscCall(DMGetNumFields(dm, &Nf));
2769566063dSJacob Faibussowitsch     PetscCall(PetscLogEventGetId("DMPlexResidualFE", &event));
2779566063dSJacob Faibussowitsch     PetscCall(PetscLogEventGetPerfInfo(stage, event, &eventInfo));
278c4762a1bSJed Brown     N        = (cEnd - cStart) * Nf * eventInfo.count;
279c4762a1bSJed Brown     flopRate = eventInfo.time != 0.0 ? eventInfo.flops / eventInfo.time : 0.0;
280c4762a1bSJed Brown     cellRate = eventInfo.time != 0.0 ? N / eventInfo.time : 0.0;
28163a3b9bcSJacob 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)));
282c4762a1bSJed Brown   }
283956f8c0dSBarry Smith #endif
284c4762a1bSJed Brown   PetscFunctionReturn(0);
285c4762a1bSJed Brown }
286c4762a1bSJed Brown 
2879371c9d4SSatish Balay int main(int argc, char **argv) {
288c4762a1bSJed Brown   DM          dm;
289c4762a1bSJed Brown   AppCtx      ctx;
290c4762a1bSJed Brown   PetscMPIInt size;
291c4762a1bSJed Brown 
292327415f7SBarry Smith   PetscFunctionBeginUser;
2939566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2949566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
29508401ef6SPierre Jolivet   PetscCheck(size <= 1, PETSC_COMM_WORLD, PETSC_ERR_SUP, "This is a uniprocessor example only.");
2969566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
2979566063dSJacob Faibussowitsch   PetscCall(PetscLogDefaultBegin());
2989566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
2999566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
3009566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
3019566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)dm, "Mesh"));
3029566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)dm, NULL, "-dm_view"));
3039566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &ctx));
3049566063dSJacob Faibussowitsch   PetscCall(TestIntegration(dm, ctx.cbs, ctx.its));
3059566063dSJacob Faibussowitsch   PetscCall(TestIntegration2(dm, ctx.cbs, ctx.its));
3069566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
3079566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
308b122ec5aSJacob Faibussowitsch   return 0;
309c4762a1bSJed Brown }
310c4762a1bSJed Brown 
311c4762a1bSJed Brown /*TEST
312c4762a1bSJed Brown   test:
313c4762a1bSJed Brown     suffix: 0
314c4762a1bSJed Brown     requires: triangle
315c4762a1bSJed Brown     args: -dm_view
316c4762a1bSJed Brown 
317c4762a1bSJed Brown   test:
318c4762a1bSJed Brown     suffix: 1
319c4762a1bSJed Brown     requires: triangle
320c4762a1bSJed Brown     args: -dm_view -potential_petscspace_degree 1
321c4762a1bSJed Brown 
322c4762a1bSJed Brown   test:
323c4762a1bSJed Brown     suffix: 2
324c4762a1bSJed Brown     requires: triangle
325c4762a1bSJed Brown     args: -dm_view -potential_petscspace_degree 2
326c4762a1bSJed Brown 
327c4762a1bSJed Brown   test:
328c4762a1bSJed Brown     suffix: 3
329c4762a1bSJed Brown     requires: triangle
330c4762a1bSJed Brown     args: -dm_view -potential_petscspace_degree 3
331c4762a1bSJed Brown TEST*/
332