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
CeedDataDestroy(CeedData * data)25d71ae5a4SJacob Faibussowitsch static PetscErrorCode CeedDataDestroy(CeedData *data)
26d71ae5a4SJacob Faibussowitsch {
27f918ec44SMatthew G. Knepley PetscFunctionBeginUser;
289566063dSJacob Faibussowitsch PetscCall(CeedVectorDestroy(&data->qdata));
299566063dSJacob Faibussowitsch PetscCall(CeedVectorDestroy(&data->uceed));
309566063dSJacob Faibussowitsch PetscCall(CeedVectorDestroy(&data->vceed));
319566063dSJacob Faibussowitsch PetscCall(CeedQFunctionDestroy(&data->qf_apply));
329566063dSJacob Faibussowitsch PetscCall(CeedOperatorDestroy(&data->op_apply));
333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
34f918ec44SMatthew G. Knepley }
35f918ec44SMatthew G. Knepley
Mass(PetscCtx ctx,const CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)36*2a8381b2SBarry Smith CEED_QFUNCTION(Mass)(PetscCtx ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)
37d71ae5a4SJacob Faibussowitsch {
38f918ec44SMatthew G. Knepley const CeedScalar *u = in[0], *qdata = in[1];
39f918ec44SMatthew G. Knepley CeedScalar *v = out[0];
40f918ec44SMatthew G. Knepley
419371c9d4SSatish Balay CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i) v[i] = qdata[i] * u[i];
42f918ec44SMatthew G. Knepley
43f918ec44SMatthew G. Knepley return 0;
44f918ec44SMatthew G. Knepley }
45f918ec44SMatthew G. Knepley
46f918ec44SMatthew G. Knepley /*
47f918ec44SMatthew G. Knepley // Reference (parent) 2D coordinates: X \in [-1, 1]^2
48f918ec44SMatthew G. Knepley //
49f918ec44SMatthew G. Knepley // Global physical coordinates given by the mesh (3D): xx \in [-l, l]^3
50f918ec44SMatthew G. Knepley //
51f918ec44SMatthew G. Knepley // Local physical coordinates on the manifold (2D): x \in [-l, l]^2
52f918ec44SMatthew G. Knepley //
53f918ec44SMatthew G. Knepley // Change of coordinates matrix computed by the library:
54f918ec44SMatthew G. Knepley // (physical 3D coords relative to reference 2D coords)
55f918ec44SMatthew G. Knepley // dxx_j/dX_i (indicial notation) [3 * 2]
56f918ec44SMatthew G. Knepley //
57f918ec44SMatthew G. Knepley // Change of coordinates x (physical 2D) relative to xx (phyisical 3D):
58f918ec44SMatthew G. Knepley // dx_i/dxx_j (indicial notation) [2 * 3]
59f918ec44SMatthew G. Knepley //
60f918ec44SMatthew G. Knepley // Change of coordinates x (physical 2D) relative to X (reference 2D):
61f918ec44SMatthew G. Knepley // (by chain rule)
62f918ec44SMatthew G. Knepley // dx_i/dX_j = dx_i/dxx_k * dxx_k/dX_j
63f918ec44SMatthew G. Knepley //
64f918ec44SMatthew G. Knepley // The quadrature data is stored in the array qdata.
65f918ec44SMatthew G. Knepley //
66f918ec44SMatthew G. Knepley // We require the determinant of the Jacobian to properly compute integrals of the form: int(u v)
67f918ec44SMatthew G. Knepley //
68f918ec44SMatthew G. Knepley // Qdata: w * det(dx_i/dX_j)
69f918ec44SMatthew G. Knepley */
SetupMassGeoCube(PetscCtx ctx,const CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)70*2a8381b2SBarry Smith CEED_QFUNCTION(SetupMassGeoCube)(PetscCtx ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)
71d71ae5a4SJacob Faibussowitsch {
72f918ec44SMatthew G. Knepley const CeedScalar *J = in[1], *w = in[2];
73f918ec44SMatthew G. Knepley CeedScalar *qdata = out[0];
74f918ec44SMatthew G. Knepley
75d71ae5a4SJacob Faibussowitsch CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i)
76d71ae5a4SJacob Faibussowitsch {
77f918ec44SMatthew G. Knepley // Read dxxdX Jacobian entries, stored as [[0 3], [1 4], [2 5]]
789371c9d4SSatish Balay const CeedScalar dxxdX[3][2] = {
799371c9d4SSatish Balay {J[i + Q * 0], J[i + Q * 3]},
80f918ec44SMatthew G. Knepley {J[i + Q * 1], J[i + Q * 4]},
819371c9d4SSatish Balay {J[i + Q * 2], J[i + Q * 5]}
829371c9d4SSatish Balay };
83f918ec44SMatthew G. Knepley // Modulus of dxxdX column vectors
84f918ec44SMatthew G. Knepley const CeedScalar modg1 = PetscSqrtReal(dxxdX[0][0] * dxxdX[0][0] + dxxdX[1][0] * dxxdX[1][0] + dxxdX[2][0] * dxxdX[2][0]);
85f918ec44SMatthew G. Knepley const CeedScalar modg2 = PetscSqrtReal(dxxdX[0][1] * dxxdX[0][1] + dxxdX[1][1] * dxxdX[1][1] + dxxdX[2][1] * dxxdX[2][1]);
86f918ec44SMatthew G. Knepley // Use normalized column vectors of dxxdX as rows of dxdxx
879371c9d4SSatish Balay const CeedScalar dxdxx[2][3] = {
889371c9d4SSatish Balay {dxxdX[0][0] / modg1, dxxdX[1][0] / modg1, dxxdX[2][0] / modg1},
899371c9d4SSatish Balay {dxxdX[0][1] / modg2, dxxdX[1][1] / modg2, dxxdX[2][1] / modg2}
909371c9d4SSatish Balay };
91f918ec44SMatthew G. Knepley
92f918ec44SMatthew G. Knepley CeedScalar dxdX[2][2];
93f918ec44SMatthew G. Knepley for (int j = 0; j < 2; ++j)
94f918ec44SMatthew G. Knepley for (int k = 0; k < 2; ++k) {
95f918ec44SMatthew G. Knepley dxdX[j][k] = 0;
969371c9d4SSatish Balay for (int l = 0; l < 3; ++l) dxdX[j][k] += dxdxx[j][l] * dxxdX[l][k];
97f918ec44SMatthew G. Knepley }
98f918ec44SMatthew G. Knepley qdata[i + Q * 0] = (dxdX[0][0] * dxdX[1][1] - dxdX[1][0] * dxdX[0][1]) * w[i]; /* det J * weight */
99f918ec44SMatthew G. Knepley }
100f918ec44SMatthew G. Knepley return 0;
101f918ec44SMatthew G. Knepley }
102f918ec44SMatthew G. Knepley
103f918ec44SMatthew G. Knepley /*
104f918ec44SMatthew G. Knepley // Reference (parent) 2D coordinates: X \in [-1, 1]^2
105f918ec44SMatthew G. Knepley //
106f918ec44SMatthew G. Knepley // Global 3D physical coordinates given by the mesh: xx \in [-R, R]^3
107f918ec44SMatthew G. Knepley // with R radius of the sphere
108f918ec44SMatthew G. Knepley //
109f918ec44SMatthew G. Knepley // Local 3D physical coordinates on the 2D manifold: x \in [-l, l]^3
110f918ec44SMatthew G. Knepley // with l half edge of the cube inscribed in the sphere
111f918ec44SMatthew G. Knepley //
112f918ec44SMatthew G. Knepley // Change of coordinates matrix computed by the library:
113f918ec44SMatthew G. Knepley // (physical 3D coords relative to reference 2D coords)
114f918ec44SMatthew G. Knepley // dxx_j/dX_i (indicial notation) [3 * 2]
115f918ec44SMatthew G. Knepley //
116f918ec44SMatthew G. Knepley // Change of coordinates x (on the 2D manifold) relative to xx (phyisical 3D):
117f918ec44SMatthew G. Knepley // dx_i/dxx_j (indicial notation) [3 * 3]
118f918ec44SMatthew G. Knepley //
119f918ec44SMatthew G. Knepley // Change of coordinates x (on the 2D manifold) relative to X (reference 2D):
120f918ec44SMatthew G. Knepley // (by chain rule)
121f918ec44SMatthew G. Knepley // dx_i/dX_j = dx_i/dxx_k * dxx_k/dX_j [3 * 2]
122f918ec44SMatthew G. Knepley //
123f918ec44SMatthew G. Knepley // modJ is given by the magnitude of the cross product of the columns of dx_i/dX_j
124f918ec44SMatthew G. Knepley //
125f918ec44SMatthew G. Knepley // The quadrature data is stored in the array qdata.
126f918ec44SMatthew G. Knepley //
127f918ec44SMatthew G. Knepley // We require the determinant of the Jacobian to properly compute integrals of
128f918ec44SMatthew G. Knepley // the form: int(u v)
129f918ec44SMatthew G. Knepley //
130f918ec44SMatthew G. Knepley // Qdata: modJ * w
131f918ec44SMatthew G. Knepley */
SetupMassGeoSphere(PetscCtx ctx,const CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)132*2a8381b2SBarry Smith CEED_QFUNCTION(SetupMassGeoSphere)(PetscCtx ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)
133d71ae5a4SJacob Faibussowitsch {
134f918ec44SMatthew G. Knepley const CeedScalar *X = in[0], *J = in[1], *w = in[2];
135f918ec44SMatthew G. Knepley CeedScalar *qdata = out[0];
136f918ec44SMatthew G. Knepley
137d71ae5a4SJacob Faibussowitsch CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i)
138d71ae5a4SJacob Faibussowitsch {
139f918ec44SMatthew G. Knepley const CeedScalar xx[3][1] = {{X[i + 0 * Q]}, {X[i + 1 * Q]}, {X[i + 2 * Q]}};
140f918ec44SMatthew G. Knepley // Read dxxdX Jacobian entries, stored as [[0 3], [1 4], [2 5]]
1419371c9d4SSatish Balay const CeedScalar dxxdX[3][2] = {
1429371c9d4SSatish Balay {J[i + Q * 0], J[i + Q * 3]},
143f918ec44SMatthew G. Knepley {J[i + Q * 1], J[i + Q * 4]},
1449371c9d4SSatish Balay {J[i + Q * 2], J[i + Q * 5]}
1459371c9d4SSatish Balay };
146f918ec44SMatthew G. Knepley // Setup
147f918ec44SMatthew G. Knepley const CeedScalar modxxsq = xx[0][0] * xx[0][0] + xx[1][0] * xx[1][0] + xx[2][0] * xx[2][0];
148f918ec44SMatthew G. Knepley CeedScalar xxsq[3][3];
149f918ec44SMatthew G. Knepley for (int j = 0; j < 3; ++j)
150f918ec44SMatthew G. Knepley for (int k = 0; k < 3; ++k) {
151f918ec44SMatthew G. Knepley xxsq[j][k] = 0.;
1529371c9d4SSatish Balay for (int l = 0; l < 1; ++l) xxsq[j][k] += xx[j][l] * xx[k][l] / (sqrt(modxxsq) * modxxsq);
153f918ec44SMatthew G. Knepley }
154f918ec44SMatthew G. Knepley
1559371c9d4SSatish Balay const CeedScalar dxdxx[3][3] = {
1569371c9d4SSatish Balay {1. / sqrt(modxxsq) - xxsq[0][0], -xxsq[0][1], -xxsq[0][2] },
157f918ec44SMatthew G. Knepley {-xxsq[1][0], 1. / sqrt(modxxsq) - xxsq[1][1], -xxsq[1][2] },
1589371c9d4SSatish Balay {-xxsq[2][0], -xxsq[2][1], 1. / sqrt(modxxsq) - xxsq[2][2]}
1599371c9d4SSatish Balay };
160f918ec44SMatthew G. Knepley
161f918ec44SMatthew G. Knepley CeedScalar dxdX[3][2];
162f918ec44SMatthew G. Knepley for (int j = 0; j < 3; ++j)
163f918ec44SMatthew G. Knepley for (int k = 0; k < 2; ++k) {
164f918ec44SMatthew G. Knepley dxdX[j][k] = 0.;
1659371c9d4SSatish Balay for (int l = 0; l < 3; ++l) dxdX[j][k] += dxdxx[j][l] * dxxdX[l][k];
166f918ec44SMatthew G. Knepley }
167f918ec44SMatthew G. Knepley // J is given by the cross product of the columns of dxdX
1689371c9d4SSatish 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]}};
169f918ec44SMatthew G. Knepley // Use the magnitude of J as our detJ (volume scaling factor)
170f918ec44SMatthew G. Knepley const CeedScalar modJ = sqrt(J[0][0] * J[0][0] + J[1][0] * J[1][0] + J[2][0] * J[2][0]);
171f918ec44SMatthew G. Knepley qdata[i + Q * 0] = modJ * w[i];
172f918ec44SMatthew G. Knepley }
173f918ec44SMatthew G. Knepley return 0;
174f918ec44SMatthew G. Knepley }
175f918ec44SMatthew G. Knepley
ProcessOptions(MPI_Comm comm,AppCtx * ctx)176d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *ctx)
177d71ae5a4SJacob Faibussowitsch {
178f918ec44SMatthew G. Knepley DMPlexShape shape = DM_SHAPE_UNKNOWN;
179f918ec44SMatthew G. Knepley
180f918ec44SMatthew G. Knepley PetscFunctionBeginUser;
181d0609cedSBarry Smith PetscOptionsBegin(comm, "", "libCEED Test Options", "DMPLEX");
182d0609cedSBarry Smith PetscOptionsEnd();
1839566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnum(NULL, NULL, "-dm_plex_shape", DMPlexShapes, (PetscEnum *)&shape, NULL));
18430602db0SMatthew G. Knepley ctx->setupgeo = NULL;
18530602db0SMatthew G. Knepley ctx->setupgeofname = NULL;
18630602db0SMatthew G. Knepley ctx->apply = Mass;
18730602db0SMatthew G. Knepley ctx->applyfname = Mass_loc;
18830602db0SMatthew G. Knepley ctx->areaExact = 0.0;
189f918ec44SMatthew G. Knepley switch (shape) {
190f918ec44SMatthew G. Knepley case DM_SHAPE_BOX_SURFACE:
191f918ec44SMatthew G. Knepley ctx->setupgeo = SetupMassGeoCube;
192f918ec44SMatthew G. Knepley ctx->setupgeofname = SetupMassGeoCube_loc;
193f918ec44SMatthew G. Knepley ctx->areaExact = 6.0;
194f918ec44SMatthew G. Knepley break;
195f918ec44SMatthew G. Knepley case DM_SHAPE_SPHERE:
196f918ec44SMatthew G. Knepley ctx->setupgeo = SetupMassGeoSphere;
197f918ec44SMatthew G. Knepley ctx->setupgeofname = SetupMassGeoSphere_loc;
198f918ec44SMatthew G. Knepley ctx->areaExact = 4.0 * M_PI;
199f918ec44SMatthew G. Knepley break;
200d71ae5a4SJacob Faibussowitsch default:
201d71ae5a4SJacob Faibussowitsch break;
202f918ec44SMatthew G. Knepley }
2033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
204f918ec44SMatthew G. Knepley }
205f918ec44SMatthew G. Knepley
CreateMesh(MPI_Comm comm,AppCtx * ctx,DM * dm)206d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm)
207d71ae5a4SJacob Faibussowitsch {
208f918ec44SMatthew G. Knepley PetscFunctionBegin;
2099566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm));
2109566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX));
2119566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm));
2129566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
21330602db0SMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED
21430602db0SMatthew G. Knepley {
21530602db0SMatthew G. Knepley Ceed ceed;
21630602db0SMatthew G. Knepley const char *usedresource;
21730602db0SMatthew G. Knepley
2189566063dSJacob Faibussowitsch PetscCall(DMGetCeed(*dm, &ceed));
2199566063dSJacob Faibussowitsch PetscCall(CeedGetResource(ceed, &usedresource));
2209566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)*dm), "libCEED Backend: %s\n", usedresource));
22130602db0SMatthew G. Knepley }
22230602db0SMatthew G. Knepley #endif
2233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
224f918ec44SMatthew G. Knepley }
225f918ec44SMatthew G. Knepley
SetupDiscretization(DM dm)226d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm)
227d71ae5a4SJacob Faibussowitsch {
228f918ec44SMatthew G. Knepley DM cdm;
229a2c9b50fSJeremy L Thompson PetscFE fe, cfe;
230a2c9b50fSJeremy L Thompson PetscInt dim, cnc;
231f918ec44SMatthew G. Knepley PetscBool simplex;
232f918ec44SMatthew G. Knepley
233f918ec44SMatthew G. Knepley PetscFunctionBeginUser;
2349566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
2359566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex));
2369566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, PETSC_DETERMINE, &fe));
2379566063dSJacob Faibussowitsch PetscCall(PetscFESetName(fe, "indicator"));
2389566063dSJacob Faibussowitsch PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
2399566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe));
2409566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm));
2419566063dSJacob Faibussowitsch PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
2429566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cnc));
2439566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, cnc, simplex, NULL, PETSC_DETERMINE, &cfe));
2444c712d99Sksagiyam PetscCall(DMSetCoordinateDisc(dm, cfe, PETSC_FALSE, PETSC_TRUE));
2459566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&cfe));
2469566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm));
2479566063dSJacob Faibussowitsch PetscCall(DMPlexSetClosurePermutationTensor(cdm, PETSC_DETERMINE, NULL));
2483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
249f918ec44SMatthew G. Knepley }
250f918ec44SMatthew G. Knepley
LibCeedSetupByDegree(DM dm,AppCtx * ctx,CeedData * data)251d71ae5a4SJacob Faibussowitsch static PetscErrorCode LibCeedSetupByDegree(DM dm, AppCtx *ctx, CeedData *data)
252d71ae5a4SJacob Faibussowitsch {
253f918ec44SMatthew G. Knepley PetscDS ds;
254f918ec44SMatthew G. Knepley PetscFE fe, cfe;
255f918ec44SMatthew G. Knepley Ceed ceed;
256f918ec44SMatthew G. Knepley CeedElemRestriction Erestrictx, Erestrictu, Erestrictq;
257f918ec44SMatthew G. Knepley CeedQFunction qf_setupgeo;
258f918ec44SMatthew G. Knepley CeedOperator op_setupgeo;
259f918ec44SMatthew G. Knepley CeedVector xcoord;
260f918ec44SMatthew G. Knepley CeedBasis basisu, basisx;
261f918ec44SMatthew G. Knepley CeedInt Nqdata = 1;
262f918ec44SMatthew G. Knepley CeedInt nqpts, nqptsx;
263f918ec44SMatthew G. Knepley DM cdm;
264f918ec44SMatthew G. Knepley Vec coords;
265f918ec44SMatthew G. Knepley const PetscScalar *coordArray;
266f918ec44SMatthew G. Knepley PetscInt dim, cdim, cStart, cEnd, Ncell;
267f918ec44SMatthew G. Knepley
268f918ec44SMatthew G. Knepley PetscFunctionBeginUser;
2699566063dSJacob Faibussowitsch PetscCall(DMGetCeed(dm, &ceed));
2709566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
2719566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim));
2729566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
273f918ec44SMatthew G. Knepley Ncell = cEnd - cStart;
274f918ec44SMatthew G. Knepley // CEED bases
2759566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds));
2769566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&fe));
2779566063dSJacob Faibussowitsch PetscCall(PetscFEGetCeedBasis(fe, &basisu));
2789566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm));
2799566063dSJacob Faibussowitsch PetscCall(DMGetDS(cdm, &ds));
2809566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&cfe));
2819566063dSJacob Faibussowitsch PetscCall(PetscFEGetCeedBasis(cfe, &basisx));
282f918ec44SMatthew G. Knepley
2839566063dSJacob Faibussowitsch PetscCall(DMPlexGetCeedRestriction(cdm, NULL, 0, 0, 0, &Erestrictx));
2849566063dSJacob Faibussowitsch PetscCall(DMPlexGetCeedRestriction(dm, NULL, 0, 0, 0, &Erestrictu));
2859566063dSJacob Faibussowitsch PetscCall(CeedBasisGetNumQuadraturePoints(basisu, &nqpts));
2869566063dSJacob Faibussowitsch PetscCall(CeedBasisGetNumQuadraturePoints(basisx, &nqptsx));
28785f1dce8SJed Brown PetscCheck(nqptsx == nqpts, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of qpoints for u %" CeedInt_FMT " != %" CeedInt_FMT " Number of qpoints for x", nqpts, nqptsx);
2889566063dSJacob Faibussowitsch PetscCall(CeedElemRestrictionCreateStrided(ceed, Ncell, nqpts, Nqdata, Nqdata * Ncell * nqpts, CEED_STRIDES_BACKEND, &Erestrictq));
289f918ec44SMatthew G. Knepley
2909566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coords));
2919566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coords, &coordArray));
2929566063dSJacob Faibussowitsch PetscCall(CeedElemRestrictionCreateVector(Erestrictx, &xcoord, NULL));
2939566063dSJacob Faibussowitsch PetscCall(CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_COPY_VALUES, (PetscScalar *)coordArray));
2949566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coords, &coordArray));
295f918ec44SMatthew G. Knepley
296f918ec44SMatthew G. Knepley // Create the vectors that will be needed in setup and apply
2979566063dSJacob Faibussowitsch PetscCall(CeedElemRestrictionCreateVector(Erestrictu, &data->uceed, NULL));
2989566063dSJacob Faibussowitsch PetscCall(CeedElemRestrictionCreateVector(Erestrictu, &data->vceed, NULL));
2999566063dSJacob Faibussowitsch PetscCall(CeedElemRestrictionCreateVector(Erestrictq, &data->qdata, NULL));
300f918ec44SMatthew G. Knepley
301f918ec44SMatthew G. Knepley // Create the Q-function that builds the operator (i.e. computes its quadrature data) and set its context data
3029566063dSJacob Faibussowitsch PetscCall(CeedQFunctionCreateInterior(ceed, 1, ctx->setupgeo, ctx->setupgeofname, &qf_setupgeo));
3039566063dSJacob Faibussowitsch PetscCall(CeedQFunctionAddInput(qf_setupgeo, "x", cdim, CEED_EVAL_INTERP));
3049566063dSJacob Faibussowitsch PetscCall(CeedQFunctionAddInput(qf_setupgeo, "dx", cdim * dim, CEED_EVAL_GRAD));
3059566063dSJacob Faibussowitsch PetscCall(CeedQFunctionAddInput(qf_setupgeo, "weight", 1, CEED_EVAL_WEIGHT));
3069566063dSJacob Faibussowitsch PetscCall(CeedQFunctionAddOutput(qf_setupgeo, "qdata", Nqdata, CEED_EVAL_NONE));
307f918ec44SMatthew G. Knepley
308f918ec44SMatthew G. Knepley // Set up the mass operator
3099566063dSJacob Faibussowitsch PetscCall(CeedQFunctionCreateInterior(ceed, 1, ctx->apply, ctx->applyfname, &data->qf_apply));
3109566063dSJacob Faibussowitsch PetscCall(CeedQFunctionAddInput(data->qf_apply, "u", 1, CEED_EVAL_INTERP));
3119566063dSJacob Faibussowitsch PetscCall(CeedQFunctionAddInput(data->qf_apply, "qdata", Nqdata, CEED_EVAL_NONE));
3129566063dSJacob Faibussowitsch PetscCall(CeedQFunctionAddOutput(data->qf_apply, "v", 1, CEED_EVAL_INTERP));
313f918ec44SMatthew G. Knepley
314f918ec44SMatthew G. Knepley // Create the operator that builds the quadrature data for the operator
3159566063dSJacob Faibussowitsch PetscCall(CeedOperatorCreate(ceed, qf_setupgeo, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setupgeo));
3169566063dSJacob Faibussowitsch PetscCall(CeedOperatorSetField(op_setupgeo, "x", Erestrictx, basisx, CEED_VECTOR_ACTIVE));
3179566063dSJacob Faibussowitsch PetscCall(CeedOperatorSetField(op_setupgeo, "dx", Erestrictx, basisx, CEED_VECTOR_ACTIVE));
3189566063dSJacob Faibussowitsch PetscCall(CeedOperatorSetField(op_setupgeo, "weight", CEED_ELEMRESTRICTION_NONE, basisx, CEED_VECTOR_NONE));
31910ae67c6SJeremy L Thompson PetscCall(CeedOperatorSetField(op_setupgeo, "qdata", Erestrictq, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
320f918ec44SMatthew G. Knepley
321f918ec44SMatthew G. Knepley // Create the mass operator
3229566063dSJacob Faibussowitsch PetscCall(CeedOperatorCreate(ceed, data->qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &data->op_apply));
3239566063dSJacob Faibussowitsch PetscCall(CeedOperatorSetField(data->op_apply, "u", Erestrictu, basisu, CEED_VECTOR_ACTIVE));
32410ae67c6SJeremy L Thompson PetscCall(CeedOperatorSetField(data->op_apply, "qdata", Erestrictq, CEED_BASIS_NONE, data->qdata));
3259566063dSJacob Faibussowitsch PetscCall(CeedOperatorSetField(data->op_apply, "v", Erestrictu, basisu, CEED_VECTOR_ACTIVE));
326f918ec44SMatthew G. Knepley
327f918ec44SMatthew G. Knepley // Setup qdata
3289566063dSJacob Faibussowitsch PetscCall(CeedOperatorApply(op_setupgeo, xcoord, data->qdata, CEED_REQUEST_IMMEDIATE));
329f918ec44SMatthew G. Knepley
3309566063dSJacob Faibussowitsch PetscCall(CeedElemRestrictionDestroy(&Erestrictq));
3319566063dSJacob Faibussowitsch PetscCall(CeedQFunctionDestroy(&qf_setupgeo));
3329566063dSJacob Faibussowitsch PetscCall(CeedOperatorDestroy(&op_setupgeo));
3339566063dSJacob Faibussowitsch PetscCall(CeedVectorDestroy(&xcoord));
3343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
335f918ec44SMatthew G. Knepley }
336f918ec44SMatthew G. Knepley
main(int argc,char ** argv)337d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
338d71ae5a4SJacob Faibussowitsch {
339f918ec44SMatthew G. Knepley MPI_Comm comm;
340f918ec44SMatthew G. Knepley DM dm;
341f918ec44SMatthew G. Knepley AppCtx ctx;
342f918ec44SMatthew G. Knepley Vec U, Uloc, V, Vloc;
343f918ec44SMatthew G. Knepley PetscScalar *v;
344f918ec44SMatthew G. Knepley PetscScalar area;
345f918ec44SMatthew G. Knepley CeedData ceeddata;
346f918ec44SMatthew G. Knepley
3479566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
348f918ec44SMatthew G. Knepley comm = PETSC_COMM_WORLD;
3499566063dSJacob Faibussowitsch PetscCall(ProcessOptions(comm, &ctx));
3509566063dSJacob Faibussowitsch PetscCall(CreateMesh(comm, &ctx, &dm));
3519566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm));
352f918ec44SMatthew G. Knepley
3539566063dSJacob Faibussowitsch PetscCall(LibCeedSetupByDegree(dm, &ctx, &ceeddata));
354f918ec44SMatthew G. Knepley
3559566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &U));
3569566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dm, &Uloc));
3579566063dSJacob Faibussowitsch PetscCall(VecDuplicate(U, &V));
3589566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Uloc, &Vloc));
359f918ec44SMatthew G. Knepley
36030602db0SMatthew G. Knepley /**/
3619566063dSJacob Faibussowitsch PetscCall(VecSet(Uloc, 1.));
3629566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(V));
3639566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Vloc));
3649566063dSJacob Faibussowitsch PetscCall(VecGetArray(Vloc, &v));
3659566063dSJacob Faibussowitsch PetscCall(CeedVectorSetArray(ceeddata.vceed, CEED_MEM_HOST, CEED_USE_POINTER, v));
3669566063dSJacob Faibussowitsch PetscCall(CeedVectorSetValue(ceeddata.uceed, 1.0));
3679566063dSJacob Faibussowitsch PetscCall(CeedOperatorApply(ceeddata.op_apply, ceeddata.uceed, ceeddata.vceed, CEED_REQUEST_IMMEDIATE));
3689566063dSJacob Faibussowitsch PetscCall(CeedVectorTakeArray(ceeddata.vceed, CEED_MEM_HOST, NULL));
3699566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Vloc, &v));
3709566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm, Vloc, ADD_VALUES, V));
3719566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm, Vloc, ADD_VALUES, V));
372f918ec44SMatthew G. Knepley
3739566063dSJacob Faibussowitsch PetscCall(VecSum(V, &area));
374f918ec44SMatthew G. Knepley if (ctx.areaExact > 0.) {
375f918ec44SMatthew G. Knepley PetscReal error = PetscAbsReal(area - ctx.areaExact);
376f918ec44SMatthew G. Knepley PetscReal tol = PETSC_SMALL;
377f918ec44SMatthew G. Knepley
37863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Exact mesh surface area : % .*f\n", PetscAbsReal(ctx.areaExact - round(ctx.areaExact)) > 1E-15 ? 14 : 1, (double)ctx.areaExact));
37963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Computed mesh surface area : % .*f\n", PetscAbsScalar(area - round(area)) > 1E-15 ? 14 : 1, (double)PetscRealPart(area)));
380f918ec44SMatthew G. Knepley if (error > tol) {
38163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Area error : % .14g\n", (double)error));
382f918ec44SMatthew G. Knepley } else {
38363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Area verifies!\n"));
384f918ec44SMatthew G. Knepley }
385f918ec44SMatthew G. Knepley }
386f918ec44SMatthew G. Knepley
3879566063dSJacob Faibussowitsch PetscCall(CeedDataDestroy(&ceeddata));
3889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U));
3899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Uloc));
3909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&V));
3919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Vloc));
3929566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
393f918ec44SMatthew G. Knepley return PetscFinalize();
394f918ec44SMatthew G. Knepley }
395f918ec44SMatthew G. Knepley
396f918ec44SMatthew G. Knepley /*TEST
397f918ec44SMatthew G. Knepley
398f918ec44SMatthew G. Knepley build:
399f918ec44SMatthew G. Knepley requires: libceed
400f918ec44SMatthew G. Knepley
401f918ec44SMatthew G. Knepley testset:
402e600fa54SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 3 -dm_view -dm_petscds_view \
403dc431b0cSMatthew G. Knepley -petscfe_default_quadrature_order 4 -cdm_default_quadrature_order 4
404d778fca9SStefano Zampini filter: sed -e "s /cpu/self/xsmm /cpu/self/opt " -e "s /cpu/self/avx /cpu/self/opt "
405f918ec44SMatthew G. Knepley
406f918ec44SMatthew G. Knepley test:
407f918ec44SMatthew G. Knepley suffix: cube_3
408f918ec44SMatthew G. Knepley args: -dm_plex_shape box_surface -dm_refine 2
409f918ec44SMatthew G. Knepley
410f918ec44SMatthew G. Knepley test:
411f918ec44SMatthew G. Knepley suffix: cube_3_p4
412f918ec44SMatthew G. Knepley nsize: 4
413fb796b39SJed Brown args: -petscpartitioner_type simple -dm_refine_pre 1 -dm_plex_shape box_surface -dm_refine 1
414f918ec44SMatthew G. Knepley
415f918ec44SMatthew G. Knepley test:
416f918ec44SMatthew G. Knepley suffix: sphere_3
417f918ec44SMatthew G. Knepley args: -dm_plex_shape sphere -dm_refine 3
418f918ec44SMatthew G. Knepley
419f918ec44SMatthew G. Knepley test:
420f918ec44SMatthew G. Knepley suffix: sphere_3_p4
421f918ec44SMatthew G. Knepley nsize: 4
422fb796b39SJed Brown args: -petscpartitioner_type simple -dm_refine_pre 1 -dm_plex_shape sphere -dm_refine 2
423f918ec44SMatthew G. Knepley
424f918ec44SMatthew G. Knepley TEST*/
425