Actual source code: ex1.c

  1: static const char help[] = "Performance Tests for FE Integration";

  3: #include <petscdmplex.h>
  4: #include <petscfe.h>
  5: #include <petscds.h>

  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;

 14: static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
 15: {
 16:   PetscFunctionBeginUser;
 17:   options->dim     = 2;
 18:   options->simplex = PETSC_TRUE;
 19:   options->its     = 1;
 20:   options->cbs     = 8;

 22:   PetscOptionsBegin(comm, "", "FE Integration Performance Options", "PETSCFE");
 23:   PetscCall(PetscOptionsInt("-dim", "The topological dimension", "ex1.c", options->dim, &options->dim, NULL));
 24:   PetscCall(PetscOptionsBool("-simplex", "Simplex or hex cells", "ex1.c", options->simplex, &options->simplex, NULL));
 25:   PetscCall(PetscOptionsInt("-its", "The number of replications for timing", "ex1.c", options->its, &options->its, NULL));
 26:   PetscCall(PetscOptionsInt("-cbs", "The number of cells in an integration block", "ex1.c", options->cbs, &options->cbs, NULL));
 27:   PetscOptionsEnd();
 28:   PetscFunctionReturn(PETSC_SUCCESS);
 29: }

 31: static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
 32: {
 33:   PetscInt d;
 34:   *u = 0.0;
 35:   for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0 * PETSC_PI * x[d]);
 36:   return PETSC_SUCCESS;
 37: }

 39: 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[])
 40: {
 41:   PetscInt d;
 42:   for (d = 0; d < dim; ++d) f0[0] += -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]);
 43: }

 45: 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[])
 46: {
 47:   PetscInt d;
 48:   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
 49: }

 51: 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[])
 52: {
 53:   PetscInt d;
 54:   for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
 55: }

 57: static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
 58: {
 59:   PetscDS        prob;
 60:   DMLabel        label;
 61:   const PetscInt id = 1;

 63:   PetscFunctionBeginUser;
 64:   PetscCall(DMGetDS(dm, &prob));
 65:   PetscCall(PetscDSSetResidual(prob, 0, f0_trig_u, f1_u));
 66:   PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu));
 67:   PetscCall(PetscDSSetExactSolution(prob, 0, trig_u, user));
 68:   PetscCall(DMGetLabel(dm, "marker", &label));
 69:   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))trig_u, NULL, user, NULL));
 70:   PetscFunctionReturn(PETSC_SUCCESS);
 71: }

 73: static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
 74: {
 75:   DM      cdm = dm;
 76:   PetscFE fe;
 77:   char    prefix[PETSC_MAX_PATH_LEN];

 79:   PetscFunctionBeginUser;
 80:   /* Create finite element */
 81:   PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
 82:   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), user->dim, 1, user->simplex, name ? prefix : NULL, -1, &fe));
 83:   PetscCall(PetscObjectSetName((PetscObject)fe, name));
 84:   /* Set discretization and boundary conditions for each mesh */
 85:   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
 86:   PetscCall(DMCreateDS(dm));
 87:   PetscCall((*setup)(dm, user));
 88:   while (cdm) {
 89:     PetscCall(DMCopyDisc(dm, cdm));
 90:     /* TODO: Check whether the boundary of coarse meshes is marked */
 91:     PetscCall(DMGetCoarseDM(cdm, &cdm));
 92:   }
 93:   PetscCall(PetscFEDestroy(&fe));
 94:   PetscFunctionReturn(PETSC_SUCCESS);
 95: }

 97: static PetscErrorCode PetscContainerUserDestroy_PetscFEGeom(void *ctx)
 98: {
 99:   PetscFEGeom *geom = (PetscFEGeom *)ctx;

101:   PetscFunctionBegin;
102:   PetscCall(PetscFEGeomDestroy(&geom));
103:   PetscFunctionReturn(PETSC_SUCCESS);
104: }

106: PetscErrorCode CellRangeGetFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
107: {
108:   char           composeStr[33] = {0};
109:   PetscObjectId  id;
110:   PetscContainer container;

112:   PetscFunctionBegin;
113:   PetscCall(PetscObjectGetId((PetscObject)quad, &id));
114:   PetscCall(PetscSNPrintf(composeStr, 32, "CellRangeGetFEGeom_%" PetscInt64_FMT "\n", id));
115:   PetscCall(PetscObjectQuery((PetscObject)cellIS, composeStr, (PetscObject *)&container));
116:   if (container) {
117:     PetscCall(PetscContainerGetPointer(container, (void **)geom));
118:   } else {
119:     PetscCall(DMFieldCreateFEGeom(coordField, cellIS, quad, faceData, geom));
120:     PetscCall(PetscObjectContainerCompose((PetscObject)cellIS, composeStr, *geom, PetscContainerUserDestroy_PetscFEGeom));
121:   }
122:   PetscFunctionReturn(PETSC_SUCCESS);
123: }

125: PetscErrorCode CellRangeRestoreFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
126: {
127:   PetscFunctionBegin;
128:   *geom = NULL;
129:   PetscFunctionReturn(PETSC_SUCCESS);
130: }

132: static PetscErrorCode CreateFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms)
133: {
134:   DMField  coordField;
135:   PetscInt Nf, f, maxDegree;

137:   PetscFunctionBeginUser;
138:   *affineQuad = NULL;
139:   *affineGeom = NULL;
140:   *quads      = NULL;
141:   *geoms      = NULL;
142:   PetscCall(PetscDSGetNumFields(ds, &Nf));
143:   PetscCall(DMGetCoordinateField(dm, &coordField));
144:   PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
145:   if (maxDegree <= 1) {
146:     PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, affineQuad));
147:     if (*affineQuad) PetscCall(CellRangeGetFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom));
148:   } else {
149:     PetscCall(PetscCalloc2(Nf, quads, Nf, geoms));
150:     for (f = 0; f < Nf; ++f) {
151:       PetscFE fe;

153:       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
154:       PetscCall(PetscFEGetQuadrature(fe, &(*quads)[f]));
155:       PetscCall(PetscObjectReference((PetscObject)(*quads)[f]));
156:       PetscCall(CellRangeGetFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f]));
157:     }
158:   }
159:   PetscFunctionReturn(PETSC_SUCCESS);
160: }

162: static PetscErrorCode DestroyFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms)
163: {
164:   DMField  coordField;
165:   PetscInt Nf, f;

167:   PetscFunctionBeginUser;
168:   PetscCall(PetscDSGetNumFields(ds, &Nf));
169:   PetscCall(DMGetCoordinateField(dm, &coordField));
170:   if (*affineQuad) {
171:     PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom));
172:     PetscCall(PetscQuadratureDestroy(affineQuad));
173:   } else {
174:     for (f = 0; f < Nf; ++f) {
175:       PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f]));
176:       PetscCall(PetscQuadratureDestroy(&(*quads)[f]));
177:     }
178:     PetscCall(PetscFree2(*quads, *geoms));
179:   }
180:   PetscFunctionReturn(PETSC_SUCCESS);
181: }

183: static PetscErrorCode TestIntegration(DM dm, PetscInt cbs, PetscInt its)
184: {
185:   PetscDS         ds;
186:   PetscFEGeom    *chunkGeom = NULL;
187:   PetscQuadrature affineQuad, *quads  = NULL;
188:   PetscFEGeom    *affineGeom, **geoms = NULL;
189:   PetscScalar    *u, *elemVec;
190:   IS              cellIS;
191:   PetscInt        depth, cStart, cEnd, cell, chunkSize = cbs, Nch = 0, Nf, f, totDim, i, k;
192:   PetscLogStage   stage;
193:   PetscLogEvent   event;

195:   PetscFunctionBeginUser;
196:   PetscCall(PetscLogStageRegister("PetscFE Residual Integration Test", &stage));
197:   PetscCall(PetscLogEventRegister("FEIntegRes", PETSCFE_CLASSID, &event));
198:   PetscCall(PetscLogStagePush(stage));
199:   PetscCall(DMPlexGetDepth(dm, &depth));
200:   PetscCall(DMGetStratumIS(dm, "depth", depth, &cellIS));
201:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
202:   PetscCall(DMGetCellDS(dm, cStart, &ds, NULL));
203:   PetscCall(PetscDSGetNumFields(ds, &Nf));
204:   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
205:   PetscCall(CreateFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms));
206:   PetscCall(PetscMalloc2(chunkSize * totDim, &u, chunkSize * totDim, &elemVec));
207:   /* Assumptions:
208:     - Single field
209:     - No input data
210:     - No auxiliary data
211:     - No time-dependence
212:   */
213:   for (i = 0; i < its; ++i) {
214:     for (cell = cStart; cell < cEnd; cell += chunkSize, ++Nch) {
215:       const PetscInt cS = cell, cE = PetscMin(cS + chunkSize, cEnd), Ne = cE - cS;

217:       PetscCall(PetscArrayzero(elemVec, chunkSize * totDim));
218:       /* TODO Replace with DMPlexGetCellFields() */
219:       for (k = 0; k < chunkSize * totDim; ++k) u[k] = 1.0;
220:       for (f = 0; f < Nf; ++f) {
221:         PetscFormKey key;
222:         PetscFEGeom *geom = affineGeom ? affineGeom : geoms[f];
223:         /* PetscQuadrature quad = affineQuad ? affineQuad : quads[f]; */

225:         key.label = NULL;
226:         key.value = 0;
227:         key.field = f;
228:         key.part  = 0;
229:         PetscCall(PetscFEGeomGetChunk(geom, cS, cE, &chunkGeom));
230:         PetscCall(PetscLogEventBegin(event, 0, 0, 0, 0));
231:         PetscCall(PetscFEIntegrateResidual(ds, key, Ne, chunkGeom, u, NULL, NULL, NULL, 0.0, elemVec));
232:         PetscCall(PetscLogEventEnd(event, 0, 0, 0, 0));
233:       }
234:     }
235:   }
236:   PetscCall(PetscFEGeomRestoreChunk(affineGeom, cStart, cEnd, &chunkGeom));
237:   PetscCall(DestroyFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms));
238:   PetscCall(ISDestroy(&cellIS));
239:   PetscCall(PetscFree2(u, elemVec));
240:   PetscCall(PetscLogStagePop());
241:   if (PetscDefined(USE_LOG)) {
242:     const char        *title = "Petsc FE Residual Integration";
243:     PetscEventPerfInfo eventInfo;
244:     PetscInt           N = (cEnd - cStart) * Nf * its;
245:     PetscReal          flopRate, cellRate;

247:     PetscCall(PetscLogEventGetPerfInfo(stage, event, &eventInfo));
248:     flopRate = eventInfo.time != 0.0 ? eventInfo.flops / eventInfo.time : 0.0;
249:     cellRate = eventInfo.time != 0.0 ? N / eventInfo.time : 0.0;
250:     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)));
251:   }
252:   PetscFunctionReturn(PETSC_SUCCESS);
253: }

255: static PetscErrorCode TestIntegration2(DM dm, PetscInt cbs, PetscInt its)
256: {
257:   Vec           X, F;
258:   PetscLogStage stage;
259:   PetscInt      i;

261:   PetscFunctionBeginUser;
262:   PetscCall(PetscLogStageRegister("DMPlex Residual Integration Test", &stage));
263:   PetscCall(PetscLogStagePush(stage));
264:   PetscCall(DMGetLocalVector(dm, &X));
265:   PetscCall(DMGetLocalVector(dm, &F));
266:   for (i = 0; i < its; ++i) PetscCall(DMPlexSNESComputeResidualFEM(dm, X, F, NULL));
267:   PetscCall(DMRestoreLocalVector(dm, &X));
268:   PetscCall(DMRestoreLocalVector(dm, &F));
269:   PetscCall(PetscLogStagePop());
270:   if (PetscDefined(USE_LOG)) {
271:     const char        *title = "DMPlex Residual Integration";
272:     PetscEventPerfInfo eventInfo;
273:     PetscReal          flopRate, cellRate;
274:     PetscInt           cStart, cEnd, Nf, N;
275:     PetscLogEvent      event;

277:     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
278:     PetscCall(DMGetNumFields(dm, &Nf));
279:     PetscCall(PetscLogEventGetId("DMPlexResidualFE", &event));
280:     PetscCall(PetscLogEventGetPerfInfo(stage, event, &eventInfo));
281:     N        = (cEnd - cStart) * Nf * eventInfo.count;
282:     flopRate = eventInfo.time != 0.0 ? eventInfo.flops / eventInfo.time : 0.0;
283:     cellRate = eventInfo.time != 0.0 ? N / eventInfo.time : 0.0;
284:     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)));
285:   }
286:   PetscFunctionReturn(PETSC_SUCCESS);
287: }

289: int main(int argc, char **argv)
290: {
291:   DM          dm;
292:   AppCtx      ctx;
293:   PetscMPIInt size;

295:   PetscFunctionBeginUser;
296:   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
297:   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
298:   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only.");
299:   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
300:   PetscCall(PetscLogDefaultBegin());
301:   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
302:   PetscCall(DMSetType(dm, DMPLEX));
303:   PetscCall(DMSetFromOptions(dm));
304:   PetscCall(PetscObjectSetName((PetscObject)dm, "Mesh"));
305:   PetscCall(PetscObjectViewFromOptions((PetscObject)dm, NULL, "-dm_view"));
306:   PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &ctx));
307:   PetscCall(TestIntegration(dm, ctx.cbs, ctx.its));
308:   PetscCall(TestIntegration2(dm, ctx.cbs, ctx.its));
309:   PetscCall(DMDestroy(&dm));
310:   PetscCall(PetscFinalize());
311:   return 0;
312: }

314: /*TEST
315:   test:
316:     suffix: 0
317:     requires: triangle
318:     args: -dm_view

320:   test:
321:     suffix: 1
322:     requires: triangle
323:     args: -dm_view -potential_petscspace_degree 1

325:   test:
326:     suffix: 2
327:     requires: triangle
328:     args: -dm_view -potential_petscspace_degree 2

330:   test:
331:     suffix: 3
332:     requires: triangle
333:     args: -dm_view -potential_petscspace_degree 3
334: TEST*/