xref: /petsc/src/dm/dt/fe/tests/ex2.c (revision 6e415bd26dfd0d33b8d7233bdee8a51de4844a26)
1 static const char help[] = "Tests for injecting basis functions";
2 
3 #include <petscdmplex.h>
4 #include <petscfe.h>
5 #include <petscds.h>
6 
7 typedef struct {
8   PetscInt its; /* Number of replications for timing */
9 } AppCtx;
10 
11 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
12 {
13   PetscFunctionBeginUser;
14   options->its = 1;
15 
16   PetscOptionsBegin(comm, "", "FE Injection Options", "PETSCFE");
17   PetscCall(PetscOptionsInt("-its", "The number of replications for timing", "ex1.c", options->its, &options->its, NULL));
18   PetscOptionsEnd();
19   PetscFunctionReturn(PETSC_SUCCESS);
20 }
21 
22 static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
23 {
24   PetscInt d;
25   *u = 0.0;
26   for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0 * PETSC_PI * x[d]);
27   return PETSC_SUCCESS;
28 }
29 
30 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[])
31 {
32   PetscInt d;
33   for (d = 0; d < dim; ++d) f0[0] += -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]);
34 }
35 
36 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[])
37 {
38   PetscInt d;
39   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
40 }
41 
42 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[])
43 {
44   PetscInt d;
45   for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
46 }
47 
48 static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
49 {
50   PetscDS        ds;
51   DMLabel        label;
52   const PetscInt id = 1;
53 
54   PetscFunctionBeginUser;
55   PetscCall(DMGetDS(dm, &ds));
56   PetscCall(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u));
57   PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
58   PetscCall(PetscDSSetExactSolution(ds, 0, trig_u, user));
59   PetscCall(DMGetLabel(dm, "marker", &label));
60   if (label) PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))trig_u, NULL, user, NULL));
61   PetscFunctionReturn(PETSC_SUCCESS);
62 }
63 
64 static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
65 {
66   DM       cdm = dm;
67   PetscFE  fe;
68   char     prefix[PETSC_MAX_PATH_LEN];
69   PetscInt dim;
70 
71   PetscFunctionBeginUser;
72   PetscCall(DMGetDimension(dm, &dim));
73   PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
74   PetscCall(DMCreateFEDefault(dm, dim, name ? prefix : NULL, -1, &fe));
75   PetscCall(PetscObjectSetName((PetscObject)fe, name));
76   /* Set discretization and boundary conditions for each mesh */
77   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
78   PetscCall(DMCreateDS(dm));
79   PetscCall((*setup)(dm, user));
80   while (cdm) {
81     PetscCall(DMCopyDisc(dm, cdm));
82     PetscCall(DMGetCoarseDM(cdm, &cdm));
83   }
84   PetscCall(PetscFEDestroy(&fe));
85   PetscFunctionReturn(PETSC_SUCCESS);
86 }
87 
88 static PetscErrorCode PetscContainerUserDestroy_PetscFEGeom(void *ctx)
89 {
90   PetscFEGeom *geom = (PetscFEGeom *)ctx;
91 
92   PetscFunctionBegin;
93   PetscCall(PetscFEGeomDestroy(&geom));
94   PetscFunctionReturn(PETSC_SUCCESS);
95 }
96 
97 PetscErrorCode CellRangeGetFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
98 {
99   char           composeStr[33] = {0};
100   PetscObjectId  id;
101   PetscContainer container;
102 
103   PetscFunctionBegin;
104   PetscCall(PetscObjectGetId((PetscObject)quad, &id));
105   PetscCall(PetscSNPrintf(composeStr, 32, "CellRangeGetFEGeom_%" PetscInt64_FMT "\n", id));
106   PetscCall(PetscObjectQuery((PetscObject)cellIS, composeStr, (PetscObject *)&container));
107   if (container) {
108     PetscCall(PetscContainerGetPointer(container, (void **)geom));
109   } else {
110     PetscCall(DMFieldCreateFEGeom(coordField, cellIS, quad, faceData, geom));
111     PetscCall(PetscObjectContainerCompose((PetscObject)cellIS, composeStr, *geom, PetscContainerUserDestroy_PetscFEGeom));
112   }
113   PetscFunctionReturn(PETSC_SUCCESS);
114 }
115 
116 PetscErrorCode CellRangeRestoreFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
117 {
118   PetscFunctionBegin;
119   *geom = NULL;
120   PetscFunctionReturn(PETSC_SUCCESS);
121 }
122 
123 static PetscErrorCode CreateFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms)
124 {
125   DMField  coordField;
126   PetscInt Nf, f, maxDegree;
127 
128   PetscFunctionBeginUser;
129   *affineQuad = NULL;
130   *affineGeom = NULL;
131   *quads      = NULL;
132   *geoms      = NULL;
133   PetscCall(PetscDSGetNumFields(ds, &Nf));
134   PetscCall(DMGetCoordinateField(dm, &coordField));
135   PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
136   if (maxDegree <= 1) {
137     PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, affineQuad));
138     if (*affineQuad) PetscCall(CellRangeGetFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom));
139   } else {
140     PetscCall(PetscCalloc2(Nf, quads, Nf, geoms));
141     for (f = 0; f < Nf; ++f) {
142       PetscFE fe;
143 
144       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
145       PetscCall(PetscFEGetQuadrature(fe, &(*quads)[f]));
146       PetscCall(PetscObjectReference((PetscObject)(*quads)[f]));
147       PetscCall(CellRangeGetFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f]));
148     }
149   }
150   PetscFunctionReturn(PETSC_SUCCESS);
151 }
152 
153 static PetscErrorCode DestroyFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms)
154 {
155   DMField  coordField;
156   PetscInt Nf, f;
157 
158   PetscFunctionBeginUser;
159   PetscCall(PetscDSGetNumFields(ds, &Nf));
160   PetscCall(DMGetCoordinateField(dm, &coordField));
161   if (*affineQuad) {
162     PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom));
163     PetscCall(PetscQuadratureDestroy(affineQuad));
164   } else {
165     for (f = 0; f < Nf; ++f) {
166       PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f]));
167       PetscCall(PetscQuadratureDestroy(&(*quads)[f]));
168     }
169     PetscCall(PetscFree2(*quads, *geoms));
170   }
171   PetscFunctionReturn(PETSC_SUCCESS);
172 }
173 
174 static PetscErrorCode TestEvaluation(DM dm)
175 {
176   PetscFE    fe;
177   PetscSpace sp;
178   PetscReal *points;
179   PetscReal *B, *D, *H;
180   PetscInt   dim, Nb, b, Nc, c, Np, p;
181 
182   PetscFunctionBeginUser;
183   PetscCall(DMGetDimension(dm, &dim));
184   PetscCall(DMGetField(dm, 0, NULL, (PetscObject *)&fe));
185   Np = 6;
186   PetscCall(PetscMalloc1(Np * dim, &points));
187   if (dim == 3) {
188     points[0]  = -1.0;
189     points[1]  = -1.0;
190     points[2]  = -1.0;
191     points[3]  = 1.0;
192     points[4]  = -1.0;
193     points[5]  = -1.0;
194     points[6]  = -1.0;
195     points[7]  = 1.0;
196     points[8]  = -1.0;
197     points[9]  = -1.0;
198     points[10] = -1.0;
199     points[11] = 1.0;
200     points[12] = 1.0;
201     points[13] = -1.0;
202     points[14] = 1.0;
203     points[15] = -1.0;
204     points[16] = 1.0;
205     points[17] = 1.0;
206   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only for 3D right now");
207   PetscCall(PetscFEGetBasisSpace(fe, &sp));
208   PetscCall(PetscSpaceGetDimension(sp, &Nb));
209   PetscCall(PetscSpaceGetNumComponents(sp, &Nc));
210   PetscCall(DMGetWorkArray(dm, Np * Nb * Nc, MPIU_REAL, &B));
211   PetscCall(DMGetWorkArray(dm, Np * Nb * Nc * dim, MPIU_REAL, &D));
212   PetscCall(DMGetWorkArray(dm, Np * Nb * Nc * dim * dim, MPIU_REAL, &H));
213   PetscCall(PetscSpaceEvaluate(sp, Np, points, B, NULL, NULL /*D, H*/));
214   for (p = 0; p < Np; ++p) {
215     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Point %" PetscInt_FMT "\n", p));
216     for (b = 0; b < Nb; ++b) {
217       PetscCall(PetscPrintf(PETSC_COMM_SELF, "B[%" PetscInt_FMT "]:", b));
218       for (c = 0; c < Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)B[(p * Nb + b) * Nc + c]));
219       PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
220 #if 0
221       for (c = 0; c < Nc; ++c) {
222         PetscCall(PetscPrintf(PETSC_COMM_SELF, " D[%" PetscInt_FMT ",%" PetscInt_FMT "]:", b, c));
223         for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double) B[((p*Nb+b)*Nc+c)*dim+d)]));
224         PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
225       }
226 #endif
227     }
228   }
229   PetscCall(DMRestoreWorkArray(dm, Np * Nb, MPIU_REAL, &B));
230   PetscCall(DMRestoreWorkArray(dm, Np * Nb * dim, MPIU_REAL, &D));
231   PetscCall(DMRestoreWorkArray(dm, Np * Nb * dim * dim, MPIU_REAL, &H));
232   PetscCall(PetscFree(points));
233   PetscFunctionReturn(PETSC_SUCCESS);
234 }
235 
236 static PetscErrorCode TestIntegration(DM dm, PetscInt cbs, PetscInt its)
237 {
238   PetscDS         ds;
239   PetscFEGeom    *chunkGeom = NULL;
240   PetscQuadrature affineQuad, *quads  = NULL;
241   PetscFEGeom    *affineGeom, **geoms = NULL;
242   PetscScalar    *u, *elemVec;
243   IS              cellIS;
244   PetscInt        depth, cStart, cEnd, cell, chunkSize = cbs, Nch = 0, Nf, f, totDim, i, k;
245 
246   PetscFunctionBeginUser;
247   PetscCall(DMPlexGetDepth(dm, &depth));
248   PetscCall(DMGetStratumIS(dm, "depth", depth, &cellIS));
249   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
250   PetscCall(DMGetCellDS(dm, cStart, &ds, NULL));
251   PetscCall(PetscDSGetNumFields(ds, &Nf));
252   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
253   PetscCall(CreateFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms));
254   PetscCall(PetscMalloc2(chunkSize * totDim, &u, chunkSize * totDim, &elemVec));
255   /* Assumptions:
256     - Single field
257     - No input data
258     - No auxiliary data
259     - No time-dependence
260   */
261   for (i = 0; i < its; ++i) {
262     for (cell = cStart; cell < cEnd; cell += chunkSize, ++Nch) {
263       const PetscInt cS = cell, cE = PetscMin(cS + chunkSize, cEnd), Ne = cE - cS;
264 
265       PetscCall(PetscArrayzero(elemVec, chunkSize * totDim));
266       /* TODO Replace with DMPlexGetCellFields() */
267       for (k = 0; k < chunkSize * totDim; ++k) u[k] = 1.0;
268       for (f = 0; f < Nf; ++f) {
269         PetscFormKey key;
270         PetscFEGeom *geom = affineGeom ? affineGeom : geoms[f];
271         /* PetscQuadrature quad = affineQuad ? affineQuad : quads[f]; */
272 
273         key.label = NULL;
274         key.value = 0;
275         key.field = f;
276         key.part  = 0;
277         PetscCall(PetscFEGeomGetChunk(geom, cS, cE, &chunkGeom));
278         PetscCall(PetscFEIntegrateResidual(ds, key, Ne, chunkGeom, u, NULL, NULL, NULL, 0.0, elemVec));
279       }
280     }
281   }
282   PetscCall(PetscFEGeomRestoreChunk(affineGeom, cStart, cEnd, &chunkGeom));
283   PetscCall(DestroyFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms));
284   PetscCall(ISDestroy(&cellIS));
285   PetscCall(PetscFree2(u, elemVec));
286   PetscFunctionReturn(PETSC_SUCCESS);
287 }
288 
289 static PetscErrorCode TestUnisolvence(DM dm)
290 {
291   Mat M;
292   Vec v;
293 
294   PetscFunctionBeginUser;
295   PetscCall(DMGetLocalVector(dm, &v));
296   PetscCall(DMRestoreLocalVector(dm, &v));
297   PetscCall(DMCreateMassMatrix(dm, dm, &M));
298   PetscCall(MatViewFromOptions(M, NULL, "-mass_view"));
299   PetscCall(MatDestroy(&M));
300   PetscFunctionReturn(PETSC_SUCCESS);
301 }
302 
303 int main(int argc, char **argv)
304 {
305   DM          dm;
306   AppCtx      ctx;
307   PetscMPIInt size;
308 
309   PetscFunctionBeginUser;
310   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
311   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
312   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only.");
313   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
314   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
315   PetscCall(DMSetType(dm, DMPLEX));
316   PetscCall(DMSetFromOptions(dm));
317   PetscCall(PetscObjectSetName((PetscObject)dm, "Mesh"));
318   PetscCall(PetscObjectViewFromOptions((PetscObject)dm, NULL, "-dm_view"));
319   PetscCall(SetupDiscretization(dm, "field", SetupPrimalProblem, &ctx));
320   PetscCall(TestEvaluation(dm));
321   PetscCall(TestIntegration(dm, 1, ctx.its));
322   PetscCall(TestUnisolvence(dm));
323   PetscCall(DMDestroy(&dm));
324   PetscCall(PetscFinalize());
325   return 0;
326 }
327 
328 /*TEST
329   test:
330     suffix: 0
331     args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -field_petscspace_degree 0
332 
333   test:
334     suffix: 2
335     args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism \
336           -field_petscspace_type sum \
337           -field_petscspace_variables 3 \
338           -field_petscspace_components 3 \
339           -field_petscspace_sum_spaces 2 \
340           -field_petscspace_sum_concatenate false \
341           -field_sumcomp_0_petscspace_variables 3 \
342           -field_sumcomp_0_petscspace_components 3 \
343           -field_sumcomp_0_petscspace_degree 1 \
344           -field_sumcomp_1_petscspace_variables 3 \
345           -field_sumcomp_1_petscspace_components 3 \
346           -field_sumcomp_1_petscspace_type wxy \
347           -field_petscdualspace_form_degree 0 \
348           -field_petscdualspace_order 1 \
349           -field_petscdualspace_components 3
350 
351 TEST*/
352