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
ProcessOptions(MPI_Comm comm,AppCtx * options)14d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
15d71ae5a4SJacob Faibussowitsch {
16c4762a1bSJed Brown PetscFunctionBeginUser;
17c4762a1bSJed Brown options->dim = 2;
18c4762a1bSJed Brown options->simplex = PETSC_TRUE;
19c4762a1bSJed Brown options->its = 1;
20c4762a1bSJed Brown options->cbs = 8;
21c4762a1bSJed Brown
22d0609cedSBarry Smith PetscOptionsBegin(comm, "", "FE Integration Performance Options", "PETSCFE");
239566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-dim", "The topological dimension", "ex1.c", options->dim, &options->dim, NULL));
249566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-simplex", "Simplex or hex cells", "ex1.c", options->simplex, &options->simplex, NULL));
259566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-its", "The number of replications for timing", "ex1.c", options->its, &options->its, NULL));
269566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-cbs", "The number of cells in an integration block", "ex1.c", options->cbs, &options->cbs, NULL));
27d0609cedSBarry Smith PetscOptionsEnd();
283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
29c4762a1bSJed Brown }
30c4762a1bSJed Brown
trig_u(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)31*2a8381b2SBarry Smith static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
32d71ae5a4SJacob Faibussowitsch {
33c4762a1bSJed Brown PetscInt d;
34c4762a1bSJed Brown *u = 0.0;
35c4762a1bSJed Brown for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0 * PETSC_PI * x[d]);
363ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
37c4762a1bSJed Brown }
38c4762a1bSJed Brown
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[])39d71ae5a4SJacob Faibussowitsch 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[])
40d71ae5a4SJacob Faibussowitsch {
41c4762a1bSJed Brown PetscInt d;
42c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[0] += -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]);
43c4762a1bSJed Brown }
44c4762a1bSJed Brown
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[])45d71ae5a4SJacob Faibussowitsch 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[])
46d71ae5a4SJacob Faibussowitsch {
47c4762a1bSJed Brown PetscInt d;
48c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u_x[d];
49c4762a1bSJed Brown }
50c4762a1bSJed Brown
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[])51d71ae5a4SJacob Faibussowitsch 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[])
52d71ae5a4SJacob Faibussowitsch {
53c4762a1bSJed Brown PetscInt d;
54c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
55c4762a1bSJed Brown }
56c4762a1bSJed Brown
SetupPrimalProblem(DM dm,AppCtx * user)57d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
58d71ae5a4SJacob Faibussowitsch {
59c4762a1bSJed Brown PetscDS prob;
6045480ffeSMatthew G. Knepley DMLabel label;
61c4762a1bSJed Brown const PetscInt id = 1;
62c4762a1bSJed Brown
63c4762a1bSJed Brown PetscFunctionBeginUser;
649566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob));
659566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_trig_u, f1_u));
669566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu));
679566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 0, trig_u, user));
689566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label));
6957d50842SBarry Smith PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)trig_u, NULL, user, NULL));
703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
71c4762a1bSJed Brown }
72c4762a1bSJed Brown
SetupDiscretization(DM dm,const char name[],PetscErrorCode (* setup)(DM,AppCtx *),AppCtx * user)73d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
74d71ae5a4SJacob Faibussowitsch {
75c4762a1bSJed Brown DM cdm = dm;
76c4762a1bSJed Brown PetscFE fe;
77c4762a1bSJed Brown char prefix[PETSC_MAX_PATH_LEN];
78c4762a1bSJed Brown
79c4762a1bSJed Brown PetscFunctionBeginUser;
80c4762a1bSJed Brown /* Create finite element */
819566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
829566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), user->dim, 1, user->simplex, name ? prefix : NULL, -1, &fe));
839566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, name));
84c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */
859566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
869566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm));
879566063dSJacob Faibussowitsch PetscCall((*setup)(dm, user));
88c4762a1bSJed Brown while (cdm) {
899566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm));
90c4762a1bSJed Brown /* TODO: Check whether the boundary of coarse meshes is marked */
919566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm));
92c4762a1bSJed Brown }
939566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe));
943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
95c4762a1bSJed Brown }
96c4762a1bSJed Brown
9749abdd8aSBarry Smith /* PetscObjectContainerCompose() compose requires void ** signature on destructor */
PetscFEGeomDestroy_Void(PetscCtxRt ctx)98*2a8381b2SBarry Smith static PetscErrorCode PetscFEGeomDestroy_Void(PetscCtxRt ctx)
99d71ae5a4SJacob Faibussowitsch {
10049abdd8aSBarry Smith return PetscFEGeomDestroy((PetscFEGeom **)ctx);
101c4762a1bSJed Brown }
102c4762a1bSJed Brown
CellRangeGetFEGeom(IS cellIS,DMField coordField,PetscQuadrature quad,PetscFEGeomMode mode,PetscFEGeom ** geom)103ac9d17c7SMatthew G. Knepley PetscErrorCode CellRangeGetFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscFEGeomMode mode, PetscFEGeom **geom)
104d71ae5a4SJacob Faibussowitsch {
105c4762a1bSJed Brown char composeStr[33] = {0};
106c4762a1bSJed Brown PetscObjectId id;
107c4762a1bSJed Brown PetscContainer container;
108c4762a1bSJed Brown
109c4762a1bSJed Brown PetscFunctionBegin;
1109566063dSJacob Faibussowitsch PetscCall(PetscObjectGetId((PetscObject)quad, &id));
11163a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(composeStr, 32, "CellRangeGetFEGeom_%" PetscInt64_FMT "\n", id));
1129566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)cellIS, composeStr, (PetscObject *)&container));
113c4762a1bSJed Brown if (container) {
114*2a8381b2SBarry Smith PetscCall(PetscContainerGetPointer(container, geom));
115c4762a1bSJed Brown } else {
116ac9d17c7SMatthew G. Knepley PetscCall(DMFieldCreateFEGeom(coordField, cellIS, quad, mode, geom));
11749abdd8aSBarry Smith PetscCall(PetscObjectContainerCompose((PetscObject)cellIS, composeStr, *geom, PetscFEGeomDestroy_Void));
118c4762a1bSJed Brown }
1193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
120c4762a1bSJed Brown }
121c4762a1bSJed Brown
CellRangeRestoreFEGeom(IS cellIS,DMField coordField,PetscQuadrature quad,PetscBool faceData,PetscFEGeom ** geom)122d71ae5a4SJacob Faibussowitsch PetscErrorCode CellRangeRestoreFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
123d71ae5a4SJacob Faibussowitsch {
124c4762a1bSJed Brown PetscFunctionBegin;
125c4762a1bSJed Brown *geom = NULL;
1263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
127c4762a1bSJed Brown }
128c4762a1bSJed Brown
CreateFEGeometry(DM dm,PetscDS ds,IS cellIS,PetscQuadrature * affineQuad,PetscFEGeom ** affineGeom,PetscQuadrature ** quads,PetscFEGeom *** geoms)129d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms)
130d71ae5a4SJacob Faibussowitsch {
131c4762a1bSJed Brown DMField coordField;
132c4762a1bSJed Brown PetscInt Nf, f, maxDegree;
133c4762a1bSJed Brown
134c4762a1bSJed Brown PetscFunctionBeginUser;
135c4762a1bSJed Brown *affineQuad = NULL;
136c4762a1bSJed Brown *affineGeom = NULL;
137c4762a1bSJed Brown *quads = NULL;
138c4762a1bSJed Brown *geoms = NULL;
1399566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf));
1409566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateField(dm, &coordField));
1419566063dSJacob Faibussowitsch PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
142c4762a1bSJed Brown if (maxDegree <= 1) {
1439566063dSJacob Faibussowitsch PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, affineQuad));
144ac9d17c7SMatthew G. Knepley if (*affineQuad) PetscCall(CellRangeGetFEGeom(cellIS, coordField, *affineQuad, PETSC_FEGEOM_BASIC, affineGeom));
145c4762a1bSJed Brown } else {
1469566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(Nf, quads, Nf, geoms));
147c4762a1bSJed Brown for (f = 0; f < Nf; ++f) {
148c4762a1bSJed Brown PetscFE fe;
149c4762a1bSJed Brown
1509566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
1519566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &(*quads)[f]));
1529566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(*quads)[f]));
153ac9d17c7SMatthew G. Knepley PetscCall(CellRangeGetFEGeom(cellIS, coordField, (*quads)[f], PETSC_FEGEOM_BASIC, &(*geoms)[f]));
154c4762a1bSJed Brown }
155c4762a1bSJed Brown }
1563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
157c4762a1bSJed Brown }
158c4762a1bSJed Brown
DestroyFEGeometry(DM dm,PetscDS ds,IS cellIS,PetscQuadrature * affineQuad,PetscFEGeom ** affineGeom,PetscQuadrature ** quads,PetscFEGeom *** geoms)159d71ae5a4SJacob Faibussowitsch static PetscErrorCode DestroyFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms)
160d71ae5a4SJacob Faibussowitsch {
161c4762a1bSJed Brown DMField coordField;
162c4762a1bSJed Brown PetscInt Nf, f;
163c4762a1bSJed Brown
164c4762a1bSJed Brown PetscFunctionBeginUser;
1659566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf));
1669566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateField(dm, &coordField));
167c4762a1bSJed Brown if (*affineQuad) {
1689566063dSJacob Faibussowitsch PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom));
1699566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(affineQuad));
170c4762a1bSJed Brown } else {
171c4762a1bSJed Brown for (f = 0; f < Nf; ++f) {
1729566063dSJacob Faibussowitsch PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f]));
1739566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&(*quads)[f]));
174c4762a1bSJed Brown }
1759566063dSJacob Faibussowitsch PetscCall(PetscFree2(*quads, *geoms));
176c4762a1bSJed Brown }
1773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
178c4762a1bSJed Brown }
179c4762a1bSJed Brown
TestIntegration(DM dm,PetscInt cbs,PetscInt its)180d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestIntegration(DM dm, PetscInt cbs, PetscInt its)
181d71ae5a4SJacob Faibussowitsch {
182c4762a1bSJed Brown PetscDS ds;
183c4762a1bSJed Brown PetscFEGeom *chunkGeom = NULL;
184c4762a1bSJed Brown PetscQuadrature affineQuad, *quads = NULL;
185c4762a1bSJed Brown PetscFEGeom *affineGeom, **geoms = NULL;
186c4762a1bSJed Brown PetscScalar *u, *elemVec;
187c4762a1bSJed Brown IS cellIS;
188c4762a1bSJed Brown PetscInt depth, cStart, cEnd, cell, chunkSize = cbs, Nch = 0, Nf, f, totDim, i, k;
189c4762a1bSJed Brown PetscLogStage stage;
190c4762a1bSJed Brown PetscLogEvent event;
191c4762a1bSJed Brown
192c4762a1bSJed Brown PetscFunctionBeginUser;
1939566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("PetscFE Residual Integration Test", &stage));
1949566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("FEIntegRes", PETSCFE_CLASSID, &event));
1959566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage));
1969566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth));
1979566063dSJacob Faibussowitsch PetscCall(DMGetStratumIS(dm, "depth", depth, &cellIS));
1989566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
19907218a29SMatthew G. Knepley PetscCall(DMGetCellDS(dm, cStart, &ds, NULL));
2009566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf));
2019566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim));
2029566063dSJacob Faibussowitsch PetscCall(CreateFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms));
2039566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(chunkSize * totDim, &u, chunkSize * totDim, &elemVec));
204c4762a1bSJed Brown /* Assumptions:
205c4762a1bSJed Brown - Single field
206c4762a1bSJed Brown - No input data
207c4762a1bSJed Brown - No auxiliary data
208c4762a1bSJed Brown - No time-dependence
209c4762a1bSJed Brown */
210c4762a1bSJed Brown for (i = 0; i < its; ++i) {
211c4762a1bSJed Brown for (cell = cStart; cell < cEnd; cell += chunkSize, ++Nch) {
212c4762a1bSJed Brown const PetscInt cS = cell, cE = PetscMin(cS + chunkSize, cEnd), Ne = cE - cS;
213c4762a1bSJed Brown
2149566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(elemVec, chunkSize * totDim));
215c4762a1bSJed Brown /* TODO Replace with DMPlexGetCellFields() */
216c4762a1bSJed Brown for (k = 0; k < chunkSize * totDim; ++k) u[k] = 1.0;
217c4762a1bSJed Brown for (f = 0; f < Nf; ++f) {
21806ad1575SMatthew G. Knepley PetscFormKey key;
219c4762a1bSJed Brown PetscFEGeom *geom = affineGeom ? affineGeom : geoms[f];
220c4762a1bSJed Brown /* PetscQuadrature quad = affineQuad ? affineQuad : quads[f]; */
221c4762a1bSJed Brown
2229371c9d4SSatish Balay key.label = NULL;
2239371c9d4SSatish Balay key.value = 0;
2249371c9d4SSatish Balay key.field = f;
2259371c9d4SSatish Balay key.part = 0;
2269566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetChunk(geom, cS, cE, &chunkGeom));
2279566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(event, 0, 0, 0, 0));
2289566063dSJacob Faibussowitsch PetscCall(PetscFEIntegrateResidual(ds, key, Ne, chunkGeom, u, NULL, NULL, NULL, 0.0, elemVec));
2299566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(event, 0, 0, 0, 0));
230c4762a1bSJed Brown }
231c4762a1bSJed Brown }
232c4762a1bSJed Brown }
2339566063dSJacob Faibussowitsch PetscCall(PetscFEGeomRestoreChunk(affineGeom, cStart, cEnd, &chunkGeom));
2349566063dSJacob Faibussowitsch PetscCall(DestroyFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms));
2359566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cellIS));
2369566063dSJacob Faibussowitsch PetscCall(PetscFree2(u, elemVec));
2379566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop());
2382611ad71SToby Isaac if (PetscDefined(USE_LOG)) {
239f0b74427SPierre Jolivet const char *title = "PETSc FE Residual Integration";
240c4762a1bSJed Brown PetscEventPerfInfo eventInfo;
241c4762a1bSJed Brown PetscInt N = (cEnd - cStart) * Nf * its;
242c4762a1bSJed Brown PetscReal flopRate, cellRate;
243c4762a1bSJed Brown
2449566063dSJacob Faibussowitsch PetscCall(PetscLogEventGetPerfInfo(stage, event, &eventInfo));
245c4762a1bSJed Brown flopRate = eventInfo.time != 0.0 ? eventInfo.flops / eventInfo.time : 0.0;
246c4762a1bSJed Brown cellRate = eventInfo.time != 0.0 ? N / eventInfo.time : 0.0;
24763a3b9bcSJacob 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)));
248c4762a1bSJed Brown }
2493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
250c4762a1bSJed Brown }
251c4762a1bSJed Brown
TestIntegration2(DM dm,PetscInt cbs,PetscInt its)252d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestIntegration2(DM dm, PetscInt cbs, PetscInt its)
253d71ae5a4SJacob Faibussowitsch {
254c4762a1bSJed Brown Vec X, F;
255c4762a1bSJed Brown PetscLogStage stage;
256c4762a1bSJed Brown PetscInt i;
257c4762a1bSJed Brown
258c4762a1bSJed Brown PetscFunctionBeginUser;
2599566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("DMPlex Residual Integration Test", &stage));
2609566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage));
2619566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &X));
2629566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &F));
26348a46eb9SPierre Jolivet for (i = 0; i < its; ++i) PetscCall(DMPlexSNESComputeResidualFEM(dm, X, F, NULL));
2649566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &X));
2659566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &F));
2669566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop());
2672611ad71SToby Isaac if (PetscDefined(USE_LOG)) {
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 }
2833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
284c4762a1bSJed Brown }
285c4762a1bSJed Brown
main(int argc,char ** argv)286d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
287d71ae5a4SJacob Faibussowitsch {
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));
295bd158744SPierre Jolivet PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "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