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