xref: /petsc/src/dm/dt/fe/tests/ex2.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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 
11*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
12dc2eb70dSMatthew G. Knepley   PetscFunctionBeginUser;
13dc2eb70dSMatthew G. Knepley   options->its = 1;
14dc2eb70dSMatthew G. Knepley 
15d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "FE Injection Options", "PETSCFE");
169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-its", "The number of replications for timing", "ex1.c", options->its, &options->its, NULL));
17d0609cedSBarry Smith   PetscOptionsEnd();
18dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
19dc2eb70dSMatthew G. Knepley }
20dc2eb70dSMatthew G. Knepley 
21*9371c9d4SSatish Balay static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
22dc2eb70dSMatthew G. Knepley   PetscInt d;
23dc2eb70dSMatthew G. Knepley   *u = 0.0;
24dc2eb70dSMatthew G. Knepley   for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0 * PETSC_PI * x[d]);
25dc2eb70dSMatthew G. Knepley   return 0;
26dc2eb70dSMatthew G. Knepley }
27dc2eb70dSMatthew G. Knepley 
28*9371c9d4SSatish Balay 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[]) {
29dc2eb70dSMatthew G. Knepley   PetscInt d;
30dc2eb70dSMatthew G. Knepley   for (d = 0; d < dim; ++d) f0[0] += -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]);
31dc2eb70dSMatthew G. Knepley }
32dc2eb70dSMatthew G. Knepley 
33*9371c9d4SSatish Balay 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[]) {
34dc2eb70dSMatthew G. Knepley   PetscInt d;
35dc2eb70dSMatthew G. Knepley   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
36dc2eb70dSMatthew G. Knepley }
37dc2eb70dSMatthew G. Knepley 
38*9371c9d4SSatish Balay 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[]) {
39dc2eb70dSMatthew G. Knepley   PetscInt d;
40dc2eb70dSMatthew G. Knepley   for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
41dc2eb70dSMatthew G. Knepley }
42dc2eb70dSMatthew G. Knepley 
43*9371c9d4SSatish Balay static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) {
44dc2eb70dSMatthew G. Knepley   PetscDS        ds;
45dc2eb70dSMatthew G. Knepley   DMLabel        label;
46dc2eb70dSMatthew G. Knepley   const PetscInt id = 1;
47dc2eb70dSMatthew G. Knepley 
48dc2eb70dSMatthew G. Knepley   PetscFunctionBeginUser;
499566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
509566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u));
519566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
529566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(ds, 0, trig_u, user));
539566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
549566063dSJacob Faibussowitsch   if (label) PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))trig_u, NULL, user, NULL));
55dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
56dc2eb70dSMatthew G. Knepley }
57dc2eb70dSMatthew G. Knepley 
58*9371c9d4SSatish Balay static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) {
59dc2eb70dSMatthew G. Knepley   DM       cdm = dm;
60dc2eb70dSMatthew G. Knepley   PetscFE  fe;
61dc2eb70dSMatthew G. Knepley   char     prefix[PETSC_MAX_PATH_LEN];
62dc2eb70dSMatthew G. Knepley   PetscInt dim;
63dc2eb70dSMatthew G. Knepley 
64dc2eb70dSMatthew G. Knepley   PetscFunctionBeginUser;
659566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
669566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
679566063dSJacob Faibussowitsch   PetscCall(DMCreateFEDefault(dm, dim, name ? prefix : NULL, -1, &fe));
689566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe, name));
69dc2eb70dSMatthew G. Knepley   /* Set discretization and boundary conditions for each mesh */
709566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
719566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
729566063dSJacob Faibussowitsch   PetscCall((*setup)(dm, user));
73dc2eb70dSMatthew G. Knepley   while (cdm) {
749566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm, cdm));
759566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
76dc2eb70dSMatthew G. Knepley   }
779566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
78dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
79dc2eb70dSMatthew G. Knepley }
80dc2eb70dSMatthew G. Knepley 
81*9371c9d4SSatish Balay static PetscErrorCode PetscContainerUserDestroy_PetscFEGeom(void *ctx) {
82dc2eb70dSMatthew G. Knepley   PetscFEGeom *geom = (PetscFEGeom *)ctx;
83dc2eb70dSMatthew G. Knepley 
84dc2eb70dSMatthew G. Knepley   PetscFunctionBegin;
859566063dSJacob Faibussowitsch   PetscCall(PetscFEGeomDestroy(&geom));
86dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
87dc2eb70dSMatthew G. Knepley }
88dc2eb70dSMatthew G. Knepley 
89*9371c9d4SSatish Balay PetscErrorCode CellRangeGetFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) {
90dc2eb70dSMatthew G. Knepley   char           composeStr[33] = {0};
91dc2eb70dSMatthew G. Knepley   PetscObjectId  id;
92dc2eb70dSMatthew G. Knepley   PetscContainer container;
93dc2eb70dSMatthew G. Knepley 
94dc2eb70dSMatthew G. Knepley   PetscFunctionBegin;
959566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetId((PetscObject)quad, &id));
9663a3b9bcSJacob Faibussowitsch   PetscCall(PetscSNPrintf(composeStr, 32, "CellRangeGetFEGeom_%" PetscInt64_FMT "\n", id));
979566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)cellIS, composeStr, (PetscObject *)&container));
98dc2eb70dSMatthew G. Knepley   if (container) {
999566063dSJacob Faibussowitsch     PetscCall(PetscContainerGetPointer(container, (void **)geom));
100dc2eb70dSMatthew G. Knepley   } else {
1019566063dSJacob Faibussowitsch     PetscCall(DMFieldCreateFEGeom(coordField, cellIS, quad, faceData, geom));
1029566063dSJacob Faibussowitsch     PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container));
1039566063dSJacob Faibussowitsch     PetscCall(PetscContainerSetPointer(container, (void *)*geom));
1049566063dSJacob Faibussowitsch     PetscCall(PetscContainerSetUserDestroy(container, PetscContainerUserDestroy_PetscFEGeom));
1059566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)cellIS, composeStr, (PetscObject)container));
1069566063dSJacob Faibussowitsch     PetscCall(PetscContainerDestroy(&container));
107dc2eb70dSMatthew G. Knepley   }
108dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
109dc2eb70dSMatthew G. Knepley }
110dc2eb70dSMatthew G. Knepley 
111*9371c9d4SSatish Balay PetscErrorCode CellRangeRestoreFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) {
112dc2eb70dSMatthew G. Knepley   PetscFunctionBegin;
113dc2eb70dSMatthew G. Knepley   *geom = NULL;
114dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
115dc2eb70dSMatthew G. Knepley }
116dc2eb70dSMatthew G. Knepley 
117*9371c9d4SSatish Balay static PetscErrorCode CreateFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms) {
118dc2eb70dSMatthew G. Knepley   DMField  coordField;
119dc2eb70dSMatthew G. Knepley   PetscInt Nf, f, maxDegree;
120dc2eb70dSMatthew G. Knepley 
121dc2eb70dSMatthew G. Knepley   PetscFunctionBeginUser;
122dc2eb70dSMatthew G. Knepley   *affineQuad = NULL;
123dc2eb70dSMatthew G. Knepley   *affineGeom = NULL;
124dc2eb70dSMatthew G. Knepley   *quads      = NULL;
125dc2eb70dSMatthew G. Knepley   *geoms      = NULL;
1269566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
1279566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateField(dm, &coordField));
1289566063dSJacob Faibussowitsch   PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
129dc2eb70dSMatthew G. Knepley   if (maxDegree <= 1) {
1309566063dSJacob Faibussowitsch     PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, affineQuad));
1319566063dSJacob Faibussowitsch     if (*affineQuad) PetscCall(CellRangeGetFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom));
132dc2eb70dSMatthew G. Knepley   } else {
1339566063dSJacob Faibussowitsch     PetscCall(PetscCalloc2(Nf, quads, Nf, geoms));
134dc2eb70dSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
135dc2eb70dSMatthew G. Knepley       PetscFE fe;
136dc2eb70dSMatthew G. Knepley 
1379566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
1389566063dSJacob Faibussowitsch       PetscCall(PetscFEGetQuadrature(fe, &(*quads)[f]));
1399566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)(*quads)[f]));
1409566063dSJacob Faibussowitsch       PetscCall(CellRangeGetFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f]));
141dc2eb70dSMatthew G. Knepley     }
142dc2eb70dSMatthew G. Knepley   }
143dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
144dc2eb70dSMatthew G. Knepley }
145dc2eb70dSMatthew G. Knepley 
146*9371c9d4SSatish Balay static PetscErrorCode DestroyFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms) {
147dc2eb70dSMatthew G. Knepley   DMField  coordField;
148dc2eb70dSMatthew G. Knepley   PetscInt Nf, f;
149dc2eb70dSMatthew G. Knepley 
150dc2eb70dSMatthew G. Knepley   PetscFunctionBeginUser;
1519566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
1529566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateField(dm, &coordField));
153dc2eb70dSMatthew G. Knepley   if (*affineQuad) {
1549566063dSJacob Faibussowitsch     PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom));
1559566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureDestroy(affineQuad));
156dc2eb70dSMatthew G. Knepley   } else {
157dc2eb70dSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
1589566063dSJacob Faibussowitsch       PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f]));
1599566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureDestroy(&(*quads)[f]));
160dc2eb70dSMatthew G. Knepley     }
1619566063dSJacob Faibussowitsch     PetscCall(PetscFree2(*quads, *geoms));
162dc2eb70dSMatthew G. Knepley   }
163dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
164dc2eb70dSMatthew G. Knepley }
165dc2eb70dSMatthew G. Knepley 
166*9371c9d4SSatish Balay static PetscErrorCode TestEvaluation(DM dm) {
167dc2eb70dSMatthew G. Knepley   PetscFE    fe;
168dc2eb70dSMatthew G. Knepley   PetscSpace sp;
169dc2eb70dSMatthew G. Knepley   PetscReal *points;
170dc2eb70dSMatthew G. Knepley   PetscReal *B, *D, *H;
171dc2eb70dSMatthew G. Knepley   PetscInt   dim, Nb, b, Nc, c, Np, p;
172dc2eb70dSMatthew G. Knepley 
173dc2eb70dSMatthew G. Knepley   PetscFunctionBeginUser;
1749566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1759566063dSJacob Faibussowitsch   PetscCall(DMGetField(dm, 0, NULL, (PetscObject *)&fe));
176dc2eb70dSMatthew G. Knepley   Np = 6;
1779566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Np * dim, &points));
178dc2eb70dSMatthew G. Knepley   if (dim == 3) {
179*9371c9d4SSatish Balay     points[0]  = -1.0;
180*9371c9d4SSatish Balay     points[1]  = -1.0;
181*9371c9d4SSatish Balay     points[2]  = -1.0;
182*9371c9d4SSatish Balay     points[3]  = 1.0;
183*9371c9d4SSatish Balay     points[4]  = -1.0;
184*9371c9d4SSatish Balay     points[5]  = -1.0;
185*9371c9d4SSatish Balay     points[6]  = -1.0;
186*9371c9d4SSatish Balay     points[7]  = 1.0;
187*9371c9d4SSatish Balay     points[8]  = -1.0;
188*9371c9d4SSatish Balay     points[9]  = -1.0;
189*9371c9d4SSatish Balay     points[10] = -1.0;
190*9371c9d4SSatish Balay     points[11] = 1.0;
191*9371c9d4SSatish Balay     points[12] = 1.0;
192*9371c9d4SSatish Balay     points[13] = -1.0;
193*9371c9d4SSatish Balay     points[14] = 1.0;
194*9371c9d4SSatish Balay     points[15] = -1.0;
195*9371c9d4SSatish Balay     points[16] = 1.0;
196*9371c9d4SSatish Balay     points[17] = 1.0;
197dc2eb70dSMatthew G. Knepley   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only for 3D right now");
1989566063dSJacob Faibussowitsch   PetscCall(PetscFEGetBasisSpace(fe, &sp));
1999566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetDimension(sp, &Nb));
2009566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumComponents(sp, &Nc));
2019566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, Np * Nb * Nc, MPIU_REAL, &B));
2029566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, Np * Nb * Nc * dim, MPIU_REAL, &D));
2039566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, Np * Nb * Nc * dim * dim, MPIU_REAL, &H));
2049566063dSJacob Faibussowitsch   PetscCall(PetscSpaceEvaluate(sp, Np, points, B, NULL, NULL /*D, H*/));
205dc2eb70dSMatthew G. Knepley   for (p = 0; p < Np; ++p) {
2069566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Point %" PetscInt_FMT "\n", p));
207dc2eb70dSMatthew G. Knepley     for (b = 0; b < Nb; ++b) {
2089566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "B[%" PetscInt_FMT "]:", b));
2099566063dSJacob Faibussowitsch       for (c = 0; c < Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)B[(p * Nb + b) * Nc + c]));
2109566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
211dc2eb70dSMatthew G. Knepley #if 0
212dc2eb70dSMatthew G. Knepley       for (c = 0; c < Nc; ++c) {
2139566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, " D[%" PetscInt_FMT ",%" PetscInt_FMT "]:", b, c));
2149566063dSJacob Faibussowitsch         for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double) B[((p*Nb+b)*Nc+c)*dim+d)]));
2159566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
216dc2eb70dSMatthew G. Knepley       }
217dc2eb70dSMatthew G. Knepley #endif
218dc2eb70dSMatthew G. Knepley     }
219dc2eb70dSMatthew G. Knepley   }
2209566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm, Np * Nb, MPIU_REAL, &B));
2219566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm, Np * Nb * dim, MPIU_REAL, &D));
2229566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm, Np * Nb * dim * dim, MPIU_REAL, &H));
2239566063dSJacob Faibussowitsch   PetscCall(PetscFree(points));
224dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
225dc2eb70dSMatthew G. Knepley }
226dc2eb70dSMatthew G. Knepley 
227*9371c9d4SSatish Balay static PetscErrorCode TestIntegration(DM dm, PetscInt cbs, PetscInt its) {
228dc2eb70dSMatthew G. Knepley   PetscDS         ds;
229dc2eb70dSMatthew G. Knepley   PetscFEGeom    *chunkGeom = NULL;
230dc2eb70dSMatthew G. Knepley   PetscQuadrature affineQuad, *quads  = NULL;
231dc2eb70dSMatthew G. Knepley   PetscFEGeom    *affineGeom, **geoms = NULL;
232dc2eb70dSMatthew G. Knepley   PetscScalar    *u, *elemVec;
233dc2eb70dSMatthew G. Knepley   IS              cellIS;
234dc2eb70dSMatthew G. Knepley   PetscInt        depth, cStart, cEnd, cell, chunkSize = cbs, Nch = 0, Nf, f, totDim, i, k;
235dc2eb70dSMatthew G. Knepley 
236dc2eb70dSMatthew G. Knepley   PetscFunctionBeginUser;
2379566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
2389566063dSJacob Faibussowitsch   PetscCall(DMGetStratumIS(dm, "depth", depth, &cellIS));
2399566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
2409566063dSJacob Faibussowitsch   PetscCall(DMGetCellDS(dm, cStart, &ds));
2419566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
2429566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
2439566063dSJacob Faibussowitsch   PetscCall(CreateFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms));
2449566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(chunkSize * totDim, &u, chunkSize * totDim, &elemVec));
245dc2eb70dSMatthew G. Knepley   /* Assumptions:
246dc2eb70dSMatthew G. Knepley     - Single field
247dc2eb70dSMatthew G. Knepley     - No input data
248dc2eb70dSMatthew G. Knepley     - No auxiliary data
249dc2eb70dSMatthew G. Knepley     - No time-dependence
250dc2eb70dSMatthew G. Knepley   */
251dc2eb70dSMatthew G. Knepley   for (i = 0; i < its; ++i) {
252dc2eb70dSMatthew G. Knepley     for (cell = cStart; cell < cEnd; cell += chunkSize, ++Nch) {
253dc2eb70dSMatthew G. Knepley       const PetscInt cS = cell, cE = PetscMin(cS + chunkSize, cEnd), Ne = cE - cS;
254dc2eb70dSMatthew G. Knepley 
2559566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(elemVec, chunkSize * totDim));
256dc2eb70dSMatthew G. Knepley       /* TODO Replace with DMPlexGetCellFields() */
257dc2eb70dSMatthew G. Knepley       for (k = 0; k < chunkSize * totDim; ++k) u[k] = 1.0;
258dc2eb70dSMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
259dc2eb70dSMatthew G. Knepley         PetscFormKey key;
260dc2eb70dSMatthew G. Knepley         PetscFEGeom *geom = affineGeom ? affineGeom : geoms[f];
261dc2eb70dSMatthew G. Knepley         /* PetscQuadrature quad = affineQuad ? affineQuad : quads[f]; */
262dc2eb70dSMatthew G. Knepley 
263*9371c9d4SSatish Balay         key.label = NULL;
264*9371c9d4SSatish Balay         key.value = 0;
265*9371c9d4SSatish Balay         key.field = f;
266*9371c9d4SSatish Balay         key.part  = 0;
2679566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomGetChunk(geom, cS, cE, &chunkGeom));
2689566063dSJacob Faibussowitsch         PetscCall(PetscFEIntegrateResidual(ds, key, Ne, chunkGeom, u, NULL, NULL, NULL, 0.0, elemVec));
269dc2eb70dSMatthew G. Knepley       }
270dc2eb70dSMatthew G. Knepley     }
271dc2eb70dSMatthew G. Knepley   }
2729566063dSJacob Faibussowitsch   PetscCall(PetscFEGeomRestoreChunk(affineGeom, cStart, cEnd, &chunkGeom));
2739566063dSJacob Faibussowitsch   PetscCall(DestroyFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms));
2749566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cellIS));
2759566063dSJacob Faibussowitsch   PetscCall(PetscFree2(u, elemVec));
276dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
277dc2eb70dSMatthew G. Knepley }
278dc2eb70dSMatthew G. Knepley 
279*9371c9d4SSatish Balay static PetscErrorCode TestUnisolvence(DM dm) {
280dc2eb70dSMatthew G. Knepley   Mat M;
281dc2eb70dSMatthew G. Knepley   Vec v;
282dc2eb70dSMatthew G. Knepley 
283dc2eb70dSMatthew G. Knepley   PetscFunctionBeginUser;
2849566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &v));
2859566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &v));
2869566063dSJacob Faibussowitsch   PetscCall(DMCreateMassMatrix(dm, dm, &M));
2879566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(M, NULL, "-mass_view"));
2889566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&M));
289dc2eb70dSMatthew G. Knepley   PetscFunctionReturn(0);
290dc2eb70dSMatthew G. Knepley }
291dc2eb70dSMatthew G. Knepley 
292*9371c9d4SSatish Balay int main(int argc, char **argv) {
293dc2eb70dSMatthew G. Knepley   DM          dm;
294dc2eb70dSMatthew G. Knepley   AppCtx      ctx;
295dc2eb70dSMatthew G. Knepley   PetscMPIInt size;
296dc2eb70dSMatthew G. Knepley 
297327415f7SBarry Smith   PetscFunctionBeginUser;
2989566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2999566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
300dc2eb70dSMatthew G. Knepley   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_SUP, "This is a uniprocessor example only.");
3019566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
3029566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
3039566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
3049566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
3059566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)dm, "Mesh"));
3069566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)dm, NULL, "-dm_view"));
3079566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, "field", SetupPrimalProblem, &ctx));
3089566063dSJacob Faibussowitsch   PetscCall(TestEvaluation(dm));
3099566063dSJacob Faibussowitsch   PetscCall(TestIntegration(dm, 1, ctx.its));
3109566063dSJacob Faibussowitsch   PetscCall(TestUnisolvence(dm));
3119566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
3129566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
313b122ec5aSJacob Faibussowitsch   return 0;
314dc2eb70dSMatthew G. Knepley }
315dc2eb70dSMatthew G. Knepley 
316dc2eb70dSMatthew G. Knepley /*TEST
317dc2eb70dSMatthew G. Knepley   test:
318dc2eb70dSMatthew G. Knepley     suffix: 0
319dc2eb70dSMatthew G. Knepley     args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -field_petscspace_degree 0
320dc2eb70dSMatthew G. Knepley 
321dc2eb70dSMatthew G. Knepley   test:
322dc2eb70dSMatthew G. Knepley     suffix: 2
323dc2eb70dSMatthew G. Knepley     args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism \
324dc2eb70dSMatthew G. Knepley           -field_petscspace_type sum \
325dc2eb70dSMatthew G. Knepley           -field_petscspace_variables 3 \
326dc2eb70dSMatthew G. Knepley           -field_petscspace_components 3 \
327dc2eb70dSMatthew G. Knepley           -field_petscspace_sum_spaces 2 \
328dc2eb70dSMatthew G. Knepley           -field_petscspace_sum_concatenate false \
329dc2eb70dSMatthew G. Knepley           -field_sumcomp_0_petscspace_variables 3 \
330dc2eb70dSMatthew G. Knepley           -field_sumcomp_0_petscspace_components 3 \
331dc2eb70dSMatthew G. Knepley           -field_sumcomp_0_petscspace_degree 1 \
332dc2eb70dSMatthew G. Knepley           -field_sumcomp_1_petscspace_variables 3 \
333dc2eb70dSMatthew G. Knepley           -field_sumcomp_1_petscspace_components 3 \
334dc2eb70dSMatthew G. Knepley           -field_sumcomp_1_petscspace_type wxy \
335dc2eb70dSMatthew G. Knepley           -field_petscdualspace_form_degree 0 \
336dc2eb70dSMatthew G. Knepley           -field_petscdualspace_order 1 \
337dc2eb70dSMatthew G. Knepley           -field_petscdualspace_components 3
338dc2eb70dSMatthew G. Knepley 
339dc2eb70dSMatthew G. Knepley TEST*/
340