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 14c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 15c4762a1bSJed Brown { 16c4762a1bSJed Brown PetscErrorCode ierr; 17c4762a1bSJed Brown 18c4762a1bSJed Brown PetscFunctionBeginUser; 19c4762a1bSJed Brown options->dim = 2; 20c4762a1bSJed Brown options->simplex = PETSC_TRUE; 21c4762a1bSJed Brown options->its = 1; 22c4762a1bSJed Brown options->cbs = 8; 23c4762a1bSJed Brown 24c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "FE Integration Performance Options", "PETSCFE");CHKERRQ(ierr); 255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-dim", "The topological dimension", "ex1.c", options->dim, &options->dim, NULL)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-simplex", "Simplex or hex cells", "ex1.c", options->simplex, &options->simplex, NULL)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-its", "The number of replications for timing", "ex1.c", options->its, &options->its, NULL)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-cbs", "The number of cells in an integration block", "ex1.c", options->cbs, &options->cbs, NULL)); 291e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 30c4762a1bSJed Brown PetscFunctionReturn(0); 31c4762a1bSJed Brown } 32c4762a1bSJed Brown 33c4762a1bSJed Brown static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 34c4762a1bSJed Brown { 35c4762a1bSJed Brown PetscInt d; 36c4762a1bSJed Brown *u = 0.0; 37c4762a1bSJed Brown for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0*PETSC_PI*x[d]); 38c4762a1bSJed Brown return 0; 39c4762a1bSJed Brown } 40c4762a1bSJed Brown 41c4762a1bSJed Brown static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 42c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 43c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 44c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 45c4762a1bSJed Brown { 46c4762a1bSJed Brown PetscInt d; 47c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[0] += -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]); 48c4762a1bSJed Brown } 49c4762a1bSJed Brown 50c4762a1bSJed Brown static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 51c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 52c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 53c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 54c4762a1bSJed Brown { 55c4762a1bSJed Brown PetscInt d; 56c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 57c4762a1bSJed Brown } 58c4762a1bSJed Brown 59c4762a1bSJed Brown static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 60c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 61c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 62c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 63c4762a1bSJed Brown { 64c4762a1bSJed Brown PetscInt d; 65c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 66c4762a1bSJed Brown } 67c4762a1bSJed Brown 68c4762a1bSJed Brown static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 69c4762a1bSJed Brown { 70c4762a1bSJed Brown PetscDS prob; 7145480ffeSMatthew G. Knepley DMLabel label; 72c4762a1bSJed Brown const PetscInt id = 1; 73c4762a1bSJed Brown 74c4762a1bSJed Brown PetscFunctionBeginUser; 755f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &prob)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(prob, 0, f0_trig_u, f1_u)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolution(prob, 0, trig_u, user)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "marker", &label)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) trig_u, NULL, user, NULL)); 81c4762a1bSJed Brown PetscFunctionReturn(0); 82c4762a1bSJed Brown } 83c4762a1bSJed Brown 84c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 85c4762a1bSJed Brown { 86c4762a1bSJed Brown DM cdm = dm; 87c4762a1bSJed Brown PetscFE fe; 88c4762a1bSJed Brown char prefix[PETSC_MAX_PATH_LEN]; 89c4762a1bSJed Brown 90c4762a1bSJed Brown PetscFunctionBeginUser; 91c4762a1bSJed Brown /* Create finite element */ 925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(PetscObjectComm((PetscObject) dm), user->dim, 1, user->simplex, name ? prefix : NULL, -1, &fe)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe, name)); 95c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 965f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dm)); 985f80ce2aSJacob Faibussowitsch CHKERRQ((*setup)(dm, user)); 99c4762a1bSJed Brown while (cdm) { 1005f80ce2aSJacob Faibussowitsch CHKERRQ(DMCopyDisc(dm,cdm)); 101c4762a1bSJed Brown /* TODO: Check whether the boundary of coarse meshes is marked */ 1025f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoarseDM(cdm, &cdm)); 103c4762a1bSJed Brown } 1045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe)); 105c4762a1bSJed Brown PetscFunctionReturn(0); 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 108c4762a1bSJed Brown static PetscErrorCode PetscContainerUserDestroy_PetscFEGeom(void *ctx) 109c4762a1bSJed Brown { 110c4762a1bSJed Brown PetscFEGeom *geom = (PetscFEGeom *) ctx; 111c4762a1bSJed Brown 112c4762a1bSJed Brown PetscFunctionBegin; 1135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEGeomDestroy(&geom)); 114c4762a1bSJed Brown PetscFunctionReturn(0); 115c4762a1bSJed Brown } 116c4762a1bSJed Brown 117c4762a1bSJed Brown PetscErrorCode CellRangeGetFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) 118c4762a1bSJed Brown { 119c4762a1bSJed Brown char composeStr[33] = {0}; 120c4762a1bSJed Brown PetscObjectId id; 121c4762a1bSJed Brown PetscContainer container; 122c4762a1bSJed Brown 123c4762a1bSJed Brown PetscFunctionBegin; 1245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetId((PetscObject) quad, &id)); 1255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(composeStr, 32, "CellRangeGetFEGeom_%x\n", id)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectQuery((PetscObject) cellIS, composeStr, (PetscObject *) &container)); 127c4762a1bSJed Brown if (container) { 1285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerGetPointer(container, (void **) geom)); 129c4762a1bSJed Brown } else { 1305f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldCreateFEGeom(coordField, cellIS, quad, faceData, geom)); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerCreate(PETSC_COMM_SELF, &container)); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerSetPointer(container, (void *) *geom)); 1335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerSetUserDestroy(container, PetscContainerUserDestroy_PetscFEGeom)); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectCompose((PetscObject) cellIS, composeStr, (PetscObject) container)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerDestroy(&container)); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown PetscFunctionReturn(0); 138c4762a1bSJed Brown } 139c4762a1bSJed Brown 140c4762a1bSJed Brown PetscErrorCode CellRangeRestoreFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) 141c4762a1bSJed Brown { 142c4762a1bSJed Brown PetscFunctionBegin; 143c4762a1bSJed Brown *geom = NULL; 144c4762a1bSJed Brown PetscFunctionReturn(0); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147c4762a1bSJed Brown static PetscErrorCode CreateFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms) 148c4762a1bSJed Brown { 149c4762a1bSJed Brown DMField coordField; 150c4762a1bSJed Brown PetscInt Nf, f, maxDegree; 151c4762a1bSJed Brown 152c4762a1bSJed Brown PetscFunctionBeginUser; 153c4762a1bSJed Brown *affineQuad = NULL; 154c4762a1bSJed Brown *affineGeom = NULL; 155c4762a1bSJed Brown *quads = NULL; 156c4762a1bSJed Brown *geoms = NULL; 1575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetNumFields(ds, &Nf)); 1585f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateField(dm, &coordField)); 1595f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree)); 160c4762a1bSJed Brown if (maxDegree <= 1) { 1615f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldCreateDefaultQuadrature(coordField, cellIS, affineQuad)); 1625f80ce2aSJacob Faibussowitsch if (*affineQuad) CHKERRQ(CellRangeGetFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom)); 163c4762a1bSJed Brown } else { 1645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc2(Nf, quads, Nf, geoms)); 165c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 166c4762a1bSJed Brown PetscFE fe; 167c4762a1bSJed Brown 1685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetDiscretization(ds, f, (PetscObject *) &fe)); 1695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEGetQuadrature(fe, &(*quads)[f])); 1705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject) (*quads)[f])); 1715f80ce2aSJacob Faibussowitsch CHKERRQ(CellRangeGetFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f])); 172c4762a1bSJed Brown } 173c4762a1bSJed Brown } 174c4762a1bSJed Brown PetscFunctionReturn(0); 175c4762a1bSJed Brown } 176c4762a1bSJed Brown 177c4762a1bSJed Brown static PetscErrorCode DestroyFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms) 178c4762a1bSJed Brown { 179c4762a1bSJed Brown DMField coordField; 180c4762a1bSJed Brown PetscInt Nf, f; 181c4762a1bSJed Brown 182c4762a1bSJed Brown PetscFunctionBeginUser; 1835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetNumFields(ds, &Nf)); 1845f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateField(dm, &coordField)); 185c4762a1bSJed Brown if (*affineQuad) { 1865f80ce2aSJacob Faibussowitsch CHKERRQ(CellRangeRestoreFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom)); 1875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscQuadratureDestroy(affineQuad)); 188c4762a1bSJed Brown } else { 189c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 1905f80ce2aSJacob Faibussowitsch CHKERRQ(CellRangeRestoreFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f])); 1915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscQuadratureDestroy(&(*quads)[f])); 192c4762a1bSJed Brown } 1935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(*quads, *geoms)); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown PetscFunctionReturn(0); 196c4762a1bSJed Brown } 197c4762a1bSJed Brown 198c4762a1bSJed Brown static PetscErrorCode TestIntegration(DM dm, PetscInt cbs, PetscInt its) 199c4762a1bSJed Brown { 200c4762a1bSJed Brown PetscDS ds; 201c4762a1bSJed Brown PetscFEGeom *chunkGeom = NULL; 202c4762a1bSJed Brown PetscQuadrature affineQuad, *quads = NULL; 203c4762a1bSJed Brown PetscFEGeom *affineGeom, **geoms = NULL; 204c4762a1bSJed Brown PetscScalar *u, *elemVec; 205c4762a1bSJed Brown IS cellIS; 206c4762a1bSJed Brown PetscInt depth, cStart, cEnd, cell, chunkSize = cbs, Nch = 0, Nf, f, totDim, i, k; 207956f8c0dSBarry Smith #if defined(PETSC_USE_LOG) 208c4762a1bSJed Brown PetscLogStage stage; 209c4762a1bSJed Brown PetscLogEvent event; 210956f8c0dSBarry Smith #endif 211c4762a1bSJed Brown 212c4762a1bSJed Brown PetscFunctionBeginUser; 2135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStageRegister("PetscFE Residual Integration Test", &stage)); 2145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventRegister("FEIntegRes", PETSCFE_CLASSID, &event)); 2155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePush(stage)); 2165f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepth(dm, &depth)); 2175f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetStratumIS(dm, "depth", depth, &cellIS)); 2185f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 2195f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCellDS(dm, cStart, &ds)); 2205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetNumFields(ds, &Nf)); 2215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetTotalDimension(ds, &totDim)); 2225f80ce2aSJacob Faibussowitsch CHKERRQ(CreateFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms)); 2235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(chunkSize*totDim, &u, chunkSize*totDim, &elemVec)); 224c4762a1bSJed Brown /* Assumptions: 225c4762a1bSJed Brown - Single field 226c4762a1bSJed Brown - No input data 227c4762a1bSJed Brown - No auxiliary data 228c4762a1bSJed Brown - No time-dependence 229c4762a1bSJed Brown */ 230c4762a1bSJed Brown for (i = 0; i < its; ++i) { 231c4762a1bSJed Brown for (cell = cStart; cell < cEnd; cell += chunkSize, ++Nch) { 232c4762a1bSJed Brown const PetscInt cS = cell, cE = PetscMin(cS + chunkSize, cEnd), Ne = cE - cS; 233c4762a1bSJed Brown 2345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArrayzero(elemVec, chunkSize*totDim)); 235c4762a1bSJed Brown /* TODO Replace with DMPlexGetCellFields() */ 236c4762a1bSJed Brown for (k = 0; k < chunkSize*totDim; ++k) u[k] = 1.0; 237c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 23806ad1575SMatthew G. Knepley PetscFormKey key; 239c4762a1bSJed Brown PetscFEGeom *geom = affineGeom ? affineGeom : geoms[f]; 240c4762a1bSJed Brown /* PetscQuadrature quad = affineQuad ? affineQuad : quads[f]; */ 241c4762a1bSJed Brown 2426528b96dSMatthew G. Knepley key.label = NULL; key.value = 0; key.field = f; 2435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEGeomGetChunk(geom, cS, cE, &chunkGeom)); 2445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventBegin(event,0,0,0,0)); 2455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEIntegrateResidual(ds, key, Ne, chunkGeom, u, NULL, NULL, NULL, 0.0, elemVec)); 2465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventEnd(event,0,0,0,0)); 247c4762a1bSJed Brown } 248c4762a1bSJed Brown } 249c4762a1bSJed Brown } 2505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEGeomRestoreChunk(affineGeom, cStart, cEnd, &chunkGeom)); 2515f80ce2aSJacob Faibussowitsch CHKERRQ(DestroyFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms)); 2525f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&cellIS)); 2535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(u, elemVec)); 2545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePop()); 255956f8c0dSBarry Smith #if defined(PETSC_USE_LOG) 256c4762a1bSJed Brown { 257c4762a1bSJed Brown const char *title = "Petsc FE Residual Integration"; 258c4762a1bSJed Brown PetscEventPerfInfo eventInfo; 259c4762a1bSJed Brown PetscInt N = (cEnd - cStart)*Nf*its; 260c4762a1bSJed Brown PetscReal flopRate, cellRate; 261c4762a1bSJed Brown 2625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventGetPerfInfo(stage, event, &eventInfo)); 263c4762a1bSJed Brown flopRate = eventInfo.time != 0.0 ? eventInfo.flops/eventInfo.time : 0.0; 264c4762a1bSJed Brown cellRate = eventInfo.time != 0.0 ? N/eventInfo.time : 0.0; 2655f80ce2aSJacob Faibussowitsch CHKERRQ(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))); 266c4762a1bSJed Brown } 267956f8c0dSBarry Smith #endif 268c4762a1bSJed Brown PetscFunctionReturn(0); 269c4762a1bSJed Brown } 270c4762a1bSJed Brown 271c4762a1bSJed Brown static PetscErrorCode TestIntegration2(DM dm, PetscInt cbs, PetscInt its) 272c4762a1bSJed Brown { 273c4762a1bSJed Brown Vec X, F; 274956f8c0dSBarry Smith #if defined(PETSC_USE_LOG) 275c4762a1bSJed Brown PetscLogStage stage; 276956f8c0dSBarry Smith #endif 277c4762a1bSJed Brown PetscInt i; 278c4762a1bSJed Brown 279c4762a1bSJed Brown PetscFunctionBeginUser; 2805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStageRegister("DMPlex Residual Integration Test", &stage)); 2815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePush(stage)); 2825f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dm, &X)); 2835f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dm, &F)); 284c4762a1bSJed Brown for (i = 0; i < its; ++i) { 2855f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSNESComputeResidualFEM(dm, X, F, NULL)); 286c4762a1bSJed Brown } 2875f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dm, &X)); 2885f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dm, &F)); 2895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePop()); 290956f8c0dSBarry Smith #if defined(PETSC_USE_LOG) 291c4762a1bSJed Brown { 292c4762a1bSJed Brown const char *title = "DMPlex Residual Integration"; 293c4762a1bSJed Brown PetscEventPerfInfo eventInfo; 294c4762a1bSJed Brown PetscReal flopRate, cellRate; 295c4762a1bSJed Brown PetscInt cStart, cEnd, Nf, N; 296956f8c0dSBarry Smith PetscLogEvent event; 297c4762a1bSJed Brown 2985f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 2995f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetNumFields(dm, &Nf)); 3005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventGetId("DMPlexResidualFE", &event)); 3015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventGetPerfInfo(stage, event, &eventInfo)); 302c4762a1bSJed Brown N = (cEnd - cStart)*Nf*eventInfo.count; 303c4762a1bSJed Brown flopRate = eventInfo.time != 0.0 ? eventInfo.flops/eventInfo.time : 0.0; 304c4762a1bSJed Brown cellRate = eventInfo.time != 0.0 ? N/eventInfo.time : 0.0; 3055f80ce2aSJacob Faibussowitsch CHKERRQ(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))); 306c4762a1bSJed Brown } 307956f8c0dSBarry Smith #endif 308c4762a1bSJed Brown PetscFunctionReturn(0); 309c4762a1bSJed Brown } 310c4762a1bSJed Brown 311c4762a1bSJed Brown int main(int argc, char **argv) 312c4762a1bSJed Brown { 313c4762a1bSJed Brown DM dm; 314c4762a1bSJed Brown AppCtx ctx; 315c4762a1bSJed Brown PetscMPIInt size; 316c4762a1bSJed Brown 317*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, NULL, help)); 3185f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 3192c71b3e2SJacob Faibussowitsch PetscCheckFalse(size > 1,PETSC_COMM_WORLD, PETSC_ERR_SUP, "This is a uniprocessor example only."); 3205f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 3215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogDefaultBegin()); 3225f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(PETSC_COMM_WORLD, &dm)); 3235f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(dm, DMPLEX)); 3245f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(dm)); 3255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) dm, "Mesh")); 3265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectViewFromOptions((PetscObject) dm, NULL, "-dm_view")); 3275f80ce2aSJacob Faibussowitsch CHKERRQ(SetupDiscretization(dm, "potential", SetupPrimalProblem, &ctx)); 3285f80ce2aSJacob Faibussowitsch CHKERRQ(TestIntegration(dm, ctx.cbs, ctx.its)); 3295f80ce2aSJacob Faibussowitsch CHKERRQ(TestIntegration2(dm, ctx.cbs, ctx.its)); 3305f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 331*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 332*b122ec5aSJacob Faibussowitsch return 0; 333c4762a1bSJed Brown } 334c4762a1bSJed Brown 335c4762a1bSJed Brown /*TEST 336c4762a1bSJed Brown test: 337c4762a1bSJed Brown suffix: 0 338c4762a1bSJed Brown requires: triangle 339c4762a1bSJed Brown args: -dm_view 340c4762a1bSJed Brown 341c4762a1bSJed Brown test: 342c4762a1bSJed Brown suffix: 1 343c4762a1bSJed Brown requires: triangle 344c4762a1bSJed Brown args: -dm_view -potential_petscspace_degree 1 345c4762a1bSJed Brown 346c4762a1bSJed Brown test: 347c4762a1bSJed Brown suffix: 2 348c4762a1bSJed Brown requires: triangle 349c4762a1bSJed Brown args: -dm_view -potential_petscspace_degree 2 350c4762a1bSJed Brown 351c4762a1bSJed Brown test: 352c4762a1bSJed Brown suffix: 3 353c4762a1bSJed Brown requires: triangle 354c4762a1bSJed Brown args: -dm_view -potential_petscspace_degree 3 355c4762a1bSJed Brown TEST*/ 356