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