xref: /libCEED/examples/petsc/src/petscutils.c (revision b6972d7456611f84b0e462eb1490bcb662442e6a)
1e83e87a5Sjeremylt #include "../include/petscutils.h"
2e83e87a5Sjeremylt 
3e83e87a5Sjeremylt // -----------------------------------------------------------------------------
4e83e87a5Sjeremylt // Convert PETSc MemType to libCEED MemType
5e83e87a5Sjeremylt // -----------------------------------------------------------------------------
62b730f8bSJeremy L Thompson CeedMemType MemTypeP2C(PetscMemType mem_type) { return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST; }
7e83e87a5Sjeremylt 
8*b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
9*b6972d74SZach Atkins // PETSc-libCEED memory space utilities
10*b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
11*b6972d74SZach Atkins PetscErrorCode VecP2C(Vec X_petsc, PetscMemType *mem_type, CeedVector x_ceed) {
12*b6972d74SZach Atkins   PetscScalar *x;
13*b6972d74SZach Atkins 
14*b6972d74SZach Atkins   PetscFunctionBeginUser;
15*b6972d74SZach Atkins   PetscCall(VecGetArrayAndMemType(X_petsc, &x, mem_type));
16*b6972d74SZach Atkins   CeedVectorSetArray(x_ceed, MemTypeP2C(*mem_type), CEED_USE_POINTER, x);
17*b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
18*b6972d74SZach Atkins }
19*b6972d74SZach Atkins 
20*b6972d74SZach Atkins PetscErrorCode VecC2P(CeedVector x_ceed, PetscMemType mem_type, Vec X_petsc) {
21*b6972d74SZach Atkins   PetscScalar *x;
22*b6972d74SZach Atkins 
23*b6972d74SZach Atkins   PetscFunctionBeginUser;
24*b6972d74SZach Atkins   CeedVectorTakeArray(x_ceed, MemTypeP2C(mem_type), &x);
25*b6972d74SZach Atkins   PetscCall(VecRestoreArrayAndMemType(X_petsc, &x));
26*b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
27*b6972d74SZach Atkins }
28*b6972d74SZach Atkins 
29*b6972d74SZach Atkins PetscErrorCode VecReadP2C(Vec X_petsc, PetscMemType *mem_type, CeedVector x_ceed) {
30*b6972d74SZach Atkins   PetscScalar *x;
31*b6972d74SZach Atkins 
32*b6972d74SZach Atkins   PetscFunctionBeginUser;
33*b6972d74SZach Atkins   PetscCall(VecGetArrayReadAndMemType(X_petsc, (const PetscScalar **)&x, mem_type));
34*b6972d74SZach Atkins   CeedVectorSetArray(x_ceed, MemTypeP2C(*mem_type), CEED_USE_POINTER, x);
35*b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
36*b6972d74SZach Atkins }
37*b6972d74SZach Atkins 
38*b6972d74SZach Atkins PetscErrorCode VecReadC2P(CeedVector x_ceed, PetscMemType mem_type, Vec X_petsc) {
39*b6972d74SZach Atkins   PetscScalar *x;
40*b6972d74SZach Atkins 
41*b6972d74SZach Atkins   PetscFunctionBeginUser;
42*b6972d74SZach Atkins   CeedVectorTakeArray(x_ceed, MemTypeP2C(mem_type), &x);
43*b6972d74SZach Atkins   PetscCall(VecRestoreArrayReadAndMemType(X_petsc, (const PetscScalar **)&x));
44*b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
45*b6972d74SZach Atkins }
46*b6972d74SZach Atkins 
47e83e87a5Sjeremylt // -----------------------------------------------------------------------------
482c58efb6Sjeremylt // Apply 3D Kershaw mesh transformation
492c58efb6Sjeremylt // -----------------------------------------------------------------------------
502c58efb6Sjeremylt // Transition from a value of "a" for x=0, to a value of "b" for x=1.  Optionally
512c58efb6Sjeremylt // smooth -- see the commented versions at the end.
522c58efb6Sjeremylt static double step(const double a, const double b, double x) {
532c58efb6Sjeremylt   if (x <= 0) return a;
542c58efb6Sjeremylt   if (x >= 1) return b;
552c58efb6Sjeremylt   return a + (b - a) * (x);
562c58efb6Sjeremylt }
572c58efb6Sjeremylt 
582c58efb6Sjeremylt // 1D transformation at the right boundary
592b730f8bSJeremy L Thompson static double right(const double eps, const double x) { return (x <= 0.5) ? (2 - eps) * x : 1 + eps * (x - 1); }
602c58efb6Sjeremylt 
612c58efb6Sjeremylt // 1D transformation at the left boundary
622b730f8bSJeremy L Thompson static double left(const double eps, const double x) { return 1 - right(eps, 1 - x); }
632c58efb6Sjeremylt 
642c58efb6Sjeremylt // Apply 3D Kershaw mesh transformation
652c58efb6Sjeremylt // The eps parameters are in (0, 1]
662c58efb6Sjeremylt // Uniform mesh is recovered for eps=1
679b072555Sjeremylt PetscErrorCode Kershaw(DM dm_orig, PetscScalar eps) {
682c58efb6Sjeremylt   Vec          coord;
692c58efb6Sjeremylt   PetscInt     ncoord;
702c58efb6Sjeremylt   PetscScalar *c;
712c58efb6Sjeremylt 
722c58efb6Sjeremylt   PetscFunctionBeginUser;
732b730f8bSJeremy L Thompson 
742b730f8bSJeremy L Thompson   PetscCall(DMGetCoordinatesLocal(dm_orig, &coord));
752b730f8bSJeremy L Thompson   PetscCall(VecGetLocalSize(coord, &ncoord));
762b730f8bSJeremy L Thompson   PetscCall(VecGetArray(coord, &c));
772c58efb6Sjeremylt 
782c58efb6Sjeremylt   for (PetscInt i = 0; i < ncoord; i += 3) {
792c58efb6Sjeremylt     PetscScalar x = c[i], y = c[i + 1], z = c[i + 2];
802c58efb6Sjeremylt     PetscInt    layer  = x * 6;
812c58efb6Sjeremylt     PetscScalar lambda = (x - layer / 6.0) * 6;
822c58efb6Sjeremylt     c[i]               = x;
832c58efb6Sjeremylt 
842c58efb6Sjeremylt     switch (layer) {
852c58efb6Sjeremylt       case 0:
862c58efb6Sjeremylt         c[i + 1] = left(eps, y);
872c58efb6Sjeremylt         c[i + 2] = left(eps, z);
882c58efb6Sjeremylt         break;
892c58efb6Sjeremylt       case 1:
902c58efb6Sjeremylt       case 4:
912c58efb6Sjeremylt         c[i + 1] = step(left(eps, y), right(eps, y), lambda);
922c58efb6Sjeremylt         c[i + 2] = step(left(eps, z), right(eps, z), lambda);
932c58efb6Sjeremylt         break;
942c58efb6Sjeremylt       case 2:
952c58efb6Sjeremylt         c[i + 1] = step(right(eps, y), left(eps, y), lambda / 2);
962c58efb6Sjeremylt         c[i + 2] = step(right(eps, z), left(eps, z), lambda / 2);
972c58efb6Sjeremylt         break;
982c58efb6Sjeremylt       case 3:
992c58efb6Sjeremylt         c[i + 1] = step(right(eps, y), left(eps, y), (1 + lambda) / 2);
1002c58efb6Sjeremylt         c[i + 2] = step(right(eps, z), left(eps, z), (1 + lambda) / 2);
1012c58efb6Sjeremylt         break;
1022c58efb6Sjeremylt       default:
1032c58efb6Sjeremylt         c[i + 1] = right(eps, y);
1042c58efb6Sjeremylt         c[i + 2] = right(eps, z);
1052c58efb6Sjeremylt     }
1062c58efb6Sjeremylt   }
1072b730f8bSJeremy L Thompson   PetscCall(VecRestoreArray(coord, &c));
108ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1092c58efb6Sjeremylt }
1102c58efb6Sjeremylt 
1112c58efb6Sjeremylt // -----------------------------------------------------------------------------
112e83e87a5Sjeremylt // Create BC label
113e83e87a5Sjeremylt // -----------------------------------------------------------------------------
114e83e87a5Sjeremylt static PetscErrorCode CreateBCLabel(DM dm, const char name[]) {
115e83e87a5Sjeremylt   DMLabel label;
116e83e87a5Sjeremylt 
117e83e87a5Sjeremylt   PetscFunctionBeginUser;
118e83e87a5Sjeremylt 
119751eb813Srezgarshakeri   PetscCall(DMCreateLabel(dm, name));
120751eb813Srezgarshakeri   PetscCall(DMGetLabel(dm, name, &label));
121751eb813Srezgarshakeri   PetscCall(DMPlexMarkBoundaryFaces(dm, PETSC_DETERMINE, label));
122751eb813Srezgarshakeri   PetscCall(DMPlexLabelComplete(dm, label));
123e83e87a5Sjeremylt 
124ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
125fd8f24faSSebastian Grimberg }
126e83e87a5Sjeremylt 
127e83e87a5Sjeremylt // -----------------------------------------------------------------------------
128e83e87a5Sjeremylt // This function sets up a DM for a given degree
129e83e87a5Sjeremylt // -----------------------------------------------------------------------------
1302b730f8bSJeremy L Thompson PetscErrorCode SetupDMByDegree(DM dm, PetscInt p_degree, PetscInt q_extra, PetscInt num_comp_u, PetscInt dim, bool enforce_bc) {
1312b730f8bSJeremy L Thompson   PetscInt  marker_ids[1] = {1};
132de1229c5Srezgarshakeri   PetscInt  q_degree      = p_degree + q_extra;
133e83e87a5Sjeremylt   PetscFE   fe;
13465233696SJeremy L Thompson   MPI_Comm  comm;
135129d8736Srezgarshakeri   PetscBool is_simplex = PETSC_TRUE;
136e83e87a5Sjeremylt 
137e83e87a5Sjeremylt   PetscFunctionBeginUser;
138e83e87a5Sjeremylt 
139129d8736Srezgarshakeri   // Check if simplex or tensor-product mesh
1402b730f8bSJeremy L Thompson   PetscCall(DMPlexIsSimplex(dm, &is_simplex));
141e83e87a5Sjeremylt   // Setup FE
1422b730f8bSJeremy L Thompson   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1432b730f8bSJeremy L Thompson   PetscCall(PetscFECreateLagrange(comm, dim, num_comp_u, is_simplex, p_degree, q_degree, &fe));
1442b730f8bSJeremy L Thompson   PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
1452b730f8bSJeremy L Thompson   PetscCall(DMCreateDS(dm));
146129d8736Srezgarshakeri 
147129d8736Srezgarshakeri   {
148129d8736Srezgarshakeri     // create FE field for coordinates
1497ed3e4cdSJeremy L Thompson     PetscFE  fe_coords;
1507ed3e4cdSJeremy L Thompson     PetscInt num_comp_coord;
1512b730f8bSJeremy L Thompson     PetscCall(DMGetCoordinateDim(dm, &num_comp_coord));
1522b730f8bSJeremy L Thompson     PetscCall(PetscFECreateLagrange(comm, dim, num_comp_coord, is_simplex, 1, q_degree, &fe_coords));
15349a40c8aSKenneth E. Jansen     PetscCall(DMSetCoordinateDisc(dm, fe_coords, PETSC_TRUE));
1542b730f8bSJeremy L Thompson     PetscCall(PetscFEDestroy(&fe_coords));
1557ed3e4cdSJeremy L Thompson   }
156de1229c5Srezgarshakeri 
157ab1ff315Srezgarshakeri   // Setup Dirichlet BC
158ab1ff315Srezgarshakeri   // Note bp1, bp2 are projection and we don't need to apply BC
159ab1ff315Srezgarshakeri   // For bp3,bp4, the target function is zero on the boundaries
160ab1ff315Srezgarshakeri   // So we pass bcFunc = NULL in DMAddBoundary function
1619b072555Sjeremylt   if (enforce_bc) {
1629b072555Sjeremylt     PetscBool has_label;
163*b6972d74SZach Atkins     PetscCall(DMHasLabel(dm, "marker", &has_label));
1642b730f8bSJeremy L Thompson     if (!has_label) {
165*b6972d74SZach Atkins       PetscCall(CreateBCLabel(dm, "marker"));
1662b730f8bSJeremy L Thompson     }
16769f5adf1Sjeremylt     DMLabel label;
1682b730f8bSJeremy L Thompson     PetscCall(DMGetLabel(dm, "marker", &label));
1692b730f8bSJeremy L Thompson     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, marker_ids, 0, 0, NULL, NULL, NULL, NULL, NULL));
170751eb813Srezgarshakeri     PetscCall(DMSetOptionsPrefix(dm, "final_"));
171751eb813Srezgarshakeri     PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
172e83e87a5Sjeremylt   }
173129d8736Srezgarshakeri 
174129d8736Srezgarshakeri   if (!is_simplex) {
175129d8736Srezgarshakeri     DM dm_coord;
1762b730f8bSJeremy L Thompson     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
1772b730f8bSJeremy L Thompson     PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
1782b730f8bSJeremy L Thompson     PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL));
179129d8736Srezgarshakeri   }
1802b730f8bSJeremy L Thompson   PetscCall(PetscFEDestroy(&fe));
181e83e87a5Sjeremylt 
182ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
183fd8f24faSSebastian Grimberg }
184e83e87a5Sjeremylt 
185e83e87a5Sjeremylt // -----------------------------------------------------------------------------
186e83e87a5Sjeremylt // Get CEED restriction data from DMPlex
187e83e87a5Sjeremylt // -----------------------------------------------------------------------------
1882b730f8bSJeremy L Thompson PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr) {
1897ed3e4cdSJeremy L Thompson   PetscInt num_elem, elem_size, num_dof, num_comp, *elem_restr_offsets;
190e83e87a5Sjeremylt 
191e83e87a5Sjeremylt   PetscFunctionBeginUser;
192e83e87a5Sjeremylt 
1932b730f8bSJeremy L Thompson   PetscCall(DMPlexGetLocalOffsets(dm, domain_label, value, height, 0, &num_elem, &elem_size, &num_comp, &num_dof, &elem_restr_offsets));
194e83e87a5Sjeremylt 
1952b730f8bSJeremy L Thompson   CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp, 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES, elem_restr_offsets, elem_restr);
1962b730f8bSJeremy L Thompson   PetscCall(PetscFree(elem_restr_offsets));
1977ed3e4cdSJeremy L Thompson 
198ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
199fd8f24faSSebastian Grimberg }
200e83e87a5Sjeremylt 
201e83e87a5Sjeremylt // -----------------------------------------------------------------------------
202129d8736Srezgarshakeri // Utility function - convert from DMPolytopeType to CeedElemTopology
203129d8736Srezgarshakeri // -----------------------------------------------------------------------------
204de1229c5Srezgarshakeri CeedElemTopology ElemTopologyP2C(DMPolytopeType cell_type) {
205129d8736Srezgarshakeri   switch (cell_type) {
2062b730f8bSJeremy L Thompson     case DM_POLYTOPE_TRIANGLE:
2072b730f8bSJeremy L Thompson       return CEED_TOPOLOGY_TRIANGLE;
2082b730f8bSJeremy L Thompson     case DM_POLYTOPE_QUADRILATERAL:
2092b730f8bSJeremy L Thompson       return CEED_TOPOLOGY_QUAD;
2102b730f8bSJeremy L Thompson     case DM_POLYTOPE_TETRAHEDRON:
2112b730f8bSJeremy L Thompson       return CEED_TOPOLOGY_TET;
2122b730f8bSJeremy L Thompson     case DM_POLYTOPE_HEXAHEDRON:
2132b730f8bSJeremy L Thompson       return CEED_TOPOLOGY_HEX;
2142b730f8bSJeremy L Thompson     default:
2152b730f8bSJeremy L Thompson       return 0;
216129d8736Srezgarshakeri   }
217129d8736Srezgarshakeri }
218129d8736Srezgarshakeri 
219129d8736Srezgarshakeri // -----------------------------------------------------------------------------
220f755c37aSrezgarshakeri // Convert DM field to DS field
221129d8736Srezgarshakeri // -----------------------------------------------------------------------------
2222b730f8bSJeremy L Thompson PetscErrorCode DMFieldToDSField(DM dm, DMLabel domain_label, PetscInt dm_field, PetscInt *ds_field) {
223129d8736Srezgarshakeri   PetscDS         ds;
224129d8736Srezgarshakeri   IS              field_is;
225129d8736Srezgarshakeri   const PetscInt *fields;
226129d8736Srezgarshakeri   PetscInt        num_fields;
227129d8736Srezgarshakeri 
228f755c37aSrezgarshakeri   PetscFunctionBeginUser;
229f755c37aSrezgarshakeri 
230129d8736Srezgarshakeri   // Translate dm_field to ds_field
231baf96a30SJames Wright   PetscCall(DMGetRegionDS(dm, domain_label, &field_is, &ds, NULL));
232f755c37aSrezgarshakeri   PetscCall(ISGetIndices(field_is, &fields));
233f755c37aSrezgarshakeri   PetscCall(ISGetSize(field_is, &num_fields));
234129d8736Srezgarshakeri   for (PetscInt i = 0; i < num_fields; i++) {
235129d8736Srezgarshakeri     if (dm_field == fields[i]) {
236f755c37aSrezgarshakeri       *ds_field = i;
237129d8736Srezgarshakeri       break;
238129d8736Srezgarshakeri     }
239129d8736Srezgarshakeri   }
240f755c37aSrezgarshakeri   PetscCall(ISRestoreIndices(field_is, &fields));
241f755c37aSrezgarshakeri 
2422b730f8bSJeremy L Thompson   if (*ds_field == -1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Could not find dm_field %" PetscInt_FMT " in DS", dm_field);
243f755c37aSrezgarshakeri 
244ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
245129d8736Srezgarshakeri }
246129d8736Srezgarshakeri 
247f755c37aSrezgarshakeri // -----------------------------------------------------------------------------
248f755c37aSrezgarshakeri // Create libCEED Basis from PetscTabulation
249f755c37aSrezgarshakeri // -----------------------------------------------------------------------------
2502b730f8bSJeremy L Thompson PetscErrorCode BasisCreateFromTabulation(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt face, PetscFE fe,
2512b730f8bSJeremy L Thompson                                          PetscTabulation basis_tabulation, PetscQuadrature quadrature, CeedBasis *basis) {
252f755c37aSrezgarshakeri   PetscInt           first_point;
253129d8736Srezgarshakeri   PetscInt           ids[1] = {label_value};
254129d8736Srezgarshakeri   DMLabel            depth_label;
255129d8736Srezgarshakeri   DMPolytopeType     cell_type;
256129d8736Srezgarshakeri   CeedElemTopology   elem_topo;
257f755c37aSrezgarshakeri   PetscScalar       *q_points, *interp, *grad;
258f755c37aSrezgarshakeri   const PetscScalar *q_weights;
259f755c37aSrezgarshakeri   PetscDualSpace     dual_space;
260f755c37aSrezgarshakeri   PetscInt           num_dual_basis_vectors;
261f755c37aSrezgarshakeri   PetscInt           dim, num_comp, P, Q;
262f755c37aSrezgarshakeri 
263f755c37aSrezgarshakeri   PetscFunctionBeginUser;
264f755c37aSrezgarshakeri 
265f755c37aSrezgarshakeri   // General basis information
266f755c37aSrezgarshakeri   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
267f755c37aSrezgarshakeri   PetscCall(PetscFEGetNumComponents(fe, &num_comp));
268f755c37aSrezgarshakeri   PetscCall(PetscFEGetDualSpace(fe, &dual_space));
269f755c37aSrezgarshakeri   PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
270f755c37aSrezgarshakeri   P = num_dual_basis_vectors / num_comp;
271129d8736Srezgarshakeri 
272129d8736Srezgarshakeri   // Use depth label if no domain label present
273129d8736Srezgarshakeri   if (!domain_label) {
274129d8736Srezgarshakeri     PetscInt depth;
275129d8736Srezgarshakeri 
276f755c37aSrezgarshakeri     PetscCall(DMPlexGetDepth(dm, &depth));
277f755c37aSrezgarshakeri     PetscCall(DMPlexGetDepthLabel(dm, &depth_label));
278129d8736Srezgarshakeri     ids[0] = depth - height;
279129d8736Srezgarshakeri   }
280f755c37aSrezgarshakeri 
281129d8736Srezgarshakeri   // Get cell interp, grad, and quadrature data
2822b730f8bSJeremy L Thompson   PetscCall(DMGetFirstLabeledPoint(dm, dm, domain_label ? domain_label : depth_label, 1, ids, height, &first_point, NULL));
283f755c37aSrezgarshakeri   PetscCall(DMPlexGetCellType(dm, first_point, &cell_type));
284129d8736Srezgarshakeri   elem_topo = ElemTopologyP2C(cell_type);
2852b730f8bSJeremy L Thompson   if (!elem_topo) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMPlex topology not supported");
286f755c37aSrezgarshakeri   {
287f755c37aSrezgarshakeri     size_t             q_points_size;
288f755c37aSrezgarshakeri     const PetscScalar *q_points_petsc;
289f755c37aSrezgarshakeri     PetscInt           q_dim;
290f755c37aSrezgarshakeri 
2912b730f8bSJeremy L Thompson     PetscCall(PetscQuadratureGetData(quadrature, &q_dim, NULL, &Q, &q_points_petsc, &q_weights));
292f755c37aSrezgarshakeri     q_points_size = Q * dim * sizeof(CeedScalar);
293f755c37aSrezgarshakeri     PetscCall(PetscCalloc(q_points_size, &q_points));
294f755c37aSrezgarshakeri     for (PetscInt q = 0; q < Q; q++) {
2952b730f8bSJeremy L Thompson       for (PetscInt d = 0; d < q_dim; d++) q_points[q * dim + d] = q_points_petsc[q * q_dim + d];
296f755c37aSrezgarshakeri     }
297f755c37aSrezgarshakeri   }
298f755c37aSrezgarshakeri 
299129d8736Srezgarshakeri   // Convert to libCEED orientation
300f755c37aSrezgarshakeri   {
301f755c37aSrezgarshakeri     PetscBool       is_simplex  = PETSC_FALSE;
302f755c37aSrezgarshakeri     IS              permutation = NULL;
303f755c37aSrezgarshakeri     const PetscInt *permutation_indices;
304f755c37aSrezgarshakeri 
305f755c37aSrezgarshakeri     PetscCall(DMPlexIsSimplex(dm, &is_simplex));
306f755c37aSrezgarshakeri     if (!is_simplex) {
307f755c37aSrezgarshakeri       PetscSection section;
308f755c37aSrezgarshakeri 
309f755c37aSrezgarshakeri       // -- Get permutation
310f755c37aSrezgarshakeri       PetscCall(DMGetLocalSection(dm, &section));
3112b730f8bSJeremy L Thompson       PetscCall(PetscSectionGetClosurePermutation(section, (PetscObject)dm, dim, num_comp * P, &permutation));
312f755c37aSrezgarshakeri       PetscCall(ISGetIndices(permutation, &permutation_indices));
313f755c37aSrezgarshakeri     }
314f755c37aSrezgarshakeri 
315f755c37aSrezgarshakeri     // -- Copy interp, grad matrices
316f755c37aSrezgarshakeri     PetscCall(PetscCalloc(P * Q * sizeof(CeedScalar), &interp));
317f755c37aSrezgarshakeri     PetscCall(PetscCalloc(P * Q * dim * sizeof(CeedScalar), &grad));
318129d8736Srezgarshakeri     const CeedInt c = 0;
319129d8736Srezgarshakeri     for (CeedInt q = 0; q < Q; q++) {
320f755c37aSrezgarshakeri       for (CeedInt p_ceed = 0; p_ceed < P; p_ceed++) {
3212b730f8bSJeremy L Thompson         CeedInt p_petsc = is_simplex ? (p_ceed * num_comp) : permutation_indices[p_ceed * num_comp];
322f755c37aSrezgarshakeri 
3232b730f8bSJeremy L Thompson         interp[q * P + p_ceed] = basis_tabulation->T[0][((face * Q + q) * P * num_comp + p_petsc) * num_comp + c];
324129d8736Srezgarshakeri         for (CeedInt d = 0; d < dim; d++) {
3252b730f8bSJeremy L Thompson           grad[(d * Q + q) * P + p_ceed] = basis_tabulation->T[1][(((face * Q + q) * P * num_comp + p_petsc) * num_comp + c) * dim + d];
326129d8736Srezgarshakeri         }
327129d8736Srezgarshakeri       }
328129d8736Srezgarshakeri     }
329f755c37aSrezgarshakeri 
330f755c37aSrezgarshakeri     // -- Cleanup
331f755c37aSrezgarshakeri     if (permutation) PetscCall(ISRestoreIndices(permutation, &permutation_indices));
332f755c37aSrezgarshakeri     PetscCall(ISDestroy(&permutation));
333f755c37aSrezgarshakeri   }
334f755c37aSrezgarshakeri 
335f755c37aSrezgarshakeri   // Finally, create libCEED basis
3362b730f8bSJeremy L Thompson   CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad, q_points, q_weights, basis);
337f755c37aSrezgarshakeri   PetscCall(PetscFree(q_points));
338f755c37aSrezgarshakeri   PetscCall(PetscFree(interp));
339f755c37aSrezgarshakeri   PetscCall(PetscFree(grad));
340f755c37aSrezgarshakeri 
341ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
342f755c37aSrezgarshakeri }
343f755c37aSrezgarshakeri 
344f755c37aSrezgarshakeri // -----------------------------------------------------------------------------
345f755c37aSrezgarshakeri // Get CEED Basis from DMPlex
346f755c37aSrezgarshakeri // -----------------------------------------------------------------------------
3472b730f8bSJeremy L Thompson PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, BPData bp_data,
3482b730f8bSJeremy L Thompson                                    CeedBasis *basis) {
349f755c37aSrezgarshakeri   PetscDS         ds;
350f755c37aSrezgarshakeri   PetscFE         fe;
351f755c37aSrezgarshakeri   PetscQuadrature quadrature;
352f755c37aSrezgarshakeri   PetscBool       is_simplex = PETSC_TRUE;
353f755c37aSrezgarshakeri   PetscInt        ds_field   = -1;
354f755c37aSrezgarshakeri 
355f755c37aSrezgarshakeri   PetscFunctionBeginUser;
356f755c37aSrezgarshakeri 
357f755c37aSrezgarshakeri   // Get element information
358baf96a30SJames Wright   PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL));
359f755c37aSrezgarshakeri   PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field));
360f755c37aSrezgarshakeri   PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe));
361f755c37aSrezgarshakeri   PetscCall(PetscFEGetHeightSubspace(fe, height, &fe));
362f755c37aSrezgarshakeri   PetscCall(PetscFEGetQuadrature(fe, &quadrature));
363f755c37aSrezgarshakeri 
364f755c37aSrezgarshakeri   // Check if simplex or tensor-product mesh
365f755c37aSrezgarshakeri   PetscCall(DMPlexIsSimplex(dm, &is_simplex));
366f755c37aSrezgarshakeri 
367f755c37aSrezgarshakeri   // Build libCEED basis
368f755c37aSrezgarshakeri   if (is_simplex) {
369f755c37aSrezgarshakeri     PetscTabulation basis_tabulation;
370f755c37aSrezgarshakeri     PetscInt        num_derivatives = 1, face = 0;
371f755c37aSrezgarshakeri 
372f755c37aSrezgarshakeri     PetscCall(PetscFEGetCellTabulation(fe, num_derivatives, &basis_tabulation));
3732b730f8bSJeremy L Thompson     PetscCall(BasisCreateFromTabulation(ceed, dm, domain_label, label_value, height, face, fe, basis_tabulation, quadrature, basis));
374129d8736Srezgarshakeri   } else {
375f755c37aSrezgarshakeri     PetscDualSpace dual_space;
376f755c37aSrezgarshakeri     PetscInt       num_dual_basis_vectors;
377f755c37aSrezgarshakeri     PetscInt       dim, num_comp, P, Q;
378f755c37aSrezgarshakeri 
379f755c37aSrezgarshakeri     PetscCall(PetscFEGetSpatialDimension(fe, &dim));
380f755c37aSrezgarshakeri     PetscCall(PetscFEGetNumComponents(fe, &num_comp));
381f755c37aSrezgarshakeri     PetscCall(PetscFEGetDualSpace(fe, &dual_space));
382f755c37aSrezgarshakeri     PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
383f755c37aSrezgarshakeri     P = num_dual_basis_vectors / num_comp;
384f755c37aSrezgarshakeri     PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &Q, NULL, NULL));
385f755c37aSrezgarshakeri 
386129d8736Srezgarshakeri     CeedInt P_1d = (CeedInt)round(pow(P, 1.0 / dim));
387129d8736Srezgarshakeri     CeedInt Q_1d = (CeedInt)round(pow(Q, 1.0 / dim));
388129d8736Srezgarshakeri 
3892b730f8bSJeremy L Thompson     CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp, P_1d, Q_1d, bp_data.q_mode, basis);
390129d8736Srezgarshakeri   }
391129d8736Srezgarshakeri 
392ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
393f755c37aSrezgarshakeri }
394129d8736Srezgarshakeri 
395129d8736Srezgarshakeri // -----------------------------------------------------------------------------
396de1229c5Srezgarshakeri // Utilities
397de1229c5Srezgarshakeri // -----------------------------------------------------------------------------
398de1229c5Srezgarshakeri 
399de1229c5Srezgarshakeri // Utility function, compute three factors of an integer
400de1229c5Srezgarshakeri static void Split3(PetscInt size, PetscInt m[3], bool reverse) {
401de1229c5Srezgarshakeri   for (PetscInt d = 0, size_left = size; d < 3; d++) {
402de1229c5Srezgarshakeri     PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(size_left, 1. / (3 - d)));
403de1229c5Srezgarshakeri     while (try * (size_left / try) != size_left) try++;
404de1229c5Srezgarshakeri     m[reverse ? 2 - d : d] = try;
405de1229c5Srezgarshakeri     size_left /= try;
406de1229c5Srezgarshakeri   }
407de1229c5Srezgarshakeri }
408de1229c5Srezgarshakeri 
4092b730f8bSJeremy L Thompson static int Max3(const PetscInt a[3]) { return PetscMax(a[0], PetscMax(a[1], a[2])); }
410de1229c5Srezgarshakeri 
4112b730f8bSJeremy L Thompson static int Min3(const PetscInt a[3]) { return PetscMin(a[0], PetscMin(a[1], a[2])); }
412de1229c5Srezgarshakeri 
413de1229c5Srezgarshakeri // -----------------------------------------------------------------------------
414de1229c5Srezgarshakeri // Create distribute dm
415de1229c5Srezgarshakeri // -----------------------------------------------------------------------------
416de1229c5Srezgarshakeri PetscErrorCode CreateDistributedDM(RunParams rp, DM *dm) {
417de1229c5Srezgarshakeri   PetscFunctionBeginUser;
418de1229c5Srezgarshakeri   // Setup DM
419de1229c5Srezgarshakeri   if (rp->read_mesh) {
4202b730f8bSJeremy L Thompson     PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, rp->filename, NULL, PETSC_TRUE, dm));
421de1229c5Srezgarshakeri   } else {
422de1229c5Srezgarshakeri     if (rp->user_l_nodes) {
423de1229c5Srezgarshakeri       // Find a nicely composite number of elements no less than global nodes
424de1229c5Srezgarshakeri       PetscMPIInt size;
4252b730f8bSJeremy L Thompson       PetscCall(MPI_Comm_size(rp->comm, &size));
4262b730f8bSJeremy L Thompson       for (PetscInt g_elem = PetscMax(1, size * rp->local_nodes / PetscPowInt(rp->degree, rp->dim));; g_elem++) {
427de1229c5Srezgarshakeri         Split3(g_elem, rp->mesh_elem, true);
428de1229c5Srezgarshakeri         if (Max3(rp->mesh_elem) / Min3(rp->mesh_elem) <= 2) break;
429de1229c5Srezgarshakeri       }
430de1229c5Srezgarshakeri     }
431ab1ff315Srezgarshakeri 
4322b730f8bSJeremy L Thompson     PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, rp->dim, rp->simplex, rp->mesh_elem, NULL, NULL, NULL, PETSC_TRUE, dm));
433de1229c5Srezgarshakeri   }
434de1229c5Srezgarshakeri 
4352b730f8bSJeremy L Thompson   PetscCall(DMSetFromOptions(*dm));
4362b730f8bSJeremy L Thompson   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
437de1229c5Srezgarshakeri 
438ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
439de1229c5Srezgarshakeri }
440de1229c5Srezgarshakeri 
441de1229c5Srezgarshakeri // -----------------------------------------------------------------------------
442