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