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