xref: /petsc/src/dm/dt/fe/tests/ex2.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   PetscErrorCode ierr;
14dc2eb70dSMatthew G. Knepley 
15dc2eb70dSMatthew G. Knepley   PetscFunctionBeginUser;
16dc2eb70dSMatthew G. Knepley   options->its = 1;
17dc2eb70dSMatthew G. Knepley 
18dc2eb70dSMatthew G. Knepley   ierr = PetscOptionsBegin(comm, "", "FE Injection Options", "PETSCFE");CHKERRQ(ierr);
195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-its", "The number of replications for timing", "ex1.c", options->its, &options->its, NULL));
20dc2eb70dSMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
21dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
22dc2eb70dSMatthew G. Knepley }
23dc2eb70dSMatthew G. Knepley 
24dc2eb70dSMatthew G. Knepley static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
25dc2eb70dSMatthew G. Knepley {
26dc2eb70dSMatthew G. Knepley   PetscInt d;
27dc2eb70dSMatthew G. Knepley   *u = 0.0;
28dc2eb70dSMatthew G. Knepley   for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0*PETSC_PI*x[d]);
29dc2eb70dSMatthew G. Knepley   return 0;
30dc2eb70dSMatthew G. Knepley }
31dc2eb70dSMatthew G. Knepley 
32dc2eb70dSMatthew G. Knepley static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
33dc2eb70dSMatthew G. Knepley                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
34dc2eb70dSMatthew G. Knepley                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
35dc2eb70dSMatthew G. Knepley                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
36dc2eb70dSMatthew G. Knepley {
37dc2eb70dSMatthew G. Knepley   PetscInt d;
38dc2eb70dSMatthew G. Knepley   for (d = 0; d < dim; ++d) f0[0] += -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]);
39dc2eb70dSMatthew G. Knepley }
40dc2eb70dSMatthew G. Knepley 
41dc2eb70dSMatthew G. Knepley static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
42dc2eb70dSMatthew G. Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
43dc2eb70dSMatthew G. Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
44dc2eb70dSMatthew G. Knepley                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
45dc2eb70dSMatthew G. Knepley {
46dc2eb70dSMatthew G. Knepley   PetscInt d;
47dc2eb70dSMatthew G. Knepley   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
48dc2eb70dSMatthew G. Knepley }
49dc2eb70dSMatthew G. Knepley 
50dc2eb70dSMatthew G. Knepley static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
51dc2eb70dSMatthew G. Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
52dc2eb70dSMatthew G. Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
53dc2eb70dSMatthew G. Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
54dc2eb70dSMatthew G. Knepley {
55dc2eb70dSMatthew G. Knepley   PetscInt d;
56dc2eb70dSMatthew G. Knepley   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
57dc2eb70dSMatthew G. Knepley }
58dc2eb70dSMatthew G. Knepley 
59dc2eb70dSMatthew G. Knepley static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
60dc2eb70dSMatthew G. Knepley {
61dc2eb70dSMatthew G. Knepley   PetscDS        ds;
62dc2eb70dSMatthew G. Knepley   DMLabel        label;
63dc2eb70dSMatthew G. Knepley   const PetscInt id = 1;
64dc2eb70dSMatthew G. Knepley 
65dc2eb70dSMatthew G. Knepley   PetscFunctionBeginUser;
665f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDS(dm, &ds));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetExactSolution(ds, 0, trig_u, user));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLabel(dm, "marker", &label));
715f80ce2aSJacob Faibussowitsch   if (label) CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) trig_u, NULL, user, NULL));
72dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
73dc2eb70dSMatthew G. Knepley }
74dc2eb70dSMatthew G. Knepley 
75dc2eb70dSMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
76dc2eb70dSMatthew G. Knepley {
77dc2eb70dSMatthew G. Knepley   DM             cdm = dm;
78dc2eb70dSMatthew G. Knepley   PetscFE        fe;
79dc2eb70dSMatthew G. Knepley   char           prefix[PETSC_MAX_PATH_LEN];
80dc2eb70dSMatthew G. Knepley   PetscInt       dim;
81dc2eb70dSMatthew G. Knepley 
82dc2eb70dSMatthew G. Knepley   PetscFunctionBeginUser;
835f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateFEDefault(dm, dim, name ? prefix : NULL, -1, &fe));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) fe, name));
87dc2eb70dSMatthew G. Knepley   /* Set discretization and boundary conditions for each mesh */
885f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateDS(dm));
905f80ce2aSJacob Faibussowitsch   CHKERRQ((*setup)(dm, user));
91dc2eb70dSMatthew G. Knepley   while (cdm) {
925f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCopyDisc(dm,cdm));
935f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoarseDM(cdm, &cdm));
94dc2eb70dSMatthew G. Knepley   }
955f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEDestroy(&fe));
96dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
97dc2eb70dSMatthew G. Knepley }
98dc2eb70dSMatthew G. Knepley 
99dc2eb70dSMatthew G. Knepley static PetscErrorCode PetscContainerUserDestroy_PetscFEGeom(void *ctx)
100dc2eb70dSMatthew G. Knepley {
101dc2eb70dSMatthew G. Knepley   PetscFEGeom   *geom = (PetscFEGeom *) ctx;
102dc2eb70dSMatthew G. Knepley 
103dc2eb70dSMatthew G. Knepley   PetscFunctionBegin;
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGeomDestroy(&geom));
105dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
106dc2eb70dSMatthew G. Knepley }
107dc2eb70dSMatthew G. Knepley 
108dc2eb70dSMatthew G. Knepley PetscErrorCode CellRangeGetFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
109dc2eb70dSMatthew G. Knepley {
110dc2eb70dSMatthew G. Knepley   char            composeStr[33] = {0};
111dc2eb70dSMatthew G. Knepley   PetscObjectId   id;
112dc2eb70dSMatthew G. Knepley   PetscContainer  container;
113dc2eb70dSMatthew G. Knepley 
114dc2eb70dSMatthew G. Knepley   PetscFunctionBegin;
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetId((PetscObject) quad, &id));
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSNPrintf(composeStr, 32, "CellRangeGetFEGeom_%x\n", id));
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectQuery((PetscObject) cellIS, composeStr, (PetscObject *) &container));
118dc2eb70dSMatthew G. Knepley   if (container) {
1195f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscContainerGetPointer(container, (void **) geom));
120dc2eb70dSMatthew G. Knepley   } else {
1215f80ce2aSJacob Faibussowitsch     CHKERRQ(DMFieldCreateFEGeom(coordField, cellIS, quad, faceData, geom));
1225f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscContainerCreate(PETSC_COMM_SELF, &container));
1235f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscContainerSetPointer(container, (void *) *geom));
1245f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscContainerSetUserDestroy(container, PetscContainerUserDestroy_PetscFEGeom));
1255f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectCompose((PetscObject) cellIS, composeStr, (PetscObject) container));
1265f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscContainerDestroy(&container));
127dc2eb70dSMatthew G. Knepley   }
128dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
129dc2eb70dSMatthew G. Knepley }
130dc2eb70dSMatthew G. Knepley 
131dc2eb70dSMatthew G. Knepley PetscErrorCode CellRangeRestoreFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
132dc2eb70dSMatthew G. Knepley {
133dc2eb70dSMatthew G. Knepley   PetscFunctionBegin;
134dc2eb70dSMatthew G. Knepley   *geom = NULL;
135dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
136dc2eb70dSMatthew G. Knepley }
137dc2eb70dSMatthew G. Knepley 
138dc2eb70dSMatthew G. Knepley static PetscErrorCode CreateFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms)
139dc2eb70dSMatthew G. Knepley {
140dc2eb70dSMatthew G. Knepley   DMField        coordField;
141dc2eb70dSMatthew G. Knepley   PetscInt       Nf, f, maxDegree;
142dc2eb70dSMatthew G. Knepley 
143dc2eb70dSMatthew G. Knepley   PetscFunctionBeginUser;
144dc2eb70dSMatthew G. Knepley   *affineQuad = NULL;
145dc2eb70dSMatthew G. Knepley   *affineGeom = NULL;
146dc2eb70dSMatthew G. Knepley   *quads      = NULL;
147dc2eb70dSMatthew G. Knepley   *geoms      = NULL;
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSGetNumFields(ds, &Nf));
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateField(dm, &coordField));
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
151dc2eb70dSMatthew G. Knepley   if (maxDegree <= 1) {
1525f80ce2aSJacob Faibussowitsch     CHKERRQ(DMFieldCreateDefaultQuadrature(coordField, cellIS, affineQuad));
1535f80ce2aSJacob Faibussowitsch     if (*affineQuad) CHKERRQ(CellRangeGetFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom));
154dc2eb70dSMatthew G. Knepley   } else {
1555f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc2(Nf, quads, Nf, geoms));
156dc2eb70dSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
157dc2eb70dSMatthew G. Knepley       PetscFE fe;
158dc2eb70dSMatthew G. Knepley 
1595f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSGetDiscretization(ds, f, (PetscObject *) &fe));
1605f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFEGetQuadrature(fe, &(*quads)[f]));
1615f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectReference((PetscObject) (*quads)[f]));
1625f80ce2aSJacob Faibussowitsch       CHKERRQ(CellRangeGetFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f]));
163dc2eb70dSMatthew G. Knepley     }
164dc2eb70dSMatthew G. Knepley   }
165dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
166dc2eb70dSMatthew G. Knepley }
167dc2eb70dSMatthew G. Knepley 
168dc2eb70dSMatthew G. Knepley static PetscErrorCode DestroyFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms)
169dc2eb70dSMatthew G. Knepley {
170dc2eb70dSMatthew G. Knepley   DMField        coordField;
171dc2eb70dSMatthew G. Knepley   PetscInt       Nf, f;
172dc2eb70dSMatthew G. Knepley 
173dc2eb70dSMatthew G. Knepley   PetscFunctionBeginUser;
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSGetNumFields(ds, &Nf));
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateField(dm, &coordField));
176dc2eb70dSMatthew G. Knepley   if (*affineQuad) {
1775f80ce2aSJacob Faibussowitsch     CHKERRQ(CellRangeRestoreFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom));
1785f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscQuadratureDestroy(affineQuad));
179dc2eb70dSMatthew G. Knepley   } else {
180dc2eb70dSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
1815f80ce2aSJacob Faibussowitsch       CHKERRQ(CellRangeRestoreFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f]));
1825f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscQuadratureDestroy(&(*quads)[f]));
183dc2eb70dSMatthew G. Knepley     }
1845f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree2(*quads, *geoms));
185dc2eb70dSMatthew G. Knepley   }
186dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
187dc2eb70dSMatthew G. Knepley }
188dc2eb70dSMatthew G. Knepley 
189dc2eb70dSMatthew G. Knepley static PetscErrorCode TestEvaluation(DM dm)
190dc2eb70dSMatthew G. Knepley {
191dc2eb70dSMatthew G. Knepley   PetscFE        fe;
192dc2eb70dSMatthew G. Knepley   PetscSpace     sp;
193dc2eb70dSMatthew G. Knepley   PetscReal     *points;
194dc2eb70dSMatthew G. Knepley   PetscReal     *B, *D, *H;
195dc2eb70dSMatthew G. Knepley   PetscInt       dim, Nb, b, Nc, c, Np, p;
196dc2eb70dSMatthew G. Knepley 
197dc2eb70dSMatthew G. Knepley   PetscFunctionBeginUser;
1985f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
1995f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetField(dm, 0, NULL, (PetscObject *) &fe));
200dc2eb70dSMatthew G. Knepley   Np = 6;
2015f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(Np*dim, &points));
202dc2eb70dSMatthew G. Knepley   if (dim == 3) {
203dc2eb70dSMatthew G. Knepley     points[0]  = -1.0; points[1]  = -1.0; points[2]  = -1.0;
204dc2eb70dSMatthew G. Knepley     points[3]  =  1.0; points[4]  = -1.0; points[5]  = -1.0;
205dc2eb70dSMatthew G. Knepley     points[6]  = -1.0; points[7]  =  1.0; points[8]  = -1.0;
206dc2eb70dSMatthew G. Knepley     points[9]  = -1.0; points[10] = -1.0; points[11] =  1.0;
207dc2eb70dSMatthew G. Knepley     points[12] =  1.0; points[13] = -1.0; points[14] =  1.0;
208dc2eb70dSMatthew G. Knepley     points[15] = -1.0; points[16] =  1.0; points[17] =  1.0;
209dc2eb70dSMatthew G. Knepley   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only for 3D right now");
2105f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetBasisSpace(fe, &sp));
2115f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpaceGetDimension(sp, &Nb));
2125f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpaceGetNumComponents(sp, &Nc));
2135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetWorkArray(dm, Np*Nb*Nc, MPIU_REAL, &B));
2145f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetWorkArray(dm, Np*Nb*Nc*dim, MPIU_REAL, &D));
2155f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetWorkArray(dm, Np*Nb*Nc*dim*dim, MPIU_REAL, &H));
2165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpaceEvaluate(sp, Np, points, B, NULL, NULL /*D, H*/));
217dc2eb70dSMatthew G. Knepley   for (p = 0; p < Np; ++p) {
2185f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "Point %" PetscInt_FMT "\n", p));
219dc2eb70dSMatthew G. Knepley     for (b = 0; b < Nb; ++b) {
2205f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "B[%" PetscInt_FMT "]:", b));
2215f80ce2aSJacob Faibussowitsch       for (c = 0; c < Nc; ++c) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " %g", (double) B[(p*Nb+b)*Nc+c]));
2225f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "\n"));
223dc2eb70dSMatthew G. Knepley #if 0
224dc2eb70dSMatthew G. Knepley       for (c = 0; c < Nc; ++c) {
2255f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " D[%" PetscInt_FMT ",%" PetscInt_FMT "]:", b, c));
2265f80ce2aSJacob Faibussowitsch         for (d = 0; d < dim; ++d) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " %g", (double) B[((p*Nb+b)*Nc+c)*dim+d)]));
2275f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "\n"));
228dc2eb70dSMatthew G. Knepley       }
229dc2eb70dSMatthew G. Knepley #endif
230dc2eb70dSMatthew G. Knepley     }
231dc2eb70dSMatthew G. Knepley   }
2325f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreWorkArray(dm, Np*Nb, MPIU_REAL, &B));
2335f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreWorkArray(dm, Np*Nb*dim, MPIU_REAL, &D));
2345f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreWorkArray(dm, Np*Nb*dim*dim, MPIU_REAL, &H));
2355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(points));
236dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
237dc2eb70dSMatthew G. Knepley }
238dc2eb70dSMatthew G. Knepley 
239dc2eb70dSMatthew G. Knepley static PetscErrorCode TestIntegration(DM dm, PetscInt cbs, PetscInt its)
240dc2eb70dSMatthew G. Knepley {
241dc2eb70dSMatthew G. Knepley   PetscDS         ds;
242dc2eb70dSMatthew G. Knepley   PetscFEGeom    *chunkGeom = NULL;
243dc2eb70dSMatthew G. Knepley   PetscQuadrature affineQuad,  *quads = NULL;
244dc2eb70dSMatthew G. Knepley   PetscFEGeom    *affineGeom, **geoms = NULL;
245dc2eb70dSMatthew G. Knepley   PetscScalar    *u, *elemVec;
246dc2eb70dSMatthew G. Knepley   IS              cellIS;
247dc2eb70dSMatthew G. Knepley   PetscInt        depth, cStart, cEnd, cell, chunkSize = cbs, Nch = 0, Nf, f, totDim, i, k;
248dc2eb70dSMatthew G. Knepley 
249dc2eb70dSMatthew G. Knepley   PetscFunctionBeginUser;
2505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetDepth(dm, &depth));
2515f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetStratumIS(dm, "depth", depth, &cellIS));
2525f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
2535f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCellDS(dm, cStart, &ds));
2545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSGetNumFields(ds, &Nf));
2555f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSGetTotalDimension(ds, &totDim));
2565f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms));
2575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(chunkSize*totDim, &u, chunkSize*totDim, &elemVec));
258dc2eb70dSMatthew G. Knepley   /* Assumptions:
259dc2eb70dSMatthew G. Knepley     - Single field
260dc2eb70dSMatthew G. Knepley     - No input data
261dc2eb70dSMatthew G. Knepley     - No auxiliary data
262dc2eb70dSMatthew G. Knepley     - No time-dependence
263dc2eb70dSMatthew G. Knepley   */
264dc2eb70dSMatthew G. Knepley   for (i = 0; i < its; ++i) {
265dc2eb70dSMatthew G. Knepley     for (cell = cStart; cell < cEnd; cell += chunkSize, ++Nch) {
266dc2eb70dSMatthew G. Knepley       const PetscInt cS = cell, cE = PetscMin(cS + chunkSize, cEnd), Ne = cE - cS;
267dc2eb70dSMatthew G. Knepley 
2685f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(elemVec, chunkSize*totDim));
269dc2eb70dSMatthew G. Knepley       /* TODO Replace with DMPlexGetCellFields() */
270dc2eb70dSMatthew G. Knepley       for (k = 0; k < chunkSize*totDim; ++k) u[k] = 1.0;
271dc2eb70dSMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
272dc2eb70dSMatthew G. Knepley         PetscFormKey key;
273dc2eb70dSMatthew G. Knepley         PetscFEGeom     *geom = affineGeom ? affineGeom : geoms[f];
274dc2eb70dSMatthew G. Knepley         /* PetscQuadrature quad = affineQuad ? affineQuad : quads[f]; */
275dc2eb70dSMatthew G. Knepley 
276dc2eb70dSMatthew G. Knepley         key.label = NULL; key.value = 0; key.field = f;
2775f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscFEGeomGetChunk(geom, cS, cE, &chunkGeom));
2785f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscFEIntegrateResidual(ds, key, Ne, chunkGeom, u, NULL, NULL, NULL, 0.0, elemVec));
279dc2eb70dSMatthew G. Knepley       }
280dc2eb70dSMatthew G. Knepley     }
281dc2eb70dSMatthew G. Knepley   }
2825f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGeomRestoreChunk(affineGeom, cStart, cEnd, &chunkGeom));
2835f80ce2aSJacob Faibussowitsch   CHKERRQ(DestroyFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms));
2845f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&cellIS));
2855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(u, elemVec));
286dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
287dc2eb70dSMatthew G. Knepley }
288dc2eb70dSMatthew G. Knepley 
289dc2eb70dSMatthew G. Knepley static PetscErrorCode TestUnisolvence(DM dm)
290dc2eb70dSMatthew G. Knepley {
291dc2eb70dSMatthew G. Knepley   Mat M;
292dc2eb70dSMatthew G. Knepley   Vec v;
293dc2eb70dSMatthew G. Knepley 
294dc2eb70dSMatthew G. Knepley   PetscFunctionBeginUser;
2955f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(dm, &v));
2965f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(dm, &v));
2975f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateMassMatrix(dm, dm, &M));
2985f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(M, NULL, "-mass_view"));
2995f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&M));
300dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
301dc2eb70dSMatthew G. Knepley }
302dc2eb70dSMatthew G. Knepley 
303dc2eb70dSMatthew G. Knepley int main(int argc, char **argv)
304dc2eb70dSMatthew G. Knepley {
305dc2eb70dSMatthew G. Knepley   DM             dm;
306dc2eb70dSMatthew G. Knepley   AppCtx         ctx;
307dc2eb70dSMatthew G. Knepley   PetscMPIInt    size;
308dc2eb70dSMatthew G. Knepley 
309*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL, help));
3105f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
311dc2eb70dSMatthew G. Knepley   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_SUP, "This is a uniprocessor example only.");
3125f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &ctx));
3135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(PETSC_COMM_WORLD, &dm));
3145f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(dm, DMPLEX));
3155f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(dm));
3165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) dm, "Mesh"));
3175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectViewFromOptions((PetscObject) dm, NULL, "-dm_view"));
3185f80ce2aSJacob Faibussowitsch   CHKERRQ(SetupDiscretization(dm, "field", SetupPrimalProblem, &ctx));
3195f80ce2aSJacob Faibussowitsch   CHKERRQ(TestEvaluation(dm));
3205f80ce2aSJacob Faibussowitsch   CHKERRQ(TestIntegration(dm, 1, ctx.its));
3215f80ce2aSJacob Faibussowitsch   CHKERRQ(TestUnisolvence(dm));
3225f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
323*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
324*b122ec5aSJacob Faibussowitsch   return 0;
325dc2eb70dSMatthew G. Knepley }
326dc2eb70dSMatthew G. Knepley 
327dc2eb70dSMatthew G. Knepley /*TEST
328dc2eb70dSMatthew G. Knepley   test:
329dc2eb70dSMatthew G. Knepley     suffix: 0
330dc2eb70dSMatthew G. Knepley     args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -field_petscspace_degree 0
331dc2eb70dSMatthew G. Knepley 
332dc2eb70dSMatthew G. Knepley   test:
333dc2eb70dSMatthew G. Knepley     suffix: 2
334dc2eb70dSMatthew G. Knepley     args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism \
335dc2eb70dSMatthew G. Knepley           -field_petscspace_type sum \
336dc2eb70dSMatthew G. Knepley           -field_petscspace_variables 3 \
337dc2eb70dSMatthew G. Knepley           -field_petscspace_components 3 \
338dc2eb70dSMatthew G. Knepley           -field_petscspace_sum_spaces 2 \
339dc2eb70dSMatthew G. Knepley           -field_petscspace_sum_concatenate false \
340dc2eb70dSMatthew G. Knepley           -field_sumcomp_0_petscspace_variables 3 \
341dc2eb70dSMatthew G. Knepley           -field_sumcomp_0_petscspace_components 3 \
342dc2eb70dSMatthew G. Knepley           -field_sumcomp_0_petscspace_degree 1 \
343dc2eb70dSMatthew G. Knepley           -field_sumcomp_1_petscspace_variables 3 \
344dc2eb70dSMatthew G. Knepley           -field_sumcomp_1_petscspace_components 3 \
345dc2eb70dSMatthew G. Knepley           -field_sumcomp_1_petscspace_type wxy \
346dc2eb70dSMatthew G. Knepley           -field_petscdualspace_form_degree 0 \
347dc2eb70dSMatthew G. Knepley           -field_petscdualspace_order 1 \
348dc2eb70dSMatthew G. Knepley           -field_petscdualspace_components 3
349dc2eb70dSMatthew G. Knepley 
350dc2eb70dSMatthew G. Knepley TEST*/
351