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 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 31d71ae5a4SJacob Faibussowitsch static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *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 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 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 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 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)); 699566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))trig_u, NULL, user, NULL)); 703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 71c4762a1bSJed Brown } 72c4762a1bSJed Brown 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 97d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscContainerUserDestroy_PetscFEGeom(void *ctx) 98d71ae5a4SJacob Faibussowitsch { 99c4762a1bSJed Brown PetscFEGeom *geom = (PetscFEGeom *)ctx; 100c4762a1bSJed Brown 101c4762a1bSJed Brown PetscFunctionBegin; 1029566063dSJacob Faibussowitsch PetscCall(PetscFEGeomDestroy(&geom)); 1033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 104c4762a1bSJed Brown } 105c4762a1bSJed Brown 106d71ae5a4SJacob Faibussowitsch PetscErrorCode CellRangeGetFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) 107d71ae5a4SJacob Faibussowitsch { 108c4762a1bSJed Brown char composeStr[33] = {0}; 109c4762a1bSJed Brown PetscObjectId id; 110c4762a1bSJed Brown PetscContainer container; 111c4762a1bSJed Brown 112c4762a1bSJed Brown PetscFunctionBegin; 1139566063dSJacob Faibussowitsch PetscCall(PetscObjectGetId((PetscObject)quad, &id)); 11463a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(composeStr, 32, "CellRangeGetFEGeom_%" PetscInt64_FMT "\n", id)); 1159566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)cellIS, composeStr, (PetscObject *)&container)); 116c4762a1bSJed Brown if (container) { 1179566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(container, (void **)geom)); 118c4762a1bSJed Brown } else { 1199566063dSJacob Faibussowitsch PetscCall(DMFieldCreateFEGeom(coordField, cellIS, quad, faceData, geom)); 12003e76207SPierre Jolivet PetscCall(PetscObjectContainerCompose((PetscObject)cellIS, composeStr, *geom, PetscContainerUserDestroy_PetscFEGeom)); 121c4762a1bSJed Brown } 1223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown 125d71ae5a4SJacob Faibussowitsch PetscErrorCode CellRangeRestoreFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) 126d71ae5a4SJacob Faibussowitsch { 127c4762a1bSJed Brown PetscFunctionBegin; 128c4762a1bSJed Brown *geom = NULL; 1293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 130c4762a1bSJed Brown } 131c4762a1bSJed Brown 132d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms) 133d71ae5a4SJacob Faibussowitsch { 134c4762a1bSJed Brown DMField coordField; 135c4762a1bSJed Brown PetscInt Nf, f, maxDegree; 136c4762a1bSJed Brown 137c4762a1bSJed Brown PetscFunctionBeginUser; 138c4762a1bSJed Brown *affineQuad = NULL; 139c4762a1bSJed Brown *affineGeom = NULL; 140c4762a1bSJed Brown *quads = NULL; 141c4762a1bSJed Brown *geoms = NULL; 1429566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 1439566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateField(dm, &coordField)); 1449566063dSJacob Faibussowitsch PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree)); 145c4762a1bSJed Brown if (maxDegree <= 1) { 1469566063dSJacob Faibussowitsch PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, affineQuad)); 1479566063dSJacob Faibussowitsch if (*affineQuad) PetscCall(CellRangeGetFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom)); 148c4762a1bSJed Brown } else { 1499566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(Nf, quads, Nf, geoms)); 150c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 151c4762a1bSJed Brown PetscFE fe; 152c4762a1bSJed Brown 1539566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe)); 1549566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &(*quads)[f])); 1559566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(*quads)[f])); 1569566063dSJacob Faibussowitsch PetscCall(CellRangeGetFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f])); 157c4762a1bSJed Brown } 158c4762a1bSJed Brown } 1593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 160c4762a1bSJed Brown } 161c4762a1bSJed Brown 162d71ae5a4SJacob Faibussowitsch static PetscErrorCode DestroyFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms) 163d71ae5a4SJacob Faibussowitsch { 164c4762a1bSJed Brown DMField coordField; 165c4762a1bSJed Brown PetscInt Nf, f; 166c4762a1bSJed Brown 167c4762a1bSJed Brown PetscFunctionBeginUser; 1689566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 1699566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateField(dm, &coordField)); 170c4762a1bSJed Brown if (*affineQuad) { 1719566063dSJacob Faibussowitsch PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom)); 1729566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(affineQuad)); 173c4762a1bSJed Brown } else { 174c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 1759566063dSJacob Faibussowitsch PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f])); 1769566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&(*quads)[f])); 177c4762a1bSJed Brown } 1789566063dSJacob Faibussowitsch PetscCall(PetscFree2(*quads, *geoms)); 179c4762a1bSJed Brown } 1803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 181c4762a1bSJed Brown } 182c4762a1bSJed Brown 183d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestIntegration(DM dm, PetscInt cbs, PetscInt its) 184d71ae5a4SJacob Faibussowitsch { 185c4762a1bSJed Brown PetscDS ds; 186c4762a1bSJed Brown PetscFEGeom *chunkGeom = NULL; 187c4762a1bSJed Brown PetscQuadrature affineQuad, *quads = NULL; 188c4762a1bSJed Brown PetscFEGeom *affineGeom, **geoms = NULL; 189c4762a1bSJed Brown PetscScalar *u, *elemVec; 190c4762a1bSJed Brown IS cellIS; 191c4762a1bSJed Brown PetscInt depth, cStart, cEnd, cell, chunkSize = cbs, Nch = 0, Nf, f, totDim, i, k; 192c4762a1bSJed Brown PetscLogStage stage; 193c4762a1bSJed Brown PetscLogEvent event; 194c4762a1bSJed Brown 195c4762a1bSJed Brown PetscFunctionBeginUser; 1969566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("PetscFE Residual Integration Test", &stage)); 1979566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("FEIntegRes", PETSCFE_CLASSID, &event)); 1989566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 1999566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 2009566063dSJacob Faibussowitsch PetscCall(DMGetStratumIS(dm, "depth", depth, &cellIS)); 2019566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 20207218a29SMatthew G. Knepley PetscCall(DMGetCellDS(dm, cStart, &ds, NULL)); 2039566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 2049566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 2059566063dSJacob Faibussowitsch PetscCall(CreateFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms)); 2069566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(chunkSize * totDim, &u, chunkSize * totDim, &elemVec)); 207c4762a1bSJed Brown /* Assumptions: 208c4762a1bSJed Brown - Single field 209c4762a1bSJed Brown - No input data 210c4762a1bSJed Brown - No auxiliary data 211c4762a1bSJed Brown - No time-dependence 212c4762a1bSJed Brown */ 213c4762a1bSJed Brown for (i = 0; i < its; ++i) { 214c4762a1bSJed Brown for (cell = cStart; cell < cEnd; cell += chunkSize, ++Nch) { 215c4762a1bSJed Brown const PetscInt cS = cell, cE = PetscMin(cS + chunkSize, cEnd), Ne = cE - cS; 216c4762a1bSJed Brown 2179566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(elemVec, chunkSize * totDim)); 218c4762a1bSJed Brown /* TODO Replace with DMPlexGetCellFields() */ 219c4762a1bSJed Brown for (k = 0; k < chunkSize * totDim; ++k) u[k] = 1.0; 220c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 22106ad1575SMatthew G. Knepley PetscFormKey key; 222c4762a1bSJed Brown PetscFEGeom *geom = affineGeom ? affineGeom : geoms[f]; 223c4762a1bSJed Brown /* PetscQuadrature quad = affineQuad ? affineQuad : quads[f]; */ 224c4762a1bSJed Brown 2259371c9d4SSatish Balay key.label = NULL; 2269371c9d4SSatish Balay key.value = 0; 2279371c9d4SSatish Balay key.field = f; 2289371c9d4SSatish Balay key.part = 0; 2299566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetChunk(geom, cS, cE, &chunkGeom)); 2309566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(event, 0, 0, 0, 0)); 2319566063dSJacob Faibussowitsch PetscCall(PetscFEIntegrateResidual(ds, key, Ne, chunkGeom, u, NULL, NULL, NULL, 0.0, elemVec)); 2329566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(event, 0, 0, 0, 0)); 233c4762a1bSJed Brown } 234c4762a1bSJed Brown } 235c4762a1bSJed Brown } 2369566063dSJacob Faibussowitsch PetscCall(PetscFEGeomRestoreChunk(affineGeom, cStart, cEnd, &chunkGeom)); 2379566063dSJacob Faibussowitsch PetscCall(DestroyFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms)); 2389566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cellIS)); 2399566063dSJacob Faibussowitsch PetscCall(PetscFree2(u, elemVec)); 2409566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 2412611ad71SToby Isaac if (PetscDefined(USE_LOG)) { 242c4762a1bSJed Brown const char *title = "Petsc FE Residual Integration"; 243c4762a1bSJed Brown PetscEventPerfInfo eventInfo; 244c4762a1bSJed Brown PetscInt N = (cEnd - cStart) * Nf * its; 245c4762a1bSJed Brown PetscReal flopRate, cellRate; 246c4762a1bSJed Brown 2479566063dSJacob Faibussowitsch PetscCall(PetscLogEventGetPerfInfo(stage, event, &eventInfo)); 248c4762a1bSJed Brown flopRate = eventInfo.time != 0.0 ? eventInfo.flops / eventInfo.time : 0.0; 249c4762a1bSJed Brown cellRate = eventInfo.time != 0.0 ? N / eventInfo.time : 0.0; 25063a3b9bcSJacob 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))); 251c4762a1bSJed Brown } 2523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 253c4762a1bSJed Brown } 254c4762a1bSJed Brown 255d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestIntegration2(DM dm, PetscInt cbs, PetscInt its) 256d71ae5a4SJacob Faibussowitsch { 257c4762a1bSJed Brown Vec X, F; 258c4762a1bSJed Brown PetscLogStage stage; 259c4762a1bSJed Brown PetscInt i; 260c4762a1bSJed Brown 261c4762a1bSJed Brown PetscFunctionBeginUser; 2629566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("DMPlex Residual Integration Test", &stage)); 2639566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 2649566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &X)); 2659566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &F)); 26648a46eb9SPierre Jolivet for (i = 0; i < its; ++i) PetscCall(DMPlexSNESComputeResidualFEM(dm, X, F, NULL)); 2679566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &X)); 2689566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &F)); 2699566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 2702611ad71SToby Isaac if (PetscDefined(USE_LOG)) { 271c4762a1bSJed Brown const char *title = "DMPlex Residual Integration"; 272c4762a1bSJed Brown PetscEventPerfInfo eventInfo; 273c4762a1bSJed Brown PetscReal flopRate, cellRate; 274c4762a1bSJed Brown PetscInt cStart, cEnd, Nf, N; 275956f8c0dSBarry Smith PetscLogEvent event; 276c4762a1bSJed Brown 2779566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 2789566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &Nf)); 2799566063dSJacob Faibussowitsch PetscCall(PetscLogEventGetId("DMPlexResidualFE", &event)); 2809566063dSJacob Faibussowitsch PetscCall(PetscLogEventGetPerfInfo(stage, event, &eventInfo)); 281c4762a1bSJed Brown N = (cEnd - cStart) * Nf * eventInfo.count; 282c4762a1bSJed Brown flopRate = eventInfo.time != 0.0 ? eventInfo.flops / eventInfo.time : 0.0; 283c4762a1bSJed Brown cellRate = eventInfo.time != 0.0 ? N / eventInfo.time : 0.0; 28463a3b9bcSJacob 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))); 285c4762a1bSJed Brown } 2863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 287c4762a1bSJed Brown } 288c4762a1bSJed Brown 289d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 290d71ae5a4SJacob Faibussowitsch { 291c4762a1bSJed Brown DM dm; 292c4762a1bSJed Brown AppCtx ctx; 293c4762a1bSJed Brown PetscMPIInt size; 294c4762a1bSJed Brown 295327415f7SBarry Smith PetscFunctionBeginUser; 2969566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 2979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 298*bd158744SPierre Jolivet PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only."); 2999566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 3009566063dSJacob Faibussowitsch PetscCall(PetscLogDefaultBegin()); 3019566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &dm)); 3029566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 3039566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 3049566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dm, "Mesh")); 3059566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)dm, NULL, "-dm_view")); 3069566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &ctx)); 3079566063dSJacob Faibussowitsch PetscCall(TestIntegration(dm, ctx.cbs, ctx.its)); 3089566063dSJacob Faibussowitsch PetscCall(TestIntegration2(dm, ctx.cbs, ctx.its)); 3099566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 3109566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 311b122ec5aSJacob Faibussowitsch return 0; 312c4762a1bSJed Brown } 313c4762a1bSJed Brown 314c4762a1bSJed Brown /*TEST 315c4762a1bSJed Brown test: 316c4762a1bSJed Brown suffix: 0 317c4762a1bSJed Brown requires: triangle 318c4762a1bSJed Brown args: -dm_view 319c4762a1bSJed Brown 320c4762a1bSJed Brown test: 321c4762a1bSJed Brown suffix: 1 322c4762a1bSJed Brown requires: triangle 323c4762a1bSJed Brown args: -dm_view -potential_petscspace_degree 1 324c4762a1bSJed Brown 325c4762a1bSJed Brown test: 326c4762a1bSJed Brown suffix: 2 327c4762a1bSJed Brown requires: triangle 328c4762a1bSJed Brown args: -dm_view -potential_petscspace_degree 2 329c4762a1bSJed Brown 330c4762a1bSJed Brown test: 331c4762a1bSJed Brown suffix: 3 332c4762a1bSJed Brown requires: triangle 333c4762a1bSJed Brown args: -dm_view -potential_petscspace_degree 3 334c4762a1bSJed Brown TEST*/ 335