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