xref: /petsc/src/dm/impls/plex/tests/ex42.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1 static const char help[] = "Simple libCEED test to calculate surface area using 1^T M 1";
2 
3 /*
4   This is a recreation of libCeed Example 2: https://libceed.readthedocs.io/en/latest/examples/ceed/
5 */
6 
7 #include <petscdmceed.h>
8 #include <petscdmplexceed.h>
9 #include <petscfeceed.h>
10 #include <petscdmplex.h>
11 #include <petscds.h>
12 
13 typedef struct {
14   PetscReal         areaExact;
15   CeedQFunctionUser setupgeo, apply;
16   const char       *setupgeofname, *applyfname;
17 } AppCtx;
18 
19 typedef struct {
20   CeedQFunction qf_apply;
21   CeedOperator  op_apply;
22   CeedVector    qdata, uceed, vceed;
23 } CeedData;
24 
25 static PetscErrorCode CeedDataDestroy(CeedData *data) {
26   PetscFunctionBeginUser;
27   PetscCall(CeedVectorDestroy(&data->qdata));
28   PetscCall(CeedVectorDestroy(&data->uceed));
29   PetscCall(CeedVectorDestroy(&data->vceed));
30   PetscCall(CeedQFunctionDestroy(&data->qf_apply));
31   PetscCall(CeedOperatorDestroy(&data->op_apply));
32   PetscFunctionReturn(0);
33 }
34 
35 CEED_QFUNCTION(Mass)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
36   const CeedScalar *u = in[0], *qdata = in[1];
37   CeedScalar       *v = out[0];
38 
39   CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i) v[i] = qdata[i] * u[i];
40 
41   return 0;
42 }
43 
44 /*
45 // Reference (parent) 2D coordinates: X \in [-1, 1]^2
46 //
47 // Global physical coordinates given by the mesh (3D): xx \in [-l, l]^3
48 //
49 // Local physical coordinates on the manifold (2D): x \in [-l, l]^2
50 //
51 // Change of coordinates matrix computed by the library:
52 //   (physical 3D coords relative to reference 2D coords)
53 //   dxx_j/dX_i (indicial notation) [3 * 2]
54 //
55 // Change of coordinates x (physical 2D) relative to xx (phyisical 3D):
56 //   dx_i/dxx_j (indicial notation) [2 * 3]
57 //
58 // Change of coordinates x (physical 2D) relative to X (reference 2D):
59 //   (by chain rule)
60 //   dx_i/dX_j = dx_i/dxx_k * dxx_k/dX_j
61 //
62 // The quadrature data is stored in the array qdata.
63 //
64 // We require the determinant of the Jacobian to properly compute integrals of the form: int(u v)
65 //
66 // Qdata: w * det(dx_i/dX_j)
67 */
68 CEED_QFUNCTION(SetupMassGeoCube)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
69   const CeedScalar *J = in[1], *w = in[2];
70   CeedScalar       *qdata = out[0];
71 
72   CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i) {
73     // Read dxxdX Jacobian entries, stored as [[0 3], [1 4], [2 5]]
74     const CeedScalar dxxdX[3][2] = {
75       {J[i + Q * 0], J[i + Q * 3]},
76       {J[i + Q * 1], J[i + Q * 4]},
77       {J[i + Q * 2], J[i + Q * 5]}
78     };
79     // Modulus of dxxdX column vectors
80     const CeedScalar modg1       = PetscSqrtReal(dxxdX[0][0] * dxxdX[0][0] + dxxdX[1][0] * dxxdX[1][0] + dxxdX[2][0] * dxxdX[2][0]);
81     const CeedScalar modg2       = PetscSqrtReal(dxxdX[0][1] * dxxdX[0][1] + dxxdX[1][1] * dxxdX[1][1] + dxxdX[2][1] * dxxdX[2][1]);
82     // Use normalized column vectors of dxxdX as rows of dxdxx
83     const CeedScalar dxdxx[2][3] = {
84       {dxxdX[0][0] / modg1, dxxdX[1][0] / modg1, dxxdX[2][0] / modg1},
85       {dxxdX[0][1] / modg2, dxxdX[1][1] / modg2, dxxdX[2][1] / modg2}
86     };
87 
88     CeedScalar dxdX[2][2];
89     for (int j = 0; j < 2; ++j)
90       for (int k = 0; k < 2; ++k) {
91         dxdX[j][k] = 0;
92         for (int l = 0; l < 3; ++l) dxdX[j][k] += dxdxx[j][l] * dxxdX[l][k];
93       }
94     qdata[i + Q * 0] = (dxdX[0][0] * dxdX[1][1] - dxdX[1][0] * dxdX[0][1]) * w[i]; /* det J * weight */
95   }
96   return 0;
97 }
98 
99 /*
100 // Reference (parent) 2D coordinates: X \in [-1, 1]^2
101 //
102 // Global 3D physical coordinates given by the mesh: xx \in [-R, R]^3
103 //   with R radius of the sphere
104 //
105 // Local 3D physical coordinates on the 2D manifold: x \in [-l, l]^3
106 //   with l half edge of the cube inscribed in the sphere
107 //
108 // Change of coordinates matrix computed by the library:
109 //   (physical 3D coords relative to reference 2D coords)
110 //   dxx_j/dX_i (indicial notation) [3 * 2]
111 //
112 // Change of coordinates x (on the 2D manifold) relative to xx (phyisical 3D):
113 //   dx_i/dxx_j (indicial notation) [3 * 3]
114 //
115 // Change of coordinates x (on the 2D manifold) relative to X (reference 2D):
116 //   (by chain rule)
117 //   dx_i/dX_j = dx_i/dxx_k * dxx_k/dX_j [3 * 2]
118 //
119 // modJ is given by the magnitude of the cross product of the columns of dx_i/dX_j
120 //
121 // The quadrature data is stored in the array qdata.
122 //
123 // We require the determinant of the Jacobian to properly compute integrals of
124 //   the form: int(u v)
125 //
126 // Qdata: modJ * w
127 */
128 CEED_QFUNCTION(SetupMassGeoSphere)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
129   const CeedScalar *X = in[0], *J = in[1], *w = in[2];
130   CeedScalar       *qdata = out[0];
131 
132   CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i) {
133     const CeedScalar xx[3][1]    = {{X[i + 0 * Q]}, {X[i + 1 * Q]}, {X[i + 2 * Q]}};
134     // Read dxxdX Jacobian entries, stored as [[0 3], [1 4], [2 5]]
135     const CeedScalar dxxdX[3][2] = {
136       {J[i + Q * 0], J[i + Q * 3]},
137       {J[i + Q * 1], J[i + Q * 4]},
138       {J[i + Q * 2], J[i + Q * 5]}
139     };
140     // Setup
141     const CeedScalar modxxsq = xx[0][0] * xx[0][0] + xx[1][0] * xx[1][0] + xx[2][0] * xx[2][0];
142     CeedScalar       xxsq[3][3];
143     for (int j = 0; j < 3; ++j)
144       for (int k = 0; k < 3; ++k) {
145         xxsq[j][k] = 0.;
146         for (int l = 0; l < 1; ++l) xxsq[j][k] += xx[j][l] * xx[k][l] / (sqrt(modxxsq) * modxxsq);
147       }
148 
149     const CeedScalar dxdxx[3][3] = {
150       {1. / sqrt(modxxsq) - xxsq[0][0], -xxsq[0][1],                     -xxsq[0][2]                    },
151       {-xxsq[1][0],                     1. / sqrt(modxxsq) - xxsq[1][1], -xxsq[1][2]                    },
152       {-xxsq[2][0],                     -xxsq[2][1],                     1. / sqrt(modxxsq) - xxsq[2][2]}
153     };
154 
155     CeedScalar dxdX[3][2];
156     for (int j = 0; j < 3; ++j)
157       for (int k = 0; k < 2; ++k) {
158         dxdX[j][k] = 0.;
159         for (int l = 0; l < 3; ++l) dxdX[j][k] += dxdxx[j][l] * dxxdX[l][k];
160       }
161     // J is given by the cross product of the columns of dxdX
162     const CeedScalar J[3][1] = {{dxdX[1][0] * dxdX[2][1] - dxdX[2][0] * dxdX[1][1]}, {dxdX[2][0] * dxdX[0][1] - dxdX[0][0] * dxdX[2][1]}, {dxdX[0][0] * dxdX[1][1] - dxdX[1][0] * dxdX[0][1]}};
163     // Use the magnitude of J as our detJ (volume scaling factor)
164     const CeedScalar modJ    = sqrt(J[0][0] * J[0][0] + J[1][0] * J[1][0] + J[2][0] * J[2][0]);
165     qdata[i + Q * 0]         = modJ * w[i];
166   }
167   return 0;
168 }
169 
170 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *ctx) {
171   DMPlexShape shape = DM_SHAPE_UNKNOWN;
172 
173   PetscFunctionBeginUser;
174   PetscOptionsBegin(comm, "", "libCEED Test Options", "DMPLEX");
175   PetscOptionsEnd();
176   PetscCall(PetscOptionsGetEnum(NULL, NULL, "-dm_plex_shape", DMPlexShapes, (PetscEnum *)&shape, NULL));
177   ctx->setupgeo      = NULL;
178   ctx->setupgeofname = NULL;
179   ctx->apply         = Mass;
180   ctx->applyfname    = Mass_loc;
181   ctx->areaExact     = 0.0;
182   switch (shape) {
183   case DM_SHAPE_BOX_SURFACE:
184     ctx->setupgeo      = SetupMassGeoCube;
185     ctx->setupgeofname = SetupMassGeoCube_loc;
186     ctx->areaExact     = 6.0;
187     break;
188   case DM_SHAPE_SPHERE:
189     ctx->setupgeo      = SetupMassGeoSphere;
190     ctx->setupgeofname = SetupMassGeoSphere_loc;
191     ctx->areaExact     = 4.0 * M_PI;
192     break;
193   default: break;
194   }
195   PetscFunctionReturn(0);
196 }
197 
198 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) {
199   PetscFunctionBegin;
200   PetscCall(DMCreate(comm, dm));
201   PetscCall(DMSetType(*dm, DMPLEX));
202   PetscCall(DMSetFromOptions(*dm));
203   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
204 #ifdef PETSC_HAVE_LIBCEED
205   {
206     Ceed        ceed;
207     const char *usedresource;
208 
209     PetscCall(DMGetCeed(*dm, &ceed));
210     PetscCall(CeedGetResource(ceed, &usedresource));
211     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)*dm), "libCEED Backend: %s\n", usedresource));
212   }
213 #endif
214   PetscFunctionReturn(0);
215 }
216 
217 static PetscErrorCode SetupDiscretization(DM dm) {
218   DM        cdm;
219   PetscFE   fe, cfe;
220   PetscInt  dim, cnc;
221   PetscBool simplex;
222 
223   PetscFunctionBeginUser;
224   PetscCall(DMGetDimension(dm, &dim));
225   PetscCall(DMPlexIsSimplex(dm, &simplex));
226   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, PETSC_DETERMINE, &fe));
227   PetscCall(PetscFESetName(fe, "indicator"));
228   PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
229   PetscCall(PetscFEDestroy(&fe));
230   PetscCall(DMCreateDS(dm));
231   PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
232   PetscCall(DMGetCoordinateDim(dm, &cnc));
233   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, cnc, simplex, NULL, PETSC_DETERMINE, &cfe));
234   PetscCall(DMProjectCoordinates(dm, cfe));
235   PetscCall(PetscFEDestroy(&cfe));
236   PetscCall(DMGetCoordinateDM(dm, &cdm));
237   PetscCall(DMPlexSetClosurePermutationTensor(cdm, PETSC_DETERMINE, NULL));
238   PetscFunctionReturn(0);
239 }
240 
241 static PetscErrorCode LibCeedSetupByDegree(DM dm, AppCtx *ctx, CeedData *data) {
242   PetscDS             ds;
243   PetscFE             fe, cfe;
244   Ceed                ceed;
245   CeedElemRestriction Erestrictx, Erestrictu, Erestrictq;
246   CeedQFunction       qf_setupgeo;
247   CeedOperator        op_setupgeo;
248   CeedVector          xcoord;
249   CeedBasis           basisu, basisx;
250   CeedInt             Nqdata = 1;
251   CeedInt             nqpts, nqptsx;
252   DM                  cdm;
253   Vec                 coords;
254   const PetscScalar  *coordArray;
255   PetscInt            dim, cdim, cStart, cEnd, Ncell;
256 
257   PetscFunctionBeginUser;
258   PetscCall(DMGetCeed(dm, &ceed));
259   PetscCall(DMGetDimension(dm, &dim));
260   PetscCall(DMGetCoordinateDim(dm, &cdim));
261   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
262   Ncell = cEnd - cStart;
263   // CEED bases
264   PetscCall(DMGetDS(dm, &ds));
265   PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&fe));
266   PetscCall(PetscFEGetCeedBasis(fe, &basisu));
267   PetscCall(DMGetCoordinateDM(dm, &cdm));
268   PetscCall(DMGetDS(cdm, &ds));
269   PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&cfe));
270   PetscCall(PetscFEGetCeedBasis(cfe, &basisx));
271 
272   PetscCall(DMPlexGetCeedRestriction(cdm, NULL, 0, 0, 0, &Erestrictx));
273   PetscCall(DMPlexGetCeedRestriction(dm, NULL, 0, 0, 0, &Erestrictu));
274   PetscCall(CeedBasisGetNumQuadraturePoints(basisu, &nqpts));
275   PetscCall(CeedBasisGetNumQuadraturePoints(basisx, &nqptsx));
276   PetscCheck(nqptsx == nqpts, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of qpoints for u %" PetscInt_FMT " != %" PetscInt_FMT " Number of qpoints for x", nqpts, nqptsx);
277   PetscCall(CeedElemRestrictionCreateStrided(ceed, Ncell, nqpts, Nqdata, Nqdata * Ncell * nqpts, CEED_STRIDES_BACKEND, &Erestrictq));
278 
279   PetscCall(DMGetCoordinatesLocal(dm, &coords));
280   PetscCall(VecGetArrayRead(coords, &coordArray));
281   PetscCall(CeedElemRestrictionCreateVector(Erestrictx, &xcoord, NULL));
282   PetscCall(CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_COPY_VALUES, (PetscScalar *)coordArray));
283   PetscCall(VecRestoreArrayRead(coords, &coordArray));
284 
285   // Create the vectors that will be needed in setup and apply
286   PetscCall(CeedElemRestrictionCreateVector(Erestrictu, &data->uceed, NULL));
287   PetscCall(CeedElemRestrictionCreateVector(Erestrictu, &data->vceed, NULL));
288   PetscCall(CeedElemRestrictionCreateVector(Erestrictq, &data->qdata, NULL));
289 
290   // Create the Q-function that builds the operator (i.e. computes its quadrature data) and set its context data
291   PetscCall(CeedQFunctionCreateInterior(ceed, 1, ctx->setupgeo, ctx->setupgeofname, &qf_setupgeo));
292   PetscCall(CeedQFunctionAddInput(qf_setupgeo, "x", cdim, CEED_EVAL_INTERP));
293   PetscCall(CeedQFunctionAddInput(qf_setupgeo, "dx", cdim * dim, CEED_EVAL_GRAD));
294   PetscCall(CeedQFunctionAddInput(qf_setupgeo, "weight", 1, CEED_EVAL_WEIGHT));
295   PetscCall(CeedQFunctionAddOutput(qf_setupgeo, "qdata", Nqdata, CEED_EVAL_NONE));
296 
297   // Set up the mass operator
298   PetscCall(CeedQFunctionCreateInterior(ceed, 1, ctx->apply, ctx->applyfname, &data->qf_apply));
299   PetscCall(CeedQFunctionAddInput(data->qf_apply, "u", 1, CEED_EVAL_INTERP));
300   PetscCall(CeedQFunctionAddInput(data->qf_apply, "qdata", Nqdata, CEED_EVAL_NONE));
301   PetscCall(CeedQFunctionAddOutput(data->qf_apply, "v", 1, CEED_EVAL_INTERP));
302 
303   // Create the operator that builds the quadrature data for the operator
304   PetscCall(CeedOperatorCreate(ceed, qf_setupgeo, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setupgeo));
305   PetscCall(CeedOperatorSetField(op_setupgeo, "x", Erestrictx, basisx, CEED_VECTOR_ACTIVE));
306   PetscCall(CeedOperatorSetField(op_setupgeo, "dx", Erestrictx, basisx, CEED_VECTOR_ACTIVE));
307   PetscCall(CeedOperatorSetField(op_setupgeo, "weight", CEED_ELEMRESTRICTION_NONE, basisx, CEED_VECTOR_NONE));
308   PetscCall(CeedOperatorSetField(op_setupgeo, "qdata", Erestrictq, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
309 
310   // Create the mass operator
311   PetscCall(CeedOperatorCreate(ceed, data->qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &data->op_apply));
312   PetscCall(CeedOperatorSetField(data->op_apply, "u", Erestrictu, basisu, CEED_VECTOR_ACTIVE));
313   PetscCall(CeedOperatorSetField(data->op_apply, "qdata", Erestrictq, CEED_BASIS_COLLOCATED, data->qdata));
314   PetscCall(CeedOperatorSetField(data->op_apply, "v", Erestrictu, basisu, CEED_VECTOR_ACTIVE));
315 
316   // Setup qdata
317   PetscCall(CeedOperatorApply(op_setupgeo, xcoord, data->qdata, CEED_REQUEST_IMMEDIATE));
318 
319   PetscCall(CeedElemRestrictionDestroy(&Erestrictq));
320   PetscCall(CeedQFunctionDestroy(&qf_setupgeo));
321   PetscCall(CeedOperatorDestroy(&op_setupgeo));
322   PetscCall(CeedVectorDestroy(&xcoord));
323   PetscFunctionReturn(0);
324 }
325 
326 int main(int argc, char **argv) {
327   MPI_Comm     comm;
328   DM           dm;
329   AppCtx       ctx;
330   Vec          U, Uloc, V, Vloc;
331   PetscScalar *v;
332   PetscScalar  area;
333   CeedData     ceeddata;
334 
335   PetscFunctionBeginUser;
336   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
337   comm = PETSC_COMM_WORLD;
338   PetscCall(ProcessOptions(comm, &ctx));
339   PetscCall(CreateMesh(comm, &ctx, &dm));
340   PetscCall(SetupDiscretization(dm));
341 
342   PetscCall(LibCeedSetupByDegree(dm, &ctx, &ceeddata));
343 
344   PetscCall(DMCreateGlobalVector(dm, &U));
345   PetscCall(DMCreateLocalVector(dm, &Uloc));
346   PetscCall(VecDuplicate(U, &V));
347   PetscCall(VecDuplicate(Uloc, &Vloc));
348 
349   /**/
350   PetscCall(VecSet(Uloc, 1.));
351   PetscCall(VecZeroEntries(V));
352   PetscCall(VecZeroEntries(Vloc));
353   PetscCall(VecGetArray(Vloc, &v));
354   PetscCall(CeedVectorSetArray(ceeddata.vceed, CEED_MEM_HOST, CEED_USE_POINTER, v));
355   PetscCall(CeedVectorSetValue(ceeddata.uceed, 1.0));
356   PetscCall(CeedOperatorApply(ceeddata.op_apply, ceeddata.uceed, ceeddata.vceed, CEED_REQUEST_IMMEDIATE));
357   PetscCall(CeedVectorTakeArray(ceeddata.vceed, CEED_MEM_HOST, NULL));
358   PetscCall(VecRestoreArray(Vloc, &v));
359   PetscCall(DMLocalToGlobalBegin(dm, Vloc, ADD_VALUES, V));
360   PetscCall(DMLocalToGlobalEnd(dm, Vloc, ADD_VALUES, V));
361 
362   PetscCall(VecSum(V, &area));
363   if (ctx.areaExact > 0.) {
364     PetscReal error = PetscAbsReal(area - ctx.areaExact);
365     PetscReal tol   = PETSC_SMALL;
366 
367     PetscCall(PetscPrintf(comm, "Exact mesh surface area    : % .*f\n", PetscAbsReal(ctx.areaExact - round(ctx.areaExact)) > 1E-15 ? 14 : 1, (double)ctx.areaExact));
368     PetscCall(PetscPrintf(comm, "Computed mesh surface area : % .*f\n", PetscAbsScalar(area - round(area)) > 1E-15 ? 14 : 1, (double)PetscRealPart(area)));
369     if (error > tol) {
370       PetscCall(PetscPrintf(comm, "Area error                 : % .14g\n", (double)error));
371     } else {
372       PetscCall(PetscPrintf(comm, "Area verifies!\n"));
373     }
374   }
375 
376   PetscCall(CeedDataDestroy(&ceeddata));
377   PetscCall(VecDestroy(&U));
378   PetscCall(VecDestroy(&Uloc));
379   PetscCall(VecDestroy(&V));
380   PetscCall(VecDestroy(&Vloc));
381   PetscCall(DMDestroy(&dm));
382   return PetscFinalize();
383 }
384 
385 /*TEST
386 
387   build:
388     requires: libceed
389 
390   testset:
391     args: -dm_plex_simplex 0 -petscspace_degree 3 -dm_view -dm_petscds_view \
392           -petscfe_default_quadrature_order 4 -coord_dm_default_quadrature_order 4
393 
394     test:
395       suffix: cube_3
396       args: -dm_plex_shape box_surface -dm_refine 2
397 
398     test:
399       suffix: cube_3_p4
400       nsize: 4
401       args: -petscpartitioner_type simple -dm_refine_pre 1 -dm_plex_shape box_surface -dm_refine 1
402 
403     test:
404       suffix: sphere_3
405       args: -dm_plex_shape sphere -dm_refine 3
406 
407     test:
408       suffix: sphere_3_p4
409       nsize: 4
410       args: -petscpartitioner_type simple -dm_refine_pre 1 -dm_plex_shape sphere -dm_refine 2
411 
412 TEST*/
413