1dc2eb70dSMatthew G. Knepley static const char help[] = "Tests for injecting basis functions"; 2dc2eb70dSMatthew G. Knepley 3dc2eb70dSMatthew G. Knepley #include <petscdmplex.h> 4dc2eb70dSMatthew G. Knepley #include <petscfe.h> 5dc2eb70dSMatthew G. Knepley #include <petscds.h> 6dc2eb70dSMatthew G. Knepley 7dc2eb70dSMatthew G. Knepley typedef struct { 8dc2eb70dSMatthew G. Knepley PetscInt its; /* Number of replications for timing */ 9dc2eb70dSMatthew G. Knepley } AppCtx; 10dc2eb70dSMatthew G. Knepley 11d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 12d71ae5a4SJacob Faibussowitsch { 13dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 14dc2eb70dSMatthew G. Knepley options->its = 1; 15dc2eb70dSMatthew G. Knepley 16d0609cedSBarry Smith PetscOptionsBegin(comm, "", "FE Injection Options", "PETSCFE"); 179566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-its", "The number of replications for timing", "ex1.c", options->its, &options->its, NULL)); 18d0609cedSBarry Smith PetscOptionsEnd(); 19*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 20dc2eb70dSMatthew G. Knepley } 21dc2eb70dSMatthew G. Knepley 22d71ae5a4SJacob Faibussowitsch static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 23d71ae5a4SJacob Faibussowitsch { 24dc2eb70dSMatthew G. Knepley PetscInt d; 25dc2eb70dSMatthew G. Knepley *u = 0.0; 26dc2eb70dSMatthew G. Knepley for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0 * PETSC_PI * x[d]); 27*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 28dc2eb70dSMatthew G. Knepley } 29dc2eb70dSMatthew G. Knepley 30d71ae5a4SJacob 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[]) 31d71ae5a4SJacob Faibussowitsch { 32dc2eb70dSMatthew G. Knepley PetscInt d; 33dc2eb70dSMatthew G. Knepley for (d = 0; d < dim; ++d) f0[0] += -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]); 34dc2eb70dSMatthew G. Knepley } 35dc2eb70dSMatthew G. Knepley 36d71ae5a4SJacob 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[]) 37d71ae5a4SJacob Faibussowitsch { 38dc2eb70dSMatthew G. Knepley PetscInt d; 39dc2eb70dSMatthew G. Knepley for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 40dc2eb70dSMatthew G. Knepley } 41dc2eb70dSMatthew G. Knepley 42d71ae5a4SJacob 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[]) 43d71ae5a4SJacob Faibussowitsch { 44dc2eb70dSMatthew G. Knepley PetscInt d; 45dc2eb70dSMatthew G. Knepley for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0; 46dc2eb70dSMatthew G. Knepley } 47dc2eb70dSMatthew G. Knepley 48d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 49d71ae5a4SJacob Faibussowitsch { 50dc2eb70dSMatthew G. Knepley PetscDS ds; 51dc2eb70dSMatthew G. Knepley DMLabel label; 52dc2eb70dSMatthew G. Knepley const PetscInt id = 1; 53dc2eb70dSMatthew G. Knepley 54dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 559566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 569566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u)); 579566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 589566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, trig_u, user)); 599566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 609566063dSJacob Faibussowitsch if (label) PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))trig_u, NULL, user, NULL)); 61*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 62dc2eb70dSMatthew G. Knepley } 63dc2eb70dSMatthew G. Knepley 64d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 65d71ae5a4SJacob Faibussowitsch { 66dc2eb70dSMatthew G. Knepley DM cdm = dm; 67dc2eb70dSMatthew G. Knepley PetscFE fe; 68dc2eb70dSMatthew G. Knepley char prefix[PETSC_MAX_PATH_LEN]; 69dc2eb70dSMatthew G. Knepley PetscInt dim; 70dc2eb70dSMatthew G. Knepley 71dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 729566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 739566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 749566063dSJacob Faibussowitsch PetscCall(DMCreateFEDefault(dm, dim, name ? prefix : NULL, -1, &fe)); 759566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, name)); 76dc2eb70dSMatthew G. Knepley /* Set discretization and boundary conditions for each mesh */ 779566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 789566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 799566063dSJacob Faibussowitsch PetscCall((*setup)(dm, user)); 80dc2eb70dSMatthew G. Knepley while (cdm) { 819566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 829566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 83dc2eb70dSMatthew G. Knepley } 849566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 85*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 86dc2eb70dSMatthew G. Knepley } 87dc2eb70dSMatthew G. Knepley 88d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscContainerUserDestroy_PetscFEGeom(void *ctx) 89d71ae5a4SJacob Faibussowitsch { 90dc2eb70dSMatthew G. Knepley PetscFEGeom *geom = (PetscFEGeom *)ctx; 91dc2eb70dSMatthew G. Knepley 92dc2eb70dSMatthew G. Knepley PetscFunctionBegin; 939566063dSJacob Faibussowitsch PetscCall(PetscFEGeomDestroy(&geom)); 94*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 95dc2eb70dSMatthew G. Knepley } 96dc2eb70dSMatthew G. Knepley 97d71ae5a4SJacob Faibussowitsch PetscErrorCode CellRangeGetFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) 98d71ae5a4SJacob Faibussowitsch { 99dc2eb70dSMatthew G. Knepley char composeStr[33] = {0}; 100dc2eb70dSMatthew G. Knepley PetscObjectId id; 101dc2eb70dSMatthew G. Knepley PetscContainer container; 102dc2eb70dSMatthew G. Knepley 103dc2eb70dSMatthew G. Knepley PetscFunctionBegin; 1049566063dSJacob Faibussowitsch PetscCall(PetscObjectGetId((PetscObject)quad, &id)); 10563a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(composeStr, 32, "CellRangeGetFEGeom_%" PetscInt64_FMT "\n", id)); 1069566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)cellIS, composeStr, (PetscObject *)&container)); 107dc2eb70dSMatthew G. Knepley if (container) { 1089566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(container, (void **)geom)); 109dc2eb70dSMatthew G. Knepley } else { 1109566063dSJacob Faibussowitsch PetscCall(DMFieldCreateFEGeom(coordField, cellIS, quad, faceData, geom)); 1119566063dSJacob Faibussowitsch PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container)); 1129566063dSJacob Faibussowitsch PetscCall(PetscContainerSetPointer(container, (void *)*geom)); 1139566063dSJacob Faibussowitsch PetscCall(PetscContainerSetUserDestroy(container, PetscContainerUserDestroy_PetscFEGeom)); 1149566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)cellIS, composeStr, (PetscObject)container)); 1159566063dSJacob Faibussowitsch PetscCall(PetscContainerDestroy(&container)); 116dc2eb70dSMatthew G. Knepley } 117*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 118dc2eb70dSMatthew G. Knepley } 119dc2eb70dSMatthew G. Knepley 120d71ae5a4SJacob Faibussowitsch PetscErrorCode CellRangeRestoreFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) 121d71ae5a4SJacob Faibussowitsch { 122dc2eb70dSMatthew G. Knepley PetscFunctionBegin; 123dc2eb70dSMatthew G. Knepley *geom = NULL; 124*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 125dc2eb70dSMatthew G. Knepley } 126dc2eb70dSMatthew G. Knepley 127d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms) 128d71ae5a4SJacob Faibussowitsch { 129dc2eb70dSMatthew G. Knepley DMField coordField; 130dc2eb70dSMatthew G. Knepley PetscInt Nf, f, maxDegree; 131dc2eb70dSMatthew G. Knepley 132dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 133dc2eb70dSMatthew G. Knepley *affineQuad = NULL; 134dc2eb70dSMatthew G. Knepley *affineGeom = NULL; 135dc2eb70dSMatthew G. Knepley *quads = NULL; 136dc2eb70dSMatthew G. Knepley *geoms = NULL; 1379566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 1389566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateField(dm, &coordField)); 1399566063dSJacob Faibussowitsch PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree)); 140dc2eb70dSMatthew G. Knepley if (maxDegree <= 1) { 1419566063dSJacob Faibussowitsch PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, affineQuad)); 1429566063dSJacob Faibussowitsch if (*affineQuad) PetscCall(CellRangeGetFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom)); 143dc2eb70dSMatthew G. Knepley } else { 1449566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(Nf, quads, Nf, geoms)); 145dc2eb70dSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 146dc2eb70dSMatthew G. Knepley PetscFE fe; 147dc2eb70dSMatthew G. Knepley 1489566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe)); 1499566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &(*quads)[f])); 1509566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(*quads)[f])); 1519566063dSJacob Faibussowitsch PetscCall(CellRangeGetFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f])); 152dc2eb70dSMatthew G. Knepley } 153dc2eb70dSMatthew G. Knepley } 154*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 155dc2eb70dSMatthew G. Knepley } 156dc2eb70dSMatthew G. Knepley 157d71ae5a4SJacob Faibussowitsch static PetscErrorCode DestroyFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms) 158d71ae5a4SJacob Faibussowitsch { 159dc2eb70dSMatthew G. Knepley DMField coordField; 160dc2eb70dSMatthew G. Knepley PetscInt Nf, f; 161dc2eb70dSMatthew G. Knepley 162dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 1639566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 1649566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateField(dm, &coordField)); 165dc2eb70dSMatthew G. Knepley if (*affineQuad) { 1669566063dSJacob Faibussowitsch PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom)); 1679566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(affineQuad)); 168dc2eb70dSMatthew G. Knepley } else { 169dc2eb70dSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 1709566063dSJacob Faibussowitsch PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f])); 1719566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&(*quads)[f])); 172dc2eb70dSMatthew G. Knepley } 1739566063dSJacob Faibussowitsch PetscCall(PetscFree2(*quads, *geoms)); 174dc2eb70dSMatthew G. Knepley } 175*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 176dc2eb70dSMatthew G. Knepley } 177dc2eb70dSMatthew G. Knepley 178d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestEvaluation(DM dm) 179d71ae5a4SJacob Faibussowitsch { 180dc2eb70dSMatthew G. Knepley PetscFE fe; 181dc2eb70dSMatthew G. Knepley PetscSpace sp; 182dc2eb70dSMatthew G. Knepley PetscReal *points; 183dc2eb70dSMatthew G. Knepley PetscReal *B, *D, *H; 184dc2eb70dSMatthew G. Knepley PetscInt dim, Nb, b, Nc, c, Np, p; 185dc2eb70dSMatthew G. Knepley 186dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 1879566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1889566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, 0, NULL, (PetscObject *)&fe)); 189dc2eb70dSMatthew G. Knepley Np = 6; 1909566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Np * dim, &points)); 191dc2eb70dSMatthew G. Knepley if (dim == 3) { 1929371c9d4SSatish Balay points[0] = -1.0; 1939371c9d4SSatish Balay points[1] = -1.0; 1949371c9d4SSatish Balay points[2] = -1.0; 1959371c9d4SSatish Balay points[3] = 1.0; 1969371c9d4SSatish Balay points[4] = -1.0; 1979371c9d4SSatish Balay points[5] = -1.0; 1989371c9d4SSatish Balay points[6] = -1.0; 1999371c9d4SSatish Balay points[7] = 1.0; 2009371c9d4SSatish Balay points[8] = -1.0; 2019371c9d4SSatish Balay points[9] = -1.0; 2029371c9d4SSatish Balay points[10] = -1.0; 2039371c9d4SSatish Balay points[11] = 1.0; 2049371c9d4SSatish Balay points[12] = 1.0; 2059371c9d4SSatish Balay points[13] = -1.0; 2069371c9d4SSatish Balay points[14] = 1.0; 2079371c9d4SSatish Balay points[15] = -1.0; 2089371c9d4SSatish Balay points[16] = 1.0; 2099371c9d4SSatish Balay points[17] = 1.0; 210dc2eb70dSMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only for 3D right now"); 2119566063dSJacob Faibussowitsch PetscCall(PetscFEGetBasisSpace(fe, &sp)); 2129566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDimension(sp, &Nb)); 2139566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetNumComponents(sp, &Nc)); 2149566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, Np * Nb * Nc, MPIU_REAL, &B)); 2159566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, Np * Nb * Nc * dim, MPIU_REAL, &D)); 2169566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, Np * Nb * Nc * dim * dim, MPIU_REAL, &H)); 2179566063dSJacob Faibussowitsch PetscCall(PetscSpaceEvaluate(sp, Np, points, B, NULL, NULL /*D, H*/)); 218dc2eb70dSMatthew G. Knepley for (p = 0; p < Np; ++p) { 2199566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Point %" PetscInt_FMT "\n", p)); 220dc2eb70dSMatthew G. Knepley for (b = 0; b < Nb; ++b) { 2219566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "B[%" PetscInt_FMT "]:", b)); 2229566063dSJacob Faibussowitsch for (c = 0; c < Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)B[(p * Nb + b) * Nc + c])); 2239566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 224dc2eb70dSMatthew G. Knepley #if 0 225dc2eb70dSMatthew G. Knepley for (c = 0; c < Nc; ++c) { 2269566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " D[%" PetscInt_FMT ",%" PetscInt_FMT "]:", b, c)); 2279566063dSJacob Faibussowitsch for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double) B[((p*Nb+b)*Nc+c)*dim+d)])); 2289566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 229dc2eb70dSMatthew G. Knepley } 230dc2eb70dSMatthew G. Knepley #endif 231dc2eb70dSMatthew G. Knepley } 232dc2eb70dSMatthew G. Knepley } 2339566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, Np * Nb, MPIU_REAL, &B)); 2349566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, Np * Nb * dim, MPIU_REAL, &D)); 2359566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, Np * Nb * dim * dim, MPIU_REAL, &H)); 2369566063dSJacob Faibussowitsch PetscCall(PetscFree(points)); 237*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 238dc2eb70dSMatthew G. Knepley } 239dc2eb70dSMatthew G. Knepley 240d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestIntegration(DM dm, PetscInt cbs, PetscInt its) 241d71ae5a4SJacob Faibussowitsch { 242dc2eb70dSMatthew G. Knepley PetscDS ds; 243dc2eb70dSMatthew G. Knepley PetscFEGeom *chunkGeom = NULL; 244dc2eb70dSMatthew G. Knepley PetscQuadrature affineQuad, *quads = NULL; 245dc2eb70dSMatthew G. Knepley PetscFEGeom *affineGeom, **geoms = NULL; 246dc2eb70dSMatthew G. Knepley PetscScalar *u, *elemVec; 247dc2eb70dSMatthew G. Knepley IS cellIS; 248dc2eb70dSMatthew G. Knepley PetscInt depth, cStart, cEnd, cell, chunkSize = cbs, Nch = 0, Nf, f, totDim, i, k; 249dc2eb70dSMatthew G. Knepley 250dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 2519566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 2529566063dSJacob Faibussowitsch PetscCall(DMGetStratumIS(dm, "depth", depth, &cellIS)); 2539566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 2549566063dSJacob Faibussowitsch PetscCall(DMGetCellDS(dm, cStart, &ds)); 2559566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 2569566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 2579566063dSJacob Faibussowitsch PetscCall(CreateFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms)); 2589566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(chunkSize * totDim, &u, chunkSize * totDim, &elemVec)); 259dc2eb70dSMatthew G. Knepley /* Assumptions: 260dc2eb70dSMatthew G. Knepley - Single field 261dc2eb70dSMatthew G. Knepley - No input data 262dc2eb70dSMatthew G. Knepley - No auxiliary data 263dc2eb70dSMatthew G. Knepley - No time-dependence 264dc2eb70dSMatthew G. Knepley */ 265dc2eb70dSMatthew G. Knepley for (i = 0; i < its; ++i) { 266dc2eb70dSMatthew G. Knepley for (cell = cStart; cell < cEnd; cell += chunkSize, ++Nch) { 267dc2eb70dSMatthew G. Knepley const PetscInt cS = cell, cE = PetscMin(cS + chunkSize, cEnd), Ne = cE - cS; 268dc2eb70dSMatthew G. Knepley 2699566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(elemVec, chunkSize * totDim)); 270dc2eb70dSMatthew G. Knepley /* TODO Replace with DMPlexGetCellFields() */ 271dc2eb70dSMatthew G. Knepley for (k = 0; k < chunkSize * totDim; ++k) u[k] = 1.0; 272dc2eb70dSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 273dc2eb70dSMatthew G. Knepley PetscFormKey key; 274dc2eb70dSMatthew G. Knepley PetscFEGeom *geom = affineGeom ? affineGeom : geoms[f]; 275dc2eb70dSMatthew G. Knepley /* PetscQuadrature quad = affineQuad ? affineQuad : quads[f]; */ 276dc2eb70dSMatthew G. Knepley 2779371c9d4SSatish Balay key.label = NULL; 2789371c9d4SSatish Balay key.value = 0; 2799371c9d4SSatish Balay key.field = f; 2809371c9d4SSatish Balay key.part = 0; 2819566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetChunk(geom, cS, cE, &chunkGeom)); 2829566063dSJacob Faibussowitsch PetscCall(PetscFEIntegrateResidual(ds, key, Ne, chunkGeom, u, NULL, NULL, NULL, 0.0, elemVec)); 283dc2eb70dSMatthew G. Knepley } 284dc2eb70dSMatthew G. Knepley } 285dc2eb70dSMatthew G. Knepley } 2869566063dSJacob Faibussowitsch PetscCall(PetscFEGeomRestoreChunk(affineGeom, cStart, cEnd, &chunkGeom)); 2879566063dSJacob Faibussowitsch PetscCall(DestroyFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms)); 2889566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cellIS)); 2899566063dSJacob Faibussowitsch PetscCall(PetscFree2(u, elemVec)); 290*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 291dc2eb70dSMatthew G. Knepley } 292dc2eb70dSMatthew G. Knepley 293d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestUnisolvence(DM dm) 294d71ae5a4SJacob Faibussowitsch { 295dc2eb70dSMatthew G. Knepley Mat M; 296dc2eb70dSMatthew G. Knepley Vec v; 297dc2eb70dSMatthew G. Knepley 298dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 2999566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &v)); 3009566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &v)); 3019566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(dm, dm, &M)); 3029566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M, NULL, "-mass_view")); 3039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M)); 304*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 305dc2eb70dSMatthew G. Knepley } 306dc2eb70dSMatthew G. Knepley 307d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 308d71ae5a4SJacob Faibussowitsch { 309dc2eb70dSMatthew G. Knepley DM dm; 310dc2eb70dSMatthew G. Knepley AppCtx ctx; 311dc2eb70dSMatthew G. Knepley PetscMPIInt size; 312dc2eb70dSMatthew G. Knepley 313327415f7SBarry Smith PetscFunctionBeginUser; 3149566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 3159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 316dc2eb70dSMatthew G. Knepley PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_SUP, "This is a uniprocessor example only."); 3179566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 3189566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &dm)); 3199566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 3209566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 3219566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dm, "Mesh")); 3229566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)dm, NULL, "-dm_view")); 3239566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, "field", SetupPrimalProblem, &ctx)); 3249566063dSJacob Faibussowitsch PetscCall(TestEvaluation(dm)); 3259566063dSJacob Faibussowitsch PetscCall(TestIntegration(dm, 1, ctx.its)); 3269566063dSJacob Faibussowitsch PetscCall(TestUnisolvence(dm)); 3279566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 3289566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 329b122ec5aSJacob Faibussowitsch return 0; 330dc2eb70dSMatthew G. Knepley } 331dc2eb70dSMatthew G. Knepley 332dc2eb70dSMatthew G. Knepley /*TEST 333dc2eb70dSMatthew G. Knepley test: 334dc2eb70dSMatthew G. Knepley suffix: 0 335dc2eb70dSMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -field_petscspace_degree 0 336dc2eb70dSMatthew G. Knepley 337dc2eb70dSMatthew G. Knepley test: 338dc2eb70dSMatthew G. Knepley suffix: 2 339dc2eb70dSMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism \ 340dc2eb70dSMatthew G. Knepley -field_petscspace_type sum \ 341dc2eb70dSMatthew G. Knepley -field_petscspace_variables 3 \ 342dc2eb70dSMatthew G. Knepley -field_petscspace_components 3 \ 343dc2eb70dSMatthew G. Knepley -field_petscspace_sum_spaces 2 \ 344dc2eb70dSMatthew G. Knepley -field_petscspace_sum_concatenate false \ 345dc2eb70dSMatthew G. Knepley -field_sumcomp_0_petscspace_variables 3 \ 346dc2eb70dSMatthew G. Knepley -field_sumcomp_0_petscspace_components 3 \ 347dc2eb70dSMatthew G. Knepley -field_sumcomp_0_petscspace_degree 1 \ 348dc2eb70dSMatthew G. Knepley -field_sumcomp_1_petscspace_variables 3 \ 349dc2eb70dSMatthew G. Knepley -field_sumcomp_1_petscspace_components 3 \ 350dc2eb70dSMatthew G. Knepley -field_sumcomp_1_petscspace_type wxy \ 351dc2eb70dSMatthew G. Knepley -field_petscdualspace_form_degree 0 \ 352dc2eb70dSMatthew G. Knepley -field_petscdualspace_order 1 \ 353dc2eb70dSMatthew G. Knepley -field_petscdualspace_components 3 354dc2eb70dSMatthew G. Knepley 355dc2eb70dSMatthew G. Knepley TEST*/ 356