xref: /petsc/src/dm/impls/plex/tests/ex42.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1f918ec44SMatthew G. Knepley static const char help[] = "Simple libCEED test to calculate surface area using 1^T M 1";
2f918ec44SMatthew G. Knepley 
3f918ec44SMatthew G. Knepley /*
4f918ec44SMatthew G. Knepley   This is a recreation of libCeed Example 2: https://libceed.readthedocs.io/en/latest/examples/ceed/
5f918ec44SMatthew G. Knepley */
6f918ec44SMatthew G. Knepley 
7f918ec44SMatthew G. Knepley #include <petscdmceed.h>
8f918ec44SMatthew G. Knepley #include <petscdmplexceed.h>
9f918ec44SMatthew G. Knepley #include <petscfeceed.h>
10f918ec44SMatthew G. Knepley #include <petscdmplex.h>
11f918ec44SMatthew G. Knepley #include <petscds.h>
12f918ec44SMatthew G. Knepley 
13f918ec44SMatthew G. Knepley typedef struct {
14f918ec44SMatthew G. Knepley   PetscReal         areaExact;
15f918ec44SMatthew G. Knepley   CeedQFunctionUser setupgeo, apply;
16f918ec44SMatthew G. Knepley   const char       *setupgeofname, *applyfname;
17f918ec44SMatthew G. Knepley } AppCtx;
18f918ec44SMatthew G. Knepley 
19f918ec44SMatthew G. Knepley typedef struct {
20f918ec44SMatthew G. Knepley   CeedQFunction qf_apply;
21f918ec44SMatthew G. Knepley   CeedOperator  op_apply;
22f918ec44SMatthew G. Knepley   CeedVector    qdata, uceed, vceed;
23f918ec44SMatthew G. Knepley } CeedData;
24f918ec44SMatthew G. Knepley 
25*9371c9d4SSatish Balay static PetscErrorCode CeedDataDestroy(CeedData *data) {
26f918ec44SMatthew G. Knepley   PetscFunctionBeginUser;
279566063dSJacob Faibussowitsch   PetscCall(CeedVectorDestroy(&data->qdata));
289566063dSJacob Faibussowitsch   PetscCall(CeedVectorDestroy(&data->uceed));
299566063dSJacob Faibussowitsch   PetscCall(CeedVectorDestroy(&data->vceed));
309566063dSJacob Faibussowitsch   PetscCall(CeedQFunctionDestroy(&data->qf_apply));
319566063dSJacob Faibussowitsch   PetscCall(CeedOperatorDestroy(&data->op_apply));
32f918ec44SMatthew G. Knepley   PetscFunctionReturn(0);
33f918ec44SMatthew G. Knepley }
34f918ec44SMatthew G. Knepley 
35*9371c9d4SSatish Balay CEED_QFUNCTION(Mass)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
36f918ec44SMatthew G. Knepley   const CeedScalar *u = in[0], *qdata = in[1];
37f918ec44SMatthew G. Knepley   CeedScalar       *v = out[0];
38f918ec44SMatthew G. Knepley 
39*9371c9d4SSatish Balay   CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i) v[i] = qdata[i] * u[i];
40f918ec44SMatthew G. Knepley 
41f918ec44SMatthew G. Knepley   return 0;
42f918ec44SMatthew G. Knepley }
43f918ec44SMatthew G. Knepley 
44f918ec44SMatthew G. Knepley /*
45f918ec44SMatthew G. Knepley // Reference (parent) 2D coordinates: X \in [-1, 1]^2
46f918ec44SMatthew G. Knepley //
47f918ec44SMatthew G. Knepley // Global physical coordinates given by the mesh (3D): xx \in [-l, l]^3
48f918ec44SMatthew G. Knepley //
49f918ec44SMatthew G. Knepley // Local physical coordinates on the manifold (2D): x \in [-l, l]^2
50f918ec44SMatthew G. Knepley //
51f918ec44SMatthew G. Knepley // Change of coordinates matrix computed by the library:
52f918ec44SMatthew G. Knepley //   (physical 3D coords relative to reference 2D coords)
53f918ec44SMatthew G. Knepley //   dxx_j/dX_i (indicial notation) [3 * 2]
54f918ec44SMatthew G. Knepley //
55f918ec44SMatthew G. Knepley // Change of coordinates x (physical 2D) relative to xx (phyisical 3D):
56f918ec44SMatthew G. Knepley //   dx_i/dxx_j (indicial notation) [2 * 3]
57f918ec44SMatthew G. Knepley //
58f918ec44SMatthew G. Knepley // Change of coordinates x (physical 2D) relative to X (reference 2D):
59f918ec44SMatthew G. Knepley //   (by chain rule)
60f918ec44SMatthew G. Knepley //   dx_i/dX_j = dx_i/dxx_k * dxx_k/dX_j
61f918ec44SMatthew G. Knepley //
62f918ec44SMatthew G. Knepley // The quadrature data is stored in the array qdata.
63f918ec44SMatthew G. Knepley //
64f918ec44SMatthew G. Knepley // We require the determinant of the Jacobian to properly compute integrals of the form: int(u v)
65f918ec44SMatthew G. Knepley //
66f918ec44SMatthew G. Knepley // Qdata: w * det(dx_i/dX_j)
67f918ec44SMatthew G. Knepley */
68*9371c9d4SSatish Balay CEED_QFUNCTION(SetupMassGeoCube)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
69f918ec44SMatthew G. Knepley   const CeedScalar *J = in[1], *w = in[2];
70f918ec44SMatthew G. Knepley   CeedScalar       *qdata = out[0];
71f918ec44SMatthew G. Knepley 
72*9371c9d4SSatish Balay   CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i) {
73f918ec44SMatthew G. Knepley     // Read dxxdX Jacobian entries, stored as [[0 3], [1 4], [2 5]]
74*9371c9d4SSatish Balay     const CeedScalar dxxdX[3][2] = {
75*9371c9d4SSatish Balay       {J[i + Q * 0], J[i + Q * 3]},
76f918ec44SMatthew G. Knepley       {J[i + Q * 1], J[i + Q * 4]},
77*9371c9d4SSatish Balay       {J[i + Q * 2], J[i + Q * 5]}
78*9371c9d4SSatish Balay     };
79f918ec44SMatthew G. Knepley     // Modulus of dxxdX column vectors
80f918ec44SMatthew G. Knepley     const CeedScalar modg1       = PetscSqrtReal(dxxdX[0][0] * dxxdX[0][0] + dxxdX[1][0] * dxxdX[1][0] + dxxdX[2][0] * dxxdX[2][0]);
81f918ec44SMatthew G. Knepley     const CeedScalar modg2       = PetscSqrtReal(dxxdX[0][1] * dxxdX[0][1] + dxxdX[1][1] * dxxdX[1][1] + dxxdX[2][1] * dxxdX[2][1]);
82f918ec44SMatthew G. Knepley     // Use normalized column vectors of dxxdX as rows of dxdxx
83*9371c9d4SSatish Balay     const CeedScalar dxdxx[2][3] = {
84*9371c9d4SSatish Balay       {dxxdX[0][0] / modg1, dxxdX[1][0] / modg1, dxxdX[2][0] / modg1},
85*9371c9d4SSatish Balay       {dxxdX[0][1] / modg2, dxxdX[1][1] / modg2, dxxdX[2][1] / modg2}
86*9371c9d4SSatish Balay     };
87f918ec44SMatthew G. Knepley 
88f918ec44SMatthew G. Knepley     CeedScalar dxdX[2][2];
89f918ec44SMatthew G. Knepley     for (int j = 0; j < 2; ++j)
90f918ec44SMatthew G. Knepley       for (int k = 0; k < 2; ++k) {
91f918ec44SMatthew G. Knepley         dxdX[j][k] = 0;
92*9371c9d4SSatish Balay         for (int l = 0; l < 3; ++l) dxdX[j][k] += dxdxx[j][l] * dxxdX[l][k];
93f918ec44SMatthew G. Knepley       }
94f918ec44SMatthew G. Knepley     qdata[i + Q * 0] = (dxdX[0][0] * dxdX[1][1] - dxdX[1][0] * dxdX[0][1]) * w[i]; /* det J * weight */
95f918ec44SMatthew G. Knepley   }
96f918ec44SMatthew G. Knepley   return 0;
97f918ec44SMatthew G. Knepley }
98f918ec44SMatthew G. Knepley 
99f918ec44SMatthew G. Knepley /*
100f918ec44SMatthew G. Knepley // Reference (parent) 2D coordinates: X \in [-1, 1]^2
101f918ec44SMatthew G. Knepley //
102f918ec44SMatthew G. Knepley // Global 3D physical coordinates given by the mesh: xx \in [-R, R]^3
103f918ec44SMatthew G. Knepley //   with R radius of the sphere
104f918ec44SMatthew G. Knepley //
105f918ec44SMatthew G. Knepley // Local 3D physical coordinates on the 2D manifold: x \in [-l, l]^3
106f918ec44SMatthew G. Knepley //   with l half edge of the cube inscribed in the sphere
107f918ec44SMatthew G. Knepley //
108f918ec44SMatthew G. Knepley // Change of coordinates matrix computed by the library:
109f918ec44SMatthew G. Knepley //   (physical 3D coords relative to reference 2D coords)
110f918ec44SMatthew G. Knepley //   dxx_j/dX_i (indicial notation) [3 * 2]
111f918ec44SMatthew G. Knepley //
112f918ec44SMatthew G. Knepley // Change of coordinates x (on the 2D manifold) relative to xx (phyisical 3D):
113f918ec44SMatthew G. Knepley //   dx_i/dxx_j (indicial notation) [3 * 3]
114f918ec44SMatthew G. Knepley //
115f918ec44SMatthew G. Knepley // Change of coordinates x (on the 2D manifold) relative to X (reference 2D):
116f918ec44SMatthew G. Knepley //   (by chain rule)
117f918ec44SMatthew G. Knepley //   dx_i/dX_j = dx_i/dxx_k * dxx_k/dX_j [3 * 2]
118f918ec44SMatthew G. Knepley //
119f918ec44SMatthew G. Knepley // modJ is given by the magnitude of the cross product of the columns of dx_i/dX_j
120f918ec44SMatthew G. Knepley //
121f918ec44SMatthew G. Knepley // The quadrature data is stored in the array qdata.
122f918ec44SMatthew G. Knepley //
123f918ec44SMatthew G. Knepley // We require the determinant of the Jacobian to properly compute integrals of
124f918ec44SMatthew G. Knepley //   the form: int(u v)
125f918ec44SMatthew G. Knepley //
126f918ec44SMatthew G. Knepley // Qdata: modJ * w
127f918ec44SMatthew G. Knepley */
128f918ec44SMatthew G. Knepley CEED_QFUNCTION(SetupMassGeoSphere)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
129f918ec44SMatthew G. Knepley   const CeedScalar *X = in[0], *J = in[1], *w = in[2];
130f918ec44SMatthew G. Knepley   CeedScalar       *qdata = out[0];
131f918ec44SMatthew G. Knepley 
132*9371c9d4SSatish Balay   CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i) {
133f918ec44SMatthew G. Knepley     const CeedScalar xx[3][1]    = {{X[i + 0 * Q]}, {X[i + 1 * Q]}, {X[i + 2 * Q]}};
134f918ec44SMatthew G. Knepley     // Read dxxdX Jacobian entries, stored as [[0 3], [1 4], [2 5]]
135*9371c9d4SSatish Balay     const CeedScalar dxxdX[3][2] = {
136*9371c9d4SSatish Balay       {J[i + Q * 0], J[i + Q * 3]},
137f918ec44SMatthew G. Knepley       {J[i + Q * 1], J[i + Q * 4]},
138*9371c9d4SSatish Balay       {J[i + Q * 2], J[i + Q * 5]}
139*9371c9d4SSatish Balay     };
140f918ec44SMatthew G. Knepley     // Setup
141f918ec44SMatthew G. Knepley     const CeedScalar modxxsq = xx[0][0] * xx[0][0] + xx[1][0] * xx[1][0] + xx[2][0] * xx[2][0];
142f918ec44SMatthew G. Knepley     CeedScalar       xxsq[3][3];
143f918ec44SMatthew G. Knepley     for (int j = 0; j < 3; ++j)
144f918ec44SMatthew G. Knepley       for (int k = 0; k < 3; ++k) {
145f918ec44SMatthew G. Knepley         xxsq[j][k] = 0.;
146*9371c9d4SSatish Balay         for (int l = 0; l < 1; ++l) xxsq[j][k] += xx[j][l] * xx[k][l] / (sqrt(modxxsq) * modxxsq);
147f918ec44SMatthew G. Knepley       }
148f918ec44SMatthew G. Knepley 
149*9371c9d4SSatish Balay     const CeedScalar dxdxx[3][3] = {
150*9371c9d4SSatish Balay       {1. / sqrt(modxxsq) - xxsq[0][0], -xxsq[0][1],                     -xxsq[0][2]                    },
151f918ec44SMatthew G. Knepley       {-xxsq[1][0],                     1. / sqrt(modxxsq) - xxsq[1][1], -xxsq[1][2]                    },
152*9371c9d4SSatish Balay       {-xxsq[2][0],                     -xxsq[2][1],                     1. / sqrt(modxxsq) - xxsq[2][2]}
153*9371c9d4SSatish Balay     };
154f918ec44SMatthew G. Knepley 
155f918ec44SMatthew G. Knepley     CeedScalar dxdX[3][2];
156f918ec44SMatthew G. Knepley     for (int j = 0; j < 3; ++j)
157f918ec44SMatthew G. Knepley       for (int k = 0; k < 2; ++k) {
158f918ec44SMatthew G. Knepley         dxdX[j][k] = 0.;
159*9371c9d4SSatish Balay         for (int l = 0; l < 3; ++l) dxdX[j][k] += dxdxx[j][l] * dxxdX[l][k];
160f918ec44SMatthew G. Knepley       }
161f918ec44SMatthew G. Knepley     // J is given by the cross product of the columns of dxdX
162*9371c9d4SSatish Balay     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]}};
163f918ec44SMatthew G. Knepley     // Use the magnitude of J as our detJ (volume scaling factor)
164f918ec44SMatthew G. Knepley     const CeedScalar modJ    = sqrt(J[0][0] * J[0][0] + J[1][0] * J[1][0] + J[2][0] * J[2][0]);
165f918ec44SMatthew G. Knepley     qdata[i + Q * 0]         = modJ * w[i];
166f918ec44SMatthew G. Knepley   }
167f918ec44SMatthew G. Knepley   return 0;
168f918ec44SMatthew G. Knepley }
169f918ec44SMatthew G. Knepley 
170*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *ctx) {
171f918ec44SMatthew G. Knepley   DMPlexShape shape = DM_SHAPE_UNKNOWN;
172f918ec44SMatthew G. Knepley 
173f918ec44SMatthew G. Knepley   PetscFunctionBeginUser;
174d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "libCEED Test Options", "DMPLEX");
175d0609cedSBarry Smith   PetscOptionsEnd();
1769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetEnum(NULL, NULL, "-dm_plex_shape", DMPlexShapes, (PetscEnum *)&shape, NULL));
17730602db0SMatthew G. Knepley   ctx->setupgeo      = NULL;
17830602db0SMatthew G. Knepley   ctx->setupgeofname = NULL;
17930602db0SMatthew G. Knepley   ctx->apply         = Mass;
18030602db0SMatthew G. Knepley   ctx->applyfname    = Mass_loc;
18130602db0SMatthew G. Knepley   ctx->areaExact     = 0.0;
182f918ec44SMatthew G. Knepley   switch (shape) {
183f918ec44SMatthew G. Knepley   case DM_SHAPE_BOX_SURFACE:
184f918ec44SMatthew G. Knepley     ctx->setupgeo      = SetupMassGeoCube;
185f918ec44SMatthew G. Knepley     ctx->setupgeofname = SetupMassGeoCube_loc;
186f918ec44SMatthew G. Knepley     ctx->areaExact     = 6.0;
187f918ec44SMatthew G. Knepley     break;
188f918ec44SMatthew G. Knepley   case DM_SHAPE_SPHERE:
189f918ec44SMatthew G. Knepley     ctx->setupgeo      = SetupMassGeoSphere;
190f918ec44SMatthew G. Knepley     ctx->setupgeofname = SetupMassGeoSphere_loc;
191f918ec44SMatthew G. Knepley     ctx->areaExact     = 4.0 * M_PI;
192f918ec44SMatthew G. Knepley     break;
19330602db0SMatthew G. Knepley   default: break;
194f918ec44SMatthew G. Knepley   }
195f918ec44SMatthew G. Knepley   PetscFunctionReturn(0);
196f918ec44SMatthew G. Knepley }
197f918ec44SMatthew G. Knepley 
198*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) {
199f918ec44SMatthew G. Knepley   PetscFunctionBegin;
2009566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
2019566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
2029566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
2039566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
20430602db0SMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED
20530602db0SMatthew G. Knepley   {
20630602db0SMatthew G. Knepley     Ceed        ceed;
20730602db0SMatthew G. Knepley     const char *usedresource;
20830602db0SMatthew G. Knepley 
2099566063dSJacob Faibussowitsch     PetscCall(DMGetCeed(*dm, &ceed));
2109566063dSJacob Faibussowitsch     PetscCall(CeedGetResource(ceed, &usedresource));
2119566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)*dm), "libCEED Backend: %s\n", usedresource));
21230602db0SMatthew G. Knepley   }
21330602db0SMatthew G. Knepley #endif
214f918ec44SMatthew G. Knepley   PetscFunctionReturn(0);
215f918ec44SMatthew G. Knepley }
216f918ec44SMatthew G. Knepley 
217*9371c9d4SSatish Balay static PetscErrorCode SetupDiscretization(DM dm) {
218f918ec44SMatthew G. Knepley   DM        cdm;
219a2c9b50fSJeremy L Thompson   PetscFE   fe, cfe;
220a2c9b50fSJeremy L Thompson   PetscInt  dim, cnc;
221f918ec44SMatthew G. Knepley   PetscBool simplex;
222f918ec44SMatthew G. Knepley 
223f918ec44SMatthew G. Knepley   PetscFunctionBeginUser;
2249566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
2259566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(dm, &simplex));
2269566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, PETSC_DETERMINE, &fe));
2279566063dSJacob Faibussowitsch   PetscCall(PetscFESetName(fe, "indicator"));
2289566063dSJacob Faibussowitsch   PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
2299566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
2309566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
2319566063dSJacob Faibussowitsch   PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
2329566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &cnc));
2339566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, cnc, simplex, NULL, PETSC_DETERMINE, &cfe));
2349566063dSJacob Faibussowitsch   PetscCall(DMProjectCoordinates(dm, cfe));
2359566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&cfe));
2369566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
2379566063dSJacob Faibussowitsch   PetscCall(DMPlexSetClosurePermutationTensor(cdm, PETSC_DETERMINE, NULL));
238f918ec44SMatthew G. Knepley   PetscFunctionReturn(0);
239f918ec44SMatthew G. Knepley }
240f918ec44SMatthew G. Knepley 
241*9371c9d4SSatish Balay static PetscErrorCode LibCeedSetupByDegree(DM dm, AppCtx *ctx, CeedData *data) {
242f918ec44SMatthew G. Knepley   PetscDS             ds;
243f918ec44SMatthew G. Knepley   PetscFE             fe, cfe;
244f918ec44SMatthew G. Knepley   Ceed                ceed;
245f918ec44SMatthew G. Knepley   CeedElemRestriction Erestrictx, Erestrictu, Erestrictq;
246f918ec44SMatthew G. Knepley   CeedQFunction       qf_setupgeo;
247f918ec44SMatthew G. Knepley   CeedOperator        op_setupgeo;
248f918ec44SMatthew G. Knepley   CeedVector          xcoord;
249f918ec44SMatthew G. Knepley   CeedBasis           basisu, basisx;
250f918ec44SMatthew G. Knepley   CeedInt             Nqdata = 1;
251f918ec44SMatthew G. Knepley   CeedInt             nqpts, nqptsx;
252f918ec44SMatthew G. Knepley   DM                  cdm;
253f918ec44SMatthew G. Knepley   Vec                 coords;
254f918ec44SMatthew G. Knepley   const PetscScalar  *coordArray;
255f918ec44SMatthew G. Knepley   PetscInt            dim, cdim, cStart, cEnd, Ncell;
256f918ec44SMatthew G. Knepley 
257f918ec44SMatthew G. Knepley   PetscFunctionBeginUser;
2589566063dSJacob Faibussowitsch   PetscCall(DMGetCeed(dm, &ceed));
2599566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
2609566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &cdim));
2619566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
262f918ec44SMatthew G. Knepley   Ncell = cEnd - cStart;
263f918ec44SMatthew G. Knepley   // CEED bases
2649566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
2659566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&fe));
2669566063dSJacob Faibussowitsch   PetscCall(PetscFEGetCeedBasis(fe, &basisu));
2679566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
2689566063dSJacob Faibussowitsch   PetscCall(DMGetDS(cdm, &ds));
2699566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&cfe));
2709566063dSJacob Faibussowitsch   PetscCall(PetscFEGetCeedBasis(cfe, &basisx));
271f918ec44SMatthew G. Knepley 
2729566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCeedRestriction(cdm, NULL, 0, 0, 0, &Erestrictx));
2739566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCeedRestriction(dm, NULL, 0, 0, 0, &Erestrictu));
2749566063dSJacob Faibussowitsch   PetscCall(CeedBasisGetNumQuadraturePoints(basisu, &nqpts));
2759566063dSJacob Faibussowitsch   PetscCall(CeedBasisGetNumQuadraturePoints(basisx, &nqptsx));
27663a3b9bcSJacob Faibussowitsch   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);
2779566063dSJacob Faibussowitsch   PetscCall(CeedElemRestrictionCreateStrided(ceed, Ncell, nqpts, Nqdata, Nqdata * Ncell * nqpts, CEED_STRIDES_BACKEND, &Erestrictq));
278f918ec44SMatthew G. Knepley 
2799566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coords));
2809566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(coords, &coordArray));
2819566063dSJacob Faibussowitsch   PetscCall(CeedElemRestrictionCreateVector(Erestrictx, &xcoord, NULL));
2829566063dSJacob Faibussowitsch   PetscCall(CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_COPY_VALUES, (PetscScalar *)coordArray));
2839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(coords, &coordArray));
284f918ec44SMatthew G. Knepley 
285f918ec44SMatthew G. Knepley   // Create the vectors that will be needed in setup and apply
2869566063dSJacob Faibussowitsch   PetscCall(CeedElemRestrictionCreateVector(Erestrictu, &data->uceed, NULL));
2879566063dSJacob Faibussowitsch   PetscCall(CeedElemRestrictionCreateVector(Erestrictu, &data->vceed, NULL));
2889566063dSJacob Faibussowitsch   PetscCall(CeedElemRestrictionCreateVector(Erestrictq, &data->qdata, NULL));
289f918ec44SMatthew G. Knepley 
290f918ec44SMatthew G. Knepley   // Create the Q-function that builds the operator (i.e. computes its quadrature data) and set its context data
2919566063dSJacob Faibussowitsch   PetscCall(CeedQFunctionCreateInterior(ceed, 1, ctx->setupgeo, ctx->setupgeofname, &qf_setupgeo));
2929566063dSJacob Faibussowitsch   PetscCall(CeedQFunctionAddInput(qf_setupgeo, "x", cdim, CEED_EVAL_INTERP));
2939566063dSJacob Faibussowitsch   PetscCall(CeedQFunctionAddInput(qf_setupgeo, "dx", cdim * dim, CEED_EVAL_GRAD));
2949566063dSJacob Faibussowitsch   PetscCall(CeedQFunctionAddInput(qf_setupgeo, "weight", 1, CEED_EVAL_WEIGHT));
2959566063dSJacob Faibussowitsch   PetscCall(CeedQFunctionAddOutput(qf_setupgeo, "qdata", Nqdata, CEED_EVAL_NONE));
296f918ec44SMatthew G. Knepley 
297f918ec44SMatthew G. Knepley   // Set up the mass operator
2989566063dSJacob Faibussowitsch   PetscCall(CeedQFunctionCreateInterior(ceed, 1, ctx->apply, ctx->applyfname, &data->qf_apply));
2999566063dSJacob Faibussowitsch   PetscCall(CeedQFunctionAddInput(data->qf_apply, "u", 1, CEED_EVAL_INTERP));
3009566063dSJacob Faibussowitsch   PetscCall(CeedQFunctionAddInput(data->qf_apply, "qdata", Nqdata, CEED_EVAL_NONE));
3019566063dSJacob Faibussowitsch   PetscCall(CeedQFunctionAddOutput(data->qf_apply, "v", 1, CEED_EVAL_INTERP));
302f918ec44SMatthew G. Knepley 
303f918ec44SMatthew G. Knepley   // Create the operator that builds the quadrature data for the operator
3049566063dSJacob Faibussowitsch   PetscCall(CeedOperatorCreate(ceed, qf_setupgeo, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setupgeo));
3059566063dSJacob Faibussowitsch   PetscCall(CeedOperatorSetField(op_setupgeo, "x", Erestrictx, basisx, CEED_VECTOR_ACTIVE));
3069566063dSJacob Faibussowitsch   PetscCall(CeedOperatorSetField(op_setupgeo, "dx", Erestrictx, basisx, CEED_VECTOR_ACTIVE));
3079566063dSJacob Faibussowitsch   PetscCall(CeedOperatorSetField(op_setupgeo, "weight", CEED_ELEMRESTRICTION_NONE, basisx, CEED_VECTOR_NONE));
3089566063dSJacob Faibussowitsch   PetscCall(CeedOperatorSetField(op_setupgeo, "qdata", Erestrictq, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
309f918ec44SMatthew G. Knepley 
310f918ec44SMatthew G. Knepley   // Create the mass operator
3119566063dSJacob Faibussowitsch   PetscCall(CeedOperatorCreate(ceed, data->qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &data->op_apply));
3129566063dSJacob Faibussowitsch   PetscCall(CeedOperatorSetField(data->op_apply, "u", Erestrictu, basisu, CEED_VECTOR_ACTIVE));
3139566063dSJacob Faibussowitsch   PetscCall(CeedOperatorSetField(data->op_apply, "qdata", Erestrictq, CEED_BASIS_COLLOCATED, data->qdata));
3149566063dSJacob Faibussowitsch   PetscCall(CeedOperatorSetField(data->op_apply, "v", Erestrictu, basisu, CEED_VECTOR_ACTIVE));
315f918ec44SMatthew G. Knepley 
316f918ec44SMatthew G. Knepley   // Setup qdata
3179566063dSJacob Faibussowitsch   PetscCall(CeedOperatorApply(op_setupgeo, xcoord, data->qdata, CEED_REQUEST_IMMEDIATE));
318f918ec44SMatthew G. Knepley 
3199566063dSJacob Faibussowitsch   PetscCall(CeedElemRestrictionDestroy(&Erestrictq));
3209566063dSJacob Faibussowitsch   PetscCall(CeedQFunctionDestroy(&qf_setupgeo));
3219566063dSJacob Faibussowitsch   PetscCall(CeedOperatorDestroy(&op_setupgeo));
3229566063dSJacob Faibussowitsch   PetscCall(CeedVectorDestroy(&xcoord));
323f918ec44SMatthew G. Knepley   PetscFunctionReturn(0);
324f918ec44SMatthew G. Knepley }
325f918ec44SMatthew G. Knepley 
326*9371c9d4SSatish Balay int main(int argc, char **argv) {
327f918ec44SMatthew G. Knepley   MPI_Comm     comm;
328f918ec44SMatthew G. Knepley   DM           dm;
329f918ec44SMatthew G. Knepley   AppCtx       ctx;
330f918ec44SMatthew G. Knepley   Vec          U, Uloc, V, Vloc;
331f918ec44SMatthew G. Knepley   PetscScalar *v;
332f918ec44SMatthew G. Knepley   PetscScalar  area;
333f918ec44SMatthew G. Knepley   CeedData     ceeddata;
334f918ec44SMatthew G. Knepley 
335327415f7SBarry Smith   PetscFunctionBeginUser;
3369566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
337f918ec44SMatthew G. Knepley   comm = PETSC_COMM_WORLD;
3389566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(comm, &ctx));
3399566063dSJacob Faibussowitsch   PetscCall(CreateMesh(comm, &ctx, &dm));
3409566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm));
341f918ec44SMatthew G. Knepley 
3429566063dSJacob Faibussowitsch   PetscCall(LibCeedSetupByDegree(dm, &ctx, &ceeddata));
343f918ec44SMatthew G. Knepley 
3449566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &U));
3459566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(dm, &Uloc));
3469566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U, &V));
3479566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(Uloc, &Vloc));
348f918ec44SMatthew G. Knepley 
34930602db0SMatthew G. Knepley   /**/
3509566063dSJacob Faibussowitsch   PetscCall(VecSet(Uloc, 1.));
3519566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(V));
3529566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(Vloc));
3539566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Vloc, &v));
3549566063dSJacob Faibussowitsch   PetscCall(CeedVectorSetArray(ceeddata.vceed, CEED_MEM_HOST, CEED_USE_POINTER, v));
3559566063dSJacob Faibussowitsch   PetscCall(CeedVectorSetValue(ceeddata.uceed, 1.0));
3569566063dSJacob Faibussowitsch   PetscCall(CeedOperatorApply(ceeddata.op_apply, ceeddata.uceed, ceeddata.vceed, CEED_REQUEST_IMMEDIATE));
3579566063dSJacob Faibussowitsch   PetscCall(CeedVectorTakeArray(ceeddata.vceed, CEED_MEM_HOST, NULL));
3589566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Vloc, &v));
3599566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(dm, Vloc, ADD_VALUES, V));
3609566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(dm, Vloc, ADD_VALUES, V));
361f918ec44SMatthew G. Knepley 
3629566063dSJacob Faibussowitsch   PetscCall(VecSum(V, &area));
363f918ec44SMatthew G. Knepley   if (ctx.areaExact > 0.) {
364f918ec44SMatthew G. Knepley     PetscReal error = PetscAbsReal(area - ctx.areaExact);
365f918ec44SMatthew G. Knepley     PetscReal tol   = PETSC_SMALL;
366f918ec44SMatthew G. Knepley 
36763a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Exact mesh surface area    : % .*f\n", PetscAbsReal(ctx.areaExact - round(ctx.areaExact)) > 1E-15 ? 14 : 1, (double)ctx.areaExact));
36863a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Computed mesh surface area : % .*f\n", PetscAbsScalar(area - round(area)) > 1E-15 ? 14 : 1, (double)PetscRealPart(area)));
369f918ec44SMatthew G. Knepley     if (error > tol) {
37063a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Area error                 : % .14g\n", (double)error));
371f918ec44SMatthew G. Knepley     } else {
37263a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Area verifies!\n"));
373f918ec44SMatthew G. Knepley     }
374f918ec44SMatthew G. Knepley   }
375f918ec44SMatthew G. Knepley 
3769566063dSJacob Faibussowitsch   PetscCall(CeedDataDestroy(&ceeddata));
3779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
3789566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Uloc));
3799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&V));
3809566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Vloc));
3819566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
382f918ec44SMatthew G. Knepley   return PetscFinalize();
383f918ec44SMatthew G. Knepley }
384f918ec44SMatthew G. Knepley 
385f918ec44SMatthew G. Knepley /*TEST
386f918ec44SMatthew G. Knepley 
387f918ec44SMatthew G. Knepley   build:
388f918ec44SMatthew G. Knepley     requires: libceed
389f918ec44SMatthew G. Knepley 
390f918ec44SMatthew G. Knepley   testset:
391e600fa54SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 3 -dm_view -dm_petscds_view \
392f918ec44SMatthew G. Knepley           -petscfe_default_quadrature_order 4 -coord_dm_default_quadrature_order 4
393f918ec44SMatthew G. Knepley 
394f918ec44SMatthew G. Knepley     test:
395f918ec44SMatthew G. Knepley       suffix: cube_3
396f918ec44SMatthew G. Knepley       args: -dm_plex_shape box_surface -dm_refine 2
397f918ec44SMatthew G. Knepley 
398f918ec44SMatthew G. Knepley     test:
399f918ec44SMatthew G. Knepley       suffix: cube_3_p4
400f918ec44SMatthew G. Knepley       nsize: 4
401fb796b39SJed Brown       args: -petscpartitioner_type simple -dm_refine_pre 1 -dm_plex_shape box_surface -dm_refine 1
402f918ec44SMatthew G. Knepley 
403f918ec44SMatthew G. Knepley     test:
404f918ec44SMatthew G. Knepley       suffix: sphere_3
405f918ec44SMatthew G. Knepley       args: -dm_plex_shape sphere -dm_refine 3
406f918ec44SMatthew G. Knepley 
407f918ec44SMatthew G. Knepley     test:
408f918ec44SMatthew G. Knepley       suffix: sphere_3_p4
409f918ec44SMatthew G. Knepley       nsize: 4
410fb796b39SJed Brown       args: -petscpartitioner_type simple -dm_refine_pre 1 -dm_plex_shape sphere -dm_refine 2
411f918ec44SMatthew G. Knepley 
412f918ec44SMatthew G. Knepley TEST*/
413