xref: /petsc/src/dm/dt/fe/tests/ex2.c (revision bd1587442888ecfa2176fd832616cd1e28a092f3)
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();
193ba16761SJacob 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]);
273ba16761SJacob 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));
613ba16761SJacob 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));
853ba16761SJacob 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));
943ba16761SJacob 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));
11103e76207SPierre Jolivet     PetscCall(PetscObjectContainerCompose((PetscObject)cellIS, composeStr, *geom, PetscContainerUserDestroy_PetscFEGeom));
112dc2eb70dSMatthew G. Knepley   }
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
114dc2eb70dSMatthew G. Knepley }
115dc2eb70dSMatthew G. Knepley 
116d71ae5a4SJacob Faibussowitsch PetscErrorCode CellRangeRestoreFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
117d71ae5a4SJacob Faibussowitsch {
118dc2eb70dSMatthew G. Knepley   PetscFunctionBegin;
119dc2eb70dSMatthew G. Knepley   *geom = NULL;
1203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
121dc2eb70dSMatthew G. Knepley }
122dc2eb70dSMatthew G. Knepley 
123d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms)
124d71ae5a4SJacob Faibussowitsch {
125dc2eb70dSMatthew G. Knepley   DMField  coordField;
126dc2eb70dSMatthew G. Knepley   PetscInt Nf, f, maxDegree;
127dc2eb70dSMatthew G. Knepley 
128dc2eb70dSMatthew G. Knepley   PetscFunctionBeginUser;
129dc2eb70dSMatthew G. Knepley   *affineQuad = NULL;
130dc2eb70dSMatthew G. Knepley   *affineGeom = NULL;
131dc2eb70dSMatthew G. Knepley   *quads      = NULL;
132dc2eb70dSMatthew G. Knepley   *geoms      = NULL;
1339566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
1349566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateField(dm, &coordField));
1359566063dSJacob Faibussowitsch   PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
136dc2eb70dSMatthew G. Knepley   if (maxDegree <= 1) {
1379566063dSJacob Faibussowitsch     PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, affineQuad));
1389566063dSJacob Faibussowitsch     if (*affineQuad) PetscCall(CellRangeGetFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom));
139dc2eb70dSMatthew G. Knepley   } else {
1409566063dSJacob Faibussowitsch     PetscCall(PetscCalloc2(Nf, quads, Nf, geoms));
141dc2eb70dSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
142dc2eb70dSMatthew G. Knepley       PetscFE fe;
143dc2eb70dSMatthew G. Knepley 
1449566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
1459566063dSJacob Faibussowitsch       PetscCall(PetscFEGetQuadrature(fe, &(*quads)[f]));
1469566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)(*quads)[f]));
1479566063dSJacob Faibussowitsch       PetscCall(CellRangeGetFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f]));
148dc2eb70dSMatthew G. Knepley     }
149dc2eb70dSMatthew G. Knepley   }
1503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
151dc2eb70dSMatthew G. Knepley }
152dc2eb70dSMatthew G. Knepley 
153d71ae5a4SJacob Faibussowitsch static PetscErrorCode DestroyFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms)
154d71ae5a4SJacob Faibussowitsch {
155dc2eb70dSMatthew G. Knepley   DMField  coordField;
156dc2eb70dSMatthew G. Knepley   PetscInt Nf, f;
157dc2eb70dSMatthew G. Knepley 
158dc2eb70dSMatthew G. Knepley   PetscFunctionBeginUser;
1599566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
1609566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateField(dm, &coordField));
161dc2eb70dSMatthew G. Knepley   if (*affineQuad) {
1629566063dSJacob Faibussowitsch     PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom));
1639566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureDestroy(affineQuad));
164dc2eb70dSMatthew G. Knepley   } else {
165dc2eb70dSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
1669566063dSJacob Faibussowitsch       PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f]));
1679566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureDestroy(&(*quads)[f]));
168dc2eb70dSMatthew G. Knepley     }
1699566063dSJacob Faibussowitsch     PetscCall(PetscFree2(*quads, *geoms));
170dc2eb70dSMatthew G. Knepley   }
1713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
172dc2eb70dSMatthew G. Knepley }
173dc2eb70dSMatthew G. Knepley 
174d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestEvaluation(DM dm)
175d71ae5a4SJacob Faibussowitsch {
176dc2eb70dSMatthew G. Knepley   PetscFE    fe;
177dc2eb70dSMatthew G. Knepley   PetscSpace sp;
178dc2eb70dSMatthew G. Knepley   PetscReal *points;
179dc2eb70dSMatthew G. Knepley   PetscReal *B, *D, *H;
180dc2eb70dSMatthew G. Knepley   PetscInt   dim, Nb, b, Nc, c, Np, p;
181dc2eb70dSMatthew G. Knepley 
182dc2eb70dSMatthew G. Knepley   PetscFunctionBeginUser;
1839566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1849566063dSJacob Faibussowitsch   PetscCall(DMGetField(dm, 0, NULL, (PetscObject *)&fe));
185dc2eb70dSMatthew G. Knepley   Np = 6;
1869566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Np * dim, &points));
187dc2eb70dSMatthew G. Knepley   if (dim == 3) {
1889371c9d4SSatish Balay     points[0]  = -1.0;
1899371c9d4SSatish Balay     points[1]  = -1.0;
1909371c9d4SSatish Balay     points[2]  = -1.0;
1919371c9d4SSatish Balay     points[3]  = 1.0;
1929371c9d4SSatish Balay     points[4]  = -1.0;
1939371c9d4SSatish Balay     points[5]  = -1.0;
1949371c9d4SSatish Balay     points[6]  = -1.0;
1959371c9d4SSatish Balay     points[7]  = 1.0;
1969371c9d4SSatish Balay     points[8]  = -1.0;
1979371c9d4SSatish Balay     points[9]  = -1.0;
1989371c9d4SSatish Balay     points[10] = -1.0;
1999371c9d4SSatish Balay     points[11] = 1.0;
2009371c9d4SSatish Balay     points[12] = 1.0;
2019371c9d4SSatish Balay     points[13] = -1.0;
2029371c9d4SSatish Balay     points[14] = 1.0;
2039371c9d4SSatish Balay     points[15] = -1.0;
2049371c9d4SSatish Balay     points[16] = 1.0;
2059371c9d4SSatish Balay     points[17] = 1.0;
206dc2eb70dSMatthew G. Knepley   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only for 3D right now");
2079566063dSJacob Faibussowitsch   PetscCall(PetscFEGetBasisSpace(fe, &sp));
2089566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetDimension(sp, &Nb));
2099566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumComponents(sp, &Nc));
2109566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, Np * Nb * Nc, MPIU_REAL, &B));
2119566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, Np * Nb * Nc * dim, MPIU_REAL, &D));
2129566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, Np * Nb * Nc * dim * dim, MPIU_REAL, &H));
2139566063dSJacob Faibussowitsch   PetscCall(PetscSpaceEvaluate(sp, Np, points, B, NULL, NULL /*D, H*/));
214dc2eb70dSMatthew G. Knepley   for (p = 0; p < Np; ++p) {
2159566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Point %" PetscInt_FMT "\n", p));
216dc2eb70dSMatthew G. Knepley     for (b = 0; b < Nb; ++b) {
2179566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "B[%" PetscInt_FMT "]:", b));
2189566063dSJacob Faibussowitsch       for (c = 0; c < Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)B[(p * Nb + b) * Nc + c]));
2199566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
220dc2eb70dSMatthew G. Knepley #if 0
221dc2eb70dSMatthew G. Knepley       for (c = 0; c < Nc; ++c) {
2229566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, " D[%" PetscInt_FMT ",%" PetscInt_FMT "]:", b, c));
2239566063dSJacob Faibussowitsch         for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double) B[((p*Nb+b)*Nc+c)*dim+d)]));
2249566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
225dc2eb70dSMatthew G. Knepley       }
226dc2eb70dSMatthew G. Knepley #endif
227dc2eb70dSMatthew G. Knepley     }
228dc2eb70dSMatthew G. Knepley   }
2299566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm, Np * Nb, MPIU_REAL, &B));
2309566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm, Np * Nb * dim, MPIU_REAL, &D));
2319566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm, Np * Nb * dim * dim, MPIU_REAL, &H));
2329566063dSJacob Faibussowitsch   PetscCall(PetscFree(points));
2333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
234dc2eb70dSMatthew G. Knepley }
235dc2eb70dSMatthew G. Knepley 
236d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestIntegration(DM dm, PetscInt cbs, PetscInt its)
237d71ae5a4SJacob Faibussowitsch {
238dc2eb70dSMatthew G. Knepley   PetscDS         ds;
239dc2eb70dSMatthew G. Knepley   PetscFEGeom    *chunkGeom = NULL;
240dc2eb70dSMatthew G. Knepley   PetscQuadrature affineQuad, *quads  = NULL;
241dc2eb70dSMatthew G. Knepley   PetscFEGeom    *affineGeom, **geoms = NULL;
242dc2eb70dSMatthew G. Knepley   PetscScalar    *u, *elemVec;
243dc2eb70dSMatthew G. Knepley   IS              cellIS;
244dc2eb70dSMatthew G. Knepley   PetscInt        depth, cStart, cEnd, cell, chunkSize = cbs, Nch = 0, Nf, f, totDim, i, k;
245dc2eb70dSMatthew G. Knepley 
246dc2eb70dSMatthew G. Knepley   PetscFunctionBeginUser;
2479566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
2489566063dSJacob Faibussowitsch   PetscCall(DMGetStratumIS(dm, "depth", depth, &cellIS));
2499566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
25007218a29SMatthew G. Knepley   PetscCall(DMGetCellDS(dm, cStart, &ds, NULL));
2519566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
2529566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
2539566063dSJacob Faibussowitsch   PetscCall(CreateFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms));
2549566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(chunkSize * totDim, &u, chunkSize * totDim, &elemVec));
255dc2eb70dSMatthew G. Knepley   /* Assumptions:
256dc2eb70dSMatthew G. Knepley     - Single field
257dc2eb70dSMatthew G. Knepley     - No input data
258dc2eb70dSMatthew G. Knepley     - No auxiliary data
259dc2eb70dSMatthew G. Knepley     - No time-dependence
260dc2eb70dSMatthew G. Knepley   */
261dc2eb70dSMatthew G. Knepley   for (i = 0; i < its; ++i) {
262dc2eb70dSMatthew G. Knepley     for (cell = cStart; cell < cEnd; cell += chunkSize, ++Nch) {
263dc2eb70dSMatthew G. Knepley       const PetscInt cS = cell, cE = PetscMin(cS + chunkSize, cEnd), Ne = cE - cS;
264dc2eb70dSMatthew G. Knepley 
2659566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(elemVec, chunkSize * totDim));
266dc2eb70dSMatthew G. Knepley       /* TODO Replace with DMPlexGetCellFields() */
267dc2eb70dSMatthew G. Knepley       for (k = 0; k < chunkSize * totDim; ++k) u[k] = 1.0;
268dc2eb70dSMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
269dc2eb70dSMatthew G. Knepley         PetscFormKey key;
270dc2eb70dSMatthew G. Knepley         PetscFEGeom *geom = affineGeom ? affineGeom : geoms[f];
271dc2eb70dSMatthew G. Knepley         /* PetscQuadrature quad = affineQuad ? affineQuad : quads[f]; */
272dc2eb70dSMatthew G. Knepley 
2739371c9d4SSatish Balay         key.label = NULL;
2749371c9d4SSatish Balay         key.value = 0;
2759371c9d4SSatish Balay         key.field = f;
2769371c9d4SSatish Balay         key.part  = 0;
2779566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomGetChunk(geom, cS, cE, &chunkGeom));
2789566063dSJacob Faibussowitsch         PetscCall(PetscFEIntegrateResidual(ds, key, Ne, chunkGeom, u, NULL, NULL, NULL, 0.0, elemVec));
279dc2eb70dSMatthew G. Knepley       }
280dc2eb70dSMatthew G. Knepley     }
281dc2eb70dSMatthew G. Knepley   }
2829566063dSJacob Faibussowitsch   PetscCall(PetscFEGeomRestoreChunk(affineGeom, cStart, cEnd, &chunkGeom));
2839566063dSJacob Faibussowitsch   PetscCall(DestroyFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms));
2849566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cellIS));
2859566063dSJacob Faibussowitsch   PetscCall(PetscFree2(u, elemVec));
2863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
287dc2eb70dSMatthew G. Knepley }
288dc2eb70dSMatthew G. Knepley 
289d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestUnisolvence(DM dm)
290d71ae5a4SJacob Faibussowitsch {
291dc2eb70dSMatthew G. Knepley   Mat M;
292dc2eb70dSMatthew G. Knepley   Vec v;
293dc2eb70dSMatthew G. Knepley 
294dc2eb70dSMatthew G. Knepley   PetscFunctionBeginUser;
2959566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &v));
2969566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &v));
2979566063dSJacob Faibussowitsch   PetscCall(DMCreateMassMatrix(dm, dm, &M));
2989566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(M, NULL, "-mass_view"));
2999566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&M));
3003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
301dc2eb70dSMatthew G. Knepley }
302dc2eb70dSMatthew G. Knepley 
303d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
304d71ae5a4SJacob Faibussowitsch {
305dc2eb70dSMatthew G. Knepley   DM          dm;
306dc2eb70dSMatthew G. Knepley   AppCtx      ctx;
307dc2eb70dSMatthew G. Knepley   PetscMPIInt size;
308dc2eb70dSMatthew G. Knepley 
309327415f7SBarry Smith   PetscFunctionBeginUser;
3109566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3119566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
312*bd158744SPierre Jolivet   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only.");
3139566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
3149566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
3159566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
3169566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
3179566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)dm, "Mesh"));
3189566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)dm, NULL, "-dm_view"));
3199566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, "field", SetupPrimalProblem, &ctx));
3209566063dSJacob Faibussowitsch   PetscCall(TestEvaluation(dm));
3219566063dSJacob Faibussowitsch   PetscCall(TestIntegration(dm, 1, ctx.its));
3229566063dSJacob Faibussowitsch   PetscCall(TestUnisolvence(dm));
3239566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
3249566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
325b122ec5aSJacob Faibussowitsch   return 0;
326dc2eb70dSMatthew G. Knepley }
327dc2eb70dSMatthew G. Knepley 
328dc2eb70dSMatthew G. Knepley /*TEST
329dc2eb70dSMatthew G. Knepley   test:
330dc2eb70dSMatthew G. Knepley     suffix: 0
331dc2eb70dSMatthew G. Knepley     args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -field_petscspace_degree 0
332dc2eb70dSMatthew G. Knepley 
333dc2eb70dSMatthew G. Knepley   test:
334dc2eb70dSMatthew G. Knepley     suffix: 2
335dc2eb70dSMatthew G. Knepley     args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism \
336dc2eb70dSMatthew G. Knepley           -field_petscspace_type sum \
337dc2eb70dSMatthew G. Knepley           -field_petscspace_variables 3 \
338dc2eb70dSMatthew G. Knepley           -field_petscspace_components 3 \
339dc2eb70dSMatthew G. Knepley           -field_petscspace_sum_spaces 2 \
340dc2eb70dSMatthew G. Knepley           -field_petscspace_sum_concatenate false \
341dc2eb70dSMatthew G. Knepley           -field_sumcomp_0_petscspace_variables 3 \
342dc2eb70dSMatthew G. Knepley           -field_sumcomp_0_petscspace_components 3 \
343dc2eb70dSMatthew G. Knepley           -field_sumcomp_0_petscspace_degree 1 \
344dc2eb70dSMatthew G. Knepley           -field_sumcomp_1_petscspace_variables 3 \
345dc2eb70dSMatthew G. Knepley           -field_sumcomp_1_petscspace_components 3 \
346dc2eb70dSMatthew G. Knepley           -field_sumcomp_1_petscspace_type wxy \
347dc2eb70dSMatthew G. Knepley           -field_petscdualspace_form_degree 0 \
348dc2eb70dSMatthew G. Knepley           -field_petscdualspace_order 1 \
349dc2eb70dSMatthew G. Knepley           -field_petscdualspace_components 3
350dc2eb70dSMatthew G. Knepley 
351dc2eb70dSMatthew G. Knepley TEST*/
352