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