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 11dc2eb70dSMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 12dc2eb70dSMatthew G. Knepley { 13dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 14dc2eb70dSMatthew G. Knepley options->its = 1; 15dc2eb70dSMatthew G. Knepley 16*d0609cedSBarry 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)); 18*d0609cedSBarry Smith PetscOptionsEnd(); 19dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 20dc2eb70dSMatthew G. Knepley } 21dc2eb70dSMatthew G. Knepley 22dc2eb70dSMatthew G. Knepley static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 23dc2eb70dSMatthew G. Knepley { 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]); 27dc2eb70dSMatthew G. Knepley return 0; 28dc2eb70dSMatthew G. Knepley } 29dc2eb70dSMatthew G. Knepley 30dc2eb70dSMatthew G. Knepley static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 31dc2eb70dSMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 32dc2eb70dSMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 33dc2eb70dSMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 34dc2eb70dSMatthew G. Knepley { 35dc2eb70dSMatthew G. Knepley PetscInt d; 36dc2eb70dSMatthew G. Knepley for (d = 0; d < dim; ++d) f0[0] += -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]); 37dc2eb70dSMatthew G. Knepley } 38dc2eb70dSMatthew G. Knepley 39dc2eb70dSMatthew G. Knepley static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 40dc2eb70dSMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 41dc2eb70dSMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 42dc2eb70dSMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 43dc2eb70dSMatthew G. Knepley { 44dc2eb70dSMatthew G. Knepley PetscInt d; 45dc2eb70dSMatthew G. Knepley for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 46dc2eb70dSMatthew G. Knepley } 47dc2eb70dSMatthew G. Knepley 48dc2eb70dSMatthew G. Knepley static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 49dc2eb70dSMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 50dc2eb70dSMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 51dc2eb70dSMatthew G. Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 52dc2eb70dSMatthew G. Knepley { 53dc2eb70dSMatthew G. Knepley PetscInt d; 54dc2eb70dSMatthew G. Knepley for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 55dc2eb70dSMatthew G. Knepley } 56dc2eb70dSMatthew G. Knepley 57dc2eb70dSMatthew G. Knepley static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 58dc2eb70dSMatthew G. Knepley { 59dc2eb70dSMatthew G. Knepley PetscDS ds; 60dc2eb70dSMatthew G. Knepley DMLabel label; 61dc2eb70dSMatthew G. Knepley const PetscInt id = 1; 62dc2eb70dSMatthew G. Knepley 63dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 649566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 659566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u)); 669566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 679566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, trig_u, user)); 689566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 699566063dSJacob Faibussowitsch if (label) PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) trig_u, NULL, user, NULL)); 70dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 71dc2eb70dSMatthew G. Knepley } 72dc2eb70dSMatthew G. Knepley 73dc2eb70dSMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 74dc2eb70dSMatthew G. Knepley { 75dc2eb70dSMatthew G. Knepley DM cdm = dm; 76dc2eb70dSMatthew G. Knepley PetscFE fe; 77dc2eb70dSMatthew G. Knepley char prefix[PETSC_MAX_PATH_LEN]; 78dc2eb70dSMatthew G. Knepley PetscInt dim; 79dc2eb70dSMatthew G. Knepley 80dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 819566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 829566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 839566063dSJacob Faibussowitsch PetscCall(DMCreateFEDefault(dm, dim, name ? prefix : NULL, -1, &fe)); 849566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe, name)); 85dc2eb70dSMatthew G. Knepley /* Set discretization and boundary conditions for each mesh */ 869566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe)); 879566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 889566063dSJacob Faibussowitsch PetscCall((*setup)(dm, user)); 89dc2eb70dSMatthew G. Knepley while (cdm) { 909566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm,cdm)); 919566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 92dc2eb70dSMatthew G. Knepley } 939566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 94dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 95dc2eb70dSMatthew G. Knepley } 96dc2eb70dSMatthew G. Knepley 97dc2eb70dSMatthew G. Knepley static PetscErrorCode PetscContainerUserDestroy_PetscFEGeom(void *ctx) 98dc2eb70dSMatthew G. Knepley { 99dc2eb70dSMatthew G. Knepley PetscFEGeom *geom = (PetscFEGeom *) ctx; 100dc2eb70dSMatthew G. Knepley 101dc2eb70dSMatthew G. Knepley PetscFunctionBegin; 1029566063dSJacob Faibussowitsch PetscCall(PetscFEGeomDestroy(&geom)); 103dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 104dc2eb70dSMatthew G. Knepley } 105dc2eb70dSMatthew G. Knepley 106dc2eb70dSMatthew G. Knepley PetscErrorCode CellRangeGetFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) 107dc2eb70dSMatthew G. Knepley { 108dc2eb70dSMatthew G. Knepley char composeStr[33] = {0}; 109dc2eb70dSMatthew G. Knepley PetscObjectId id; 110dc2eb70dSMatthew G. Knepley PetscContainer container; 111dc2eb70dSMatthew G. Knepley 112dc2eb70dSMatthew G. Knepley PetscFunctionBegin; 1139566063dSJacob Faibussowitsch PetscCall(PetscObjectGetId((PetscObject) quad, &id)); 1149566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(composeStr, 32, "CellRangeGetFEGeom_%x\n", id)); 1159566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject) cellIS, composeStr, (PetscObject *) &container)); 116dc2eb70dSMatthew G. Knepley if (container) { 1179566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(container, (void **) geom)); 118dc2eb70dSMatthew G. Knepley } else { 1199566063dSJacob Faibussowitsch PetscCall(DMFieldCreateFEGeom(coordField, cellIS, quad, faceData, geom)); 1209566063dSJacob Faibussowitsch PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container)); 1219566063dSJacob Faibussowitsch PetscCall(PetscContainerSetPointer(container, (void *) *geom)); 1229566063dSJacob Faibussowitsch PetscCall(PetscContainerSetUserDestroy(container, PetscContainerUserDestroy_PetscFEGeom)); 1239566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject) cellIS, composeStr, (PetscObject) container)); 1249566063dSJacob Faibussowitsch PetscCall(PetscContainerDestroy(&container)); 125dc2eb70dSMatthew G. Knepley } 126dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 127dc2eb70dSMatthew G. Knepley } 128dc2eb70dSMatthew G. Knepley 129dc2eb70dSMatthew G. Knepley PetscErrorCode CellRangeRestoreFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) 130dc2eb70dSMatthew G. Knepley { 131dc2eb70dSMatthew G. Knepley PetscFunctionBegin; 132dc2eb70dSMatthew G. Knepley *geom = NULL; 133dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 134dc2eb70dSMatthew G. Knepley } 135dc2eb70dSMatthew G. Knepley 136dc2eb70dSMatthew G. Knepley static PetscErrorCode CreateFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms) 137dc2eb70dSMatthew G. Knepley { 138dc2eb70dSMatthew G. Knepley DMField coordField; 139dc2eb70dSMatthew G. Knepley PetscInt Nf, f, maxDegree; 140dc2eb70dSMatthew G. Knepley 141dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 142dc2eb70dSMatthew G. Knepley *affineQuad = NULL; 143dc2eb70dSMatthew G. Knepley *affineGeom = NULL; 144dc2eb70dSMatthew G. Knepley *quads = NULL; 145dc2eb70dSMatthew G. Knepley *geoms = NULL; 1469566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 1479566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateField(dm, &coordField)); 1489566063dSJacob Faibussowitsch PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree)); 149dc2eb70dSMatthew G. Knepley if (maxDegree <= 1) { 1509566063dSJacob Faibussowitsch PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, affineQuad)); 1519566063dSJacob Faibussowitsch if (*affineQuad) PetscCall(CellRangeGetFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom)); 152dc2eb70dSMatthew G. Knepley } else { 1539566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(Nf, quads, Nf, geoms)); 154dc2eb70dSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 155dc2eb70dSMatthew G. Knepley PetscFE fe; 156dc2eb70dSMatthew G. Knepley 1579566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *) &fe)); 1589566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &(*quads)[f])); 1599566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject) (*quads)[f])); 1609566063dSJacob Faibussowitsch PetscCall(CellRangeGetFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f])); 161dc2eb70dSMatthew G. Knepley } 162dc2eb70dSMatthew G. Knepley } 163dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 164dc2eb70dSMatthew G. Knepley } 165dc2eb70dSMatthew G. Knepley 166dc2eb70dSMatthew G. Knepley static PetscErrorCode DestroyFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms) 167dc2eb70dSMatthew G. Knepley { 168dc2eb70dSMatthew G. Knepley DMField coordField; 169dc2eb70dSMatthew G. Knepley PetscInt Nf, f; 170dc2eb70dSMatthew G. Knepley 171dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 1729566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 1739566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateField(dm, &coordField)); 174dc2eb70dSMatthew G. Knepley if (*affineQuad) { 1759566063dSJacob Faibussowitsch PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom)); 1769566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(affineQuad)); 177dc2eb70dSMatthew G. Knepley } else { 178dc2eb70dSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 1799566063dSJacob Faibussowitsch PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f])); 1809566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&(*quads)[f])); 181dc2eb70dSMatthew G. Knepley } 1829566063dSJacob Faibussowitsch PetscCall(PetscFree2(*quads, *geoms)); 183dc2eb70dSMatthew G. Knepley } 184dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 185dc2eb70dSMatthew G. Knepley } 186dc2eb70dSMatthew G. Knepley 187dc2eb70dSMatthew G. Knepley static PetscErrorCode TestEvaluation(DM dm) 188dc2eb70dSMatthew G. Knepley { 189dc2eb70dSMatthew G. Knepley PetscFE fe; 190dc2eb70dSMatthew G. Knepley PetscSpace sp; 191dc2eb70dSMatthew G. Knepley PetscReal *points; 192dc2eb70dSMatthew G. Knepley PetscReal *B, *D, *H; 193dc2eb70dSMatthew G. Knepley PetscInt dim, Nb, b, Nc, c, Np, p; 194dc2eb70dSMatthew G. Knepley 195dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 1969566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1979566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, 0, NULL, (PetscObject *) &fe)); 198dc2eb70dSMatthew G. Knepley Np = 6; 1999566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Np*dim, &points)); 200dc2eb70dSMatthew G. Knepley if (dim == 3) { 201dc2eb70dSMatthew G. Knepley points[0] = -1.0; points[1] = -1.0; points[2] = -1.0; 202dc2eb70dSMatthew G. Knepley points[3] = 1.0; points[4] = -1.0; points[5] = -1.0; 203dc2eb70dSMatthew G. Knepley points[6] = -1.0; points[7] = 1.0; points[8] = -1.0; 204dc2eb70dSMatthew G. Knepley points[9] = -1.0; points[10] = -1.0; points[11] = 1.0; 205dc2eb70dSMatthew G. Knepley points[12] = 1.0; points[13] = -1.0; points[14] = 1.0; 206dc2eb70dSMatthew G. Knepley points[15] = -1.0; points[16] = 1.0; points[17] = 1.0; 207dc2eb70dSMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only for 3D right now"); 2089566063dSJacob Faibussowitsch PetscCall(PetscFEGetBasisSpace(fe, &sp)); 2099566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDimension(sp, &Nb)); 2109566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetNumComponents(sp, &Nc)); 2119566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, Np*Nb*Nc, MPIU_REAL, &B)); 2129566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, Np*Nb*Nc*dim, MPIU_REAL, &D)); 2139566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, Np*Nb*Nc*dim*dim, MPIU_REAL, &H)); 2149566063dSJacob Faibussowitsch PetscCall(PetscSpaceEvaluate(sp, Np, points, B, NULL, NULL /*D, H*/)); 215dc2eb70dSMatthew G. Knepley for (p = 0; p < Np; ++p) { 2169566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Point %" PetscInt_FMT "\n", p)); 217dc2eb70dSMatthew G. Knepley for (b = 0; b < Nb; ++b) { 2189566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "B[%" PetscInt_FMT "]:", b)); 2199566063dSJacob Faibussowitsch for (c = 0; c < Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double) B[(p*Nb+b)*Nc+c])); 2209566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 221dc2eb70dSMatthew G. Knepley #if 0 222dc2eb70dSMatthew G. Knepley for (c = 0; c < Nc; ++c) { 2239566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " D[%" PetscInt_FMT ",%" PetscInt_FMT "]:", b, c)); 2249566063dSJacob Faibussowitsch for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double) B[((p*Nb+b)*Nc+c)*dim+d)])); 2259566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 226dc2eb70dSMatthew G. Knepley } 227dc2eb70dSMatthew G. Knepley #endif 228dc2eb70dSMatthew G. Knepley } 229dc2eb70dSMatthew G. Knepley } 2309566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, Np*Nb, MPIU_REAL, &B)); 2319566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, Np*Nb*dim, MPIU_REAL, &D)); 2329566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, Np*Nb*dim*dim, MPIU_REAL, &H)); 2339566063dSJacob Faibussowitsch PetscCall(PetscFree(points)); 234dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 235dc2eb70dSMatthew G. Knepley } 236dc2eb70dSMatthew G. Knepley 237dc2eb70dSMatthew G. Knepley static PetscErrorCode TestIntegration(DM dm, PetscInt cbs, PetscInt its) 238dc2eb70dSMatthew G. Knepley { 239dc2eb70dSMatthew G. Knepley PetscDS ds; 240dc2eb70dSMatthew G. Knepley PetscFEGeom *chunkGeom = NULL; 241dc2eb70dSMatthew G. Knepley PetscQuadrature affineQuad, *quads = NULL; 242dc2eb70dSMatthew G. Knepley PetscFEGeom *affineGeom, **geoms = NULL; 243dc2eb70dSMatthew G. Knepley PetscScalar *u, *elemVec; 244dc2eb70dSMatthew G. Knepley IS cellIS; 245dc2eb70dSMatthew G. Knepley PetscInt depth, cStart, cEnd, cell, chunkSize = cbs, Nch = 0, Nf, f, totDim, i, k; 246dc2eb70dSMatthew G. Knepley 247dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 2489566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 2499566063dSJacob Faibussowitsch PetscCall(DMGetStratumIS(dm, "depth", depth, &cellIS)); 2509566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 2519566063dSJacob Faibussowitsch PetscCall(DMGetCellDS(dm, cStart, &ds)); 2529566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 2539566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 2549566063dSJacob Faibussowitsch PetscCall(CreateFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms)); 2559566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(chunkSize*totDim, &u, chunkSize*totDim, &elemVec)); 256dc2eb70dSMatthew G. Knepley /* Assumptions: 257dc2eb70dSMatthew G. Knepley - Single field 258dc2eb70dSMatthew G. Knepley - No input data 259dc2eb70dSMatthew G. Knepley - No auxiliary data 260dc2eb70dSMatthew G. Knepley - No time-dependence 261dc2eb70dSMatthew G. Knepley */ 262dc2eb70dSMatthew G. Knepley for (i = 0; i < its; ++i) { 263dc2eb70dSMatthew G. Knepley for (cell = cStart; cell < cEnd; cell += chunkSize, ++Nch) { 264dc2eb70dSMatthew G. Knepley const PetscInt cS = cell, cE = PetscMin(cS + chunkSize, cEnd), Ne = cE - cS; 265dc2eb70dSMatthew G. Knepley 2669566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(elemVec, chunkSize*totDim)); 267dc2eb70dSMatthew G. Knepley /* TODO Replace with DMPlexGetCellFields() */ 268dc2eb70dSMatthew G. Knepley for (k = 0; k < chunkSize*totDim; ++k) u[k] = 1.0; 269dc2eb70dSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 270dc2eb70dSMatthew G. Knepley PetscFormKey key; 271dc2eb70dSMatthew G. Knepley PetscFEGeom *geom = affineGeom ? affineGeom : geoms[f]; 272dc2eb70dSMatthew G. Knepley /* PetscQuadrature quad = affineQuad ? affineQuad : quads[f]; */ 273dc2eb70dSMatthew G. Knepley 274dc2eb70dSMatthew G. Knepley key.label = NULL; key.value = 0; key.field = f; 2759566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetChunk(geom, cS, cE, &chunkGeom)); 2769566063dSJacob Faibussowitsch PetscCall(PetscFEIntegrateResidual(ds, key, Ne, chunkGeom, u, NULL, NULL, NULL, 0.0, elemVec)); 277dc2eb70dSMatthew G. Knepley } 278dc2eb70dSMatthew G. Knepley } 279dc2eb70dSMatthew G. Knepley } 2809566063dSJacob Faibussowitsch PetscCall(PetscFEGeomRestoreChunk(affineGeom, cStart, cEnd, &chunkGeom)); 2819566063dSJacob Faibussowitsch PetscCall(DestroyFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms)); 2829566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cellIS)); 2839566063dSJacob Faibussowitsch PetscCall(PetscFree2(u, elemVec)); 284dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 285dc2eb70dSMatthew G. Knepley } 286dc2eb70dSMatthew G. Knepley 287dc2eb70dSMatthew G. Knepley static PetscErrorCode TestUnisolvence(DM dm) 288dc2eb70dSMatthew G. Knepley { 289dc2eb70dSMatthew G. Knepley Mat M; 290dc2eb70dSMatthew G. Knepley Vec v; 291dc2eb70dSMatthew G. Knepley 292dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 2939566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &v)); 2949566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &v)); 2959566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(dm, dm, &M)); 2969566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M, NULL, "-mass_view")); 2979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M)); 298dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 299dc2eb70dSMatthew G. Knepley } 300dc2eb70dSMatthew G. Knepley 301dc2eb70dSMatthew G. Knepley int main(int argc, char **argv) 302dc2eb70dSMatthew G. Knepley { 303dc2eb70dSMatthew G. Knepley DM dm; 304dc2eb70dSMatthew G. Knepley AppCtx ctx; 305dc2eb70dSMatthew G. Knepley PetscMPIInt size; 306dc2eb70dSMatthew G. Knepley 3079566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 3089566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 309dc2eb70dSMatthew G. Knepley PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_SUP, "This is a uniprocessor example only."); 3109566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 3119566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &dm)); 3129566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 3139566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 3149566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) dm, "Mesh")); 3159566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject) dm, NULL, "-dm_view")); 3169566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, "field", SetupPrimalProblem, &ctx)); 3179566063dSJacob Faibussowitsch PetscCall(TestEvaluation(dm)); 3189566063dSJacob Faibussowitsch PetscCall(TestIntegration(dm, 1, ctx.its)); 3199566063dSJacob Faibussowitsch PetscCall(TestUnisolvence(dm)); 3209566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 3219566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 322b122ec5aSJacob Faibussowitsch return 0; 323dc2eb70dSMatthew G. Knepley } 324dc2eb70dSMatthew G. Knepley 325dc2eb70dSMatthew G. Knepley /*TEST 326dc2eb70dSMatthew G. Knepley test: 327dc2eb70dSMatthew G. Knepley suffix: 0 328dc2eb70dSMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -field_petscspace_degree 0 329dc2eb70dSMatthew G. Knepley 330dc2eb70dSMatthew G. Knepley test: 331dc2eb70dSMatthew G. Knepley suffix: 2 332dc2eb70dSMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism \ 333dc2eb70dSMatthew G. Knepley -field_petscspace_type sum \ 334dc2eb70dSMatthew G. Knepley -field_petscspace_variables 3 \ 335dc2eb70dSMatthew G. Knepley -field_petscspace_components 3 \ 336dc2eb70dSMatthew G. Knepley -field_petscspace_sum_spaces 2 \ 337dc2eb70dSMatthew G. Knepley -field_petscspace_sum_concatenate false \ 338dc2eb70dSMatthew G. Knepley -field_sumcomp_0_petscspace_variables 3 \ 339dc2eb70dSMatthew G. Knepley -field_sumcomp_0_petscspace_components 3 \ 340dc2eb70dSMatthew G. Knepley -field_sumcomp_0_petscspace_degree 1 \ 341dc2eb70dSMatthew G. Knepley -field_sumcomp_1_petscspace_variables 3 \ 342dc2eb70dSMatthew G. Knepley -field_sumcomp_1_petscspace_components 3 \ 343dc2eb70dSMatthew G. Knepley -field_sumcomp_1_petscspace_type wxy \ 344dc2eb70dSMatthew G. Knepley -field_petscdualspace_form_degree 0 \ 345dc2eb70dSMatthew G. Knepley -field_petscdualspace_order 1 \ 346dc2eb70dSMatthew G. Knepley -field_petscdualspace_components 3 347dc2eb70dSMatthew G. Knepley 348dc2eb70dSMatthew G. Knepley TEST*/ 349