xref: /libCEED/examples/petsc/src/petscutils.c (revision 129d873643707da2348b9f89b9ef383c571b232a)
1e83e87a5Sjeremylt #include "../include/petscutils.h"
2e83e87a5Sjeremylt 
3e83e87a5Sjeremylt // -----------------------------------------------------------------------------
4e83e87a5Sjeremylt // Convert PETSc MemType to libCEED MemType
5e83e87a5Sjeremylt // -----------------------------------------------------------------------------
69b072555Sjeremylt CeedMemType MemTypeP2C(PetscMemType mem_type) {
79b072555Sjeremylt   return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST;
8e83e87a5Sjeremylt }
9e83e87a5Sjeremylt 
10e83e87a5Sjeremylt // -----------------------------------------------------------------------------
112c58efb6Sjeremylt // Apply 3D Kershaw mesh transformation
122c58efb6Sjeremylt // -----------------------------------------------------------------------------
132c58efb6Sjeremylt // Transition from a value of "a" for x=0, to a value of "b" for x=1.  Optionally
142c58efb6Sjeremylt // smooth -- see the commented versions at the end.
152c58efb6Sjeremylt static double step(const double a, const double b, double x) {
162c58efb6Sjeremylt   if (x <= 0) return a;
172c58efb6Sjeremylt   if (x >= 1) return b;
182c58efb6Sjeremylt   return a + (b-a) * (x);
192c58efb6Sjeremylt }
202c58efb6Sjeremylt 
212c58efb6Sjeremylt // 1D transformation at the right boundary
222c58efb6Sjeremylt static double right(const double eps, const double x) {
232c58efb6Sjeremylt   return (x <= 0.5) ? (2-eps) * x : 1 + eps*(x-1);
242c58efb6Sjeremylt }
252c58efb6Sjeremylt 
262c58efb6Sjeremylt // 1D transformation at the left boundary
272c58efb6Sjeremylt static double left(const double eps, const double x) {
282c58efb6Sjeremylt   return 1-right(eps,1-x);
292c58efb6Sjeremylt }
302c58efb6Sjeremylt 
312c58efb6Sjeremylt // Apply 3D Kershaw mesh transformation
322c58efb6Sjeremylt // The eps parameters are in (0, 1]
332c58efb6Sjeremylt // Uniform mesh is recovered for eps=1
349b072555Sjeremylt PetscErrorCode Kershaw(DM dm_orig, PetscScalar eps) {
352c58efb6Sjeremylt   PetscErrorCode ierr;
362c58efb6Sjeremylt   Vec coord;
372c58efb6Sjeremylt   PetscInt ncoord;
382c58efb6Sjeremylt   PetscScalar *c;
392c58efb6Sjeremylt 
402c58efb6Sjeremylt   PetscFunctionBeginUser;
419b072555Sjeremylt   ierr = DMGetCoordinatesLocal(dm_orig, &coord); CHKERRQ(ierr);
422c58efb6Sjeremylt   ierr = VecGetLocalSize(coord, &ncoord); CHKERRQ(ierr);
432c58efb6Sjeremylt   ierr = VecGetArray(coord, &c); CHKERRQ(ierr);
442c58efb6Sjeremylt 
452c58efb6Sjeremylt   for (PetscInt i = 0; i < ncoord; i += 3) {
462c58efb6Sjeremylt     PetscScalar x = c[i], y = c[i+1], z = c[i+2];
472c58efb6Sjeremylt     PetscInt layer = x*6;
482c58efb6Sjeremylt     PetscScalar lambda = (x-layer/6.0)*6;
492c58efb6Sjeremylt     c[i] = x;
502c58efb6Sjeremylt 
512c58efb6Sjeremylt     switch (layer) {
522c58efb6Sjeremylt     case 0:
532c58efb6Sjeremylt       c[i+1] = left(eps, y);
542c58efb6Sjeremylt       c[i+2] = left(eps, z);
552c58efb6Sjeremylt       break;
562c58efb6Sjeremylt     case 1:
572c58efb6Sjeremylt     case 4:
582c58efb6Sjeremylt       c[i+1] = step(left(eps, y), right(eps, y), lambda);
592c58efb6Sjeremylt       c[i+2] = step(left(eps, z), right(eps, z), lambda);
602c58efb6Sjeremylt       break;
612c58efb6Sjeremylt     case 2:
622c58efb6Sjeremylt       c[i+1] = step(right(eps, y), left(eps, y), lambda/2);
632c58efb6Sjeremylt       c[i+2] = step(right(eps, z), left(eps, z), lambda/2);
642c58efb6Sjeremylt       break;
652c58efb6Sjeremylt     case 3:
662c58efb6Sjeremylt       c[i+1] = step(right(eps, y), left(eps, y), (1+lambda)/2);
672c58efb6Sjeremylt       c[i+2] = step(right(eps, z), left(eps, z), (1+lambda)/2);
682c58efb6Sjeremylt       break;
692c58efb6Sjeremylt     default:
702c58efb6Sjeremylt       c[i+1] = right(eps, y);
712c58efb6Sjeremylt       c[i+2] = right(eps, z);
722c58efb6Sjeremylt     }
732c58efb6Sjeremylt   }
742c58efb6Sjeremylt   ierr = VecRestoreArray(coord, &c); CHKERRQ(ierr);
752c58efb6Sjeremylt   PetscFunctionReturn(0);
762c58efb6Sjeremylt }
772c58efb6Sjeremylt 
782c58efb6Sjeremylt // -----------------------------------------------------------------------------
79e83e87a5Sjeremylt // Create BC label
80e83e87a5Sjeremylt // -----------------------------------------------------------------------------
81e83e87a5Sjeremylt static PetscErrorCode CreateBCLabel(DM dm, const char name[]) {
82e83e87a5Sjeremylt   int ierr;
83e83e87a5Sjeremylt   DMLabel label;
84e83e87a5Sjeremylt 
85e83e87a5Sjeremylt   PetscFunctionBeginUser;
86e83e87a5Sjeremylt 
87e83e87a5Sjeremylt   ierr = DMCreateLabel(dm, name); CHKERRQ(ierr);
88e83e87a5Sjeremylt   ierr = DMGetLabel(dm, name, &label); CHKERRQ(ierr);
89e83e87a5Sjeremylt   ierr = DMPlexMarkBoundaryFaces(dm, 1, label); CHKERRQ(ierr);
90e83e87a5Sjeremylt 
91e83e87a5Sjeremylt   PetscFunctionReturn(0);
92e83e87a5Sjeremylt };
93e83e87a5Sjeremylt 
94e83e87a5Sjeremylt // -----------------------------------------------------------------------------
95e83e87a5Sjeremylt // This function sets up a DM for a given degree
96e83e87a5Sjeremylt // -----------------------------------------------------------------------------
979b072555Sjeremylt PetscErrorCode SetupDMByDegree(DM dm, PetscInt degree, PetscInt num_comp_u,
989b072555Sjeremylt                                PetscInt dim, bool enforce_bc, BCFunction bc_func) {
99e83e87a5Sjeremylt   PetscInt ierr, marker_ids[1] = {1};
100e83e87a5Sjeremylt   PetscFE fe;
10165233696SJeremy L Thompson   MPI_Comm comm;
102*129d8736Srezgarshakeri   PetscBool      is_simplex = PETSC_TRUE;
103e83e87a5Sjeremylt 
104e83e87a5Sjeremylt   PetscFunctionBeginUser;
105e83e87a5Sjeremylt 
106*129d8736Srezgarshakeri   // Check if simplex or tensor-product mesh
107*129d8736Srezgarshakeri   ierr = DMPlexIsSimplex(dm, &is_simplex); CHKERRQ(ierr);
108e83e87a5Sjeremylt   // Setup FE
10965233696SJeremy L Thompson   ierr = PetscObjectGetComm((PetscObject)dm, &comm); CHKERRQ(ierr);
110*129d8736Srezgarshakeri   ierr = PetscFECreateLagrange(comm, dim, num_comp_u, is_simplex, degree, degree,
11165233696SJeremy L Thompson                                &fe); CHKERRQ(ierr);
112e83e87a5Sjeremylt   ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr);
113*129d8736Srezgarshakeri   ierr = DMCreateDS(dm); CHKERRQ(ierr);
1147ed3e4cdSJeremy L Thompson   {
115*129d8736Srezgarshakeri     DM             dm_coord;
116*129d8736Srezgarshakeri     PetscDS        ds_coord;
117*129d8736Srezgarshakeri     PetscFE        fe_coord_current, fe_coord_new;
118*129d8736Srezgarshakeri     PetscDualSpace fe_coord_dual_space;
119*129d8736Srezgarshakeri     PetscInt       fe_coord_order, num_comp_coord;
120*129d8736Srezgarshakeri 
121*129d8736Srezgarshakeri     ierr = DMGetCoordinateDM(dm, &dm_coord); CHKERRQ(ierr);
122*129d8736Srezgarshakeri     ierr = DMGetCoordinateDim(dm, &num_comp_coord); CHKERRQ(ierr);
123*129d8736Srezgarshakeri     ierr = DMGetRegionDS(dm_coord, NULL, NULL, &ds_coord); CHKERRQ(ierr);
124*129d8736Srezgarshakeri     ierr = PetscDSGetDiscretization(ds_coord, 0, (PetscObject *)&fe_coord_current); CHKERRQ(ierr);
125*129d8736Srezgarshakeri     ierr = PetscFEGetDualSpace(fe_coord_current, &fe_coord_dual_space); CHKERRQ(ierr);
126*129d8736Srezgarshakeri     ierr = PetscDualSpaceGetOrder(fe_coord_dual_space, &fe_coord_order); CHKERRQ(ierr);
127*129d8736Srezgarshakeri 
128*129d8736Srezgarshakeri     // Create FE for coordinates
129*129d8736Srezgarshakeri     ierr = PetscFECreateLagrange(comm, dim, num_comp_coord, is_simplex, fe_coord_order, degree, &fe_coord_new);
130*129d8736Srezgarshakeri     CHKERRQ(ierr);
131*129d8736Srezgarshakeri     ierr = DMProjectCoordinates(dm, fe_coord_new); CHKERRQ(ierr);
132*129d8736Srezgarshakeri     ierr = PetscFEDestroy(&fe_coord_new); CHKERRQ(ierr);
133*129d8736Srezgarshakeri   }
134*129d8736Srezgarshakeri   /*
135*129d8736Srezgarshakeri   {
136*129d8736Srezgarshakeri     // create FE field for coordinates
1377ed3e4cdSJeremy L Thompson     PetscFE fe_coords;
1387ed3e4cdSJeremy L Thompson     PetscInt num_comp_coord;
1397ed3e4cdSJeremy L Thompson     ierr = DMGetCoordinateDim(dm, &num_comp_coord); CHKERRQ(ierr);
140*129d8736Srezgarshakeri     ierr = PetscFECreateLagrange(comm, dim, num_comp_coord, is_simplex, 1, degree,
1417ed3e4cdSJeremy L Thompson                                  &fe_coords); CHKERRQ(ierr);
1427ed3e4cdSJeremy L Thompson     ierr = DMProjectCoordinates(dm, fe_coords); CHKERRQ(ierr);
1437ed3e4cdSJeremy L Thompson     ierr = PetscFEDestroy(&fe_coords); CHKERRQ(ierr);
1447ed3e4cdSJeremy L Thompson   }
145*129d8736Srezgarshakeri   */
146e83e87a5Sjeremylt   // Setup DM
147*129d8736Srezgarshakeri   //ierr = DMCreateDS(dm); CHKERRQ(ierr);
1489b072555Sjeremylt   if (enforce_bc) {
1499b072555Sjeremylt     PetscBool has_label;
1509b072555Sjeremylt     DMHasLabel(dm, "marker", &has_label);
1519b072555Sjeremylt     if (!has_label) {CreateBCLabel(dm, "marker");}
15269f5adf1Sjeremylt     DMLabel label;
15369f5adf1Sjeremylt     ierr = DMGetLabel(dm, "marker", &label); CHKERRQ(ierr);
154b8962995SJeremy L Thompson     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1,
155b8962995SJeremy L Thompson                          marker_ids, 0, 0, NULL, (void(*)(void))bc_func,
156b8962995SJeremy L Thompson                          NULL, NULL, NULL); CHKERRQ(ierr);
157e83e87a5Sjeremylt   }
158*129d8736Srezgarshakeri 
159*129d8736Srezgarshakeri   if (!is_simplex) {
160*129d8736Srezgarshakeri     DM dm_coord;
161*129d8736Srezgarshakeri     ierr = DMGetCoordinateDM(dm, &dm_coord); CHKERRQ(ierr);
162*129d8736Srezgarshakeri     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL); CHKERRQ(ierr);
163*129d8736Srezgarshakeri     ierr = DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL); CHKERRQ(ierr);
164*129d8736Srezgarshakeri   }
165e83e87a5Sjeremylt   ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
166e83e87a5Sjeremylt 
167e83e87a5Sjeremylt   PetscFunctionReturn(0);
168e83e87a5Sjeremylt };
169e83e87a5Sjeremylt 
170e83e87a5Sjeremylt // -----------------------------------------------------------------------------
171e83e87a5Sjeremylt // Utility function - essential BC dofs are encoded in closure indices as -(i+1)
172e83e87a5Sjeremylt // -----------------------------------------------------------------------------
173e83e87a5Sjeremylt PetscInt Involute(PetscInt i) {
174e83e87a5Sjeremylt   return i >= 0 ? i : -(i + 1);
175e83e87a5Sjeremylt };
176e83e87a5Sjeremylt 
177e83e87a5Sjeremylt // -----------------------------------------------------------------------------
178e83e87a5Sjeremylt // Get CEED restriction data from DMPlex
179e83e87a5Sjeremylt // -----------------------------------------------------------------------------
1807ed3e4cdSJeremy L Thompson PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height,
1817ed3e4cdSJeremy L Thompson     DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr) {
1827ed3e4cdSJeremy L Thompson   PetscInt num_elem, elem_size, num_dof, num_comp, *elem_restr_offsets;
183e83e87a5Sjeremylt   PetscErrorCode ierr;
184e83e87a5Sjeremylt 
185e83e87a5Sjeremylt   PetscFunctionBeginUser;
186e83e87a5Sjeremylt 
1877ed3e4cdSJeremy L Thompson   ierr = DMPlexGetLocalOffsets(dm, domain_label, value, height, 0, &num_elem,
1887ed3e4cdSJeremy L Thompson                                &elem_size, &num_comp, &num_dof, &elem_restr_offsets);
1897ed3e4cdSJeremy L Thompson   CHKERRQ(ierr);
190e83e87a5Sjeremylt 
1917ed3e4cdSJeremy L Thompson   CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp,
1927ed3e4cdSJeremy L Thompson                             1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES,
1939b072555Sjeremylt                             elem_restr_offsets, elem_restr);
1949b072555Sjeremylt   ierr = PetscFree(elem_restr_offsets); CHKERRQ(ierr);
1957ed3e4cdSJeremy L Thompson 
196e83e87a5Sjeremylt   PetscFunctionReturn(0);
197e83e87a5Sjeremylt };
198e83e87a5Sjeremylt 
199e83e87a5Sjeremylt // -----------------------------------------------------------------------------
200*129d8736Srezgarshakeri // Utility function - convert from DMPolytopeType to CeedElemTopology
201*129d8736Srezgarshakeri // -----------------------------------------------------------------------------
202*129d8736Srezgarshakeri static inline CeedElemTopology ElemTopologyP2C(DMPolytopeType cell_type) {
203*129d8736Srezgarshakeri   switch (cell_type) {
204*129d8736Srezgarshakeri   case DM_POLYTOPE_TRIANGLE:      return CEED_TOPOLOGY_TRIANGLE;
205*129d8736Srezgarshakeri   case DM_POLYTOPE_QUADRILATERAL: return CEED_TOPOLOGY_QUAD;
206*129d8736Srezgarshakeri   case DM_POLYTOPE_TETRAHEDRON:   return CEED_TOPOLOGY_TET;
207*129d8736Srezgarshakeri   case DM_POLYTOPE_HEXAHEDRON:    return CEED_TOPOLOGY_HEX;
208*129d8736Srezgarshakeri   default:                        return 0;
209*129d8736Srezgarshakeri   }
210*129d8736Srezgarshakeri }
211*129d8736Srezgarshakeri 
212*129d8736Srezgarshakeri // -----------------------------------------------------------------------------
213*129d8736Srezgarshakeri // Get CEED Basis from DMPlex
214*129d8736Srezgarshakeri // -----------------------------------------------------------------------------
215*129d8736Srezgarshakeri PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label,
216*129d8736Srezgarshakeri                                    CeedInt label_value, CeedInt height,
217*129d8736Srezgarshakeri                                    CeedInt dm_field, CeedBasis *basis) {
218*129d8736Srezgarshakeri   PetscErrorCode   ierr;
219*129d8736Srezgarshakeri   PetscDS          ds;
220*129d8736Srezgarshakeri   PetscFE          fe;
221*129d8736Srezgarshakeri   PetscQuadrature  quadrature;
222*129d8736Srezgarshakeri   PetscBool        is_simplex = PETSC_TRUE;
223*129d8736Srezgarshakeri   PetscInt         dim, ds_field = -1, num_comp, P, Q;
224*129d8736Srezgarshakeri 
225*129d8736Srezgarshakeri   PetscFunctionBeginUser;
226*129d8736Srezgarshakeri 
227*129d8736Srezgarshakeri   // Get basis information
228*129d8736Srezgarshakeri   {
229*129d8736Srezgarshakeri     IS             field_is;
230*129d8736Srezgarshakeri     const PetscInt *fields;
231*129d8736Srezgarshakeri     PetscInt       num_fields;
232*129d8736Srezgarshakeri 
233*129d8736Srezgarshakeri     ierr = DMGetRegionDS(dm, domain_label, &field_is, &ds); CHKERRQ(ierr);
234*129d8736Srezgarshakeri     // Translate dm_field to ds_field
235*129d8736Srezgarshakeri     ierr = ISGetIndices(field_is, &fields); CHKERRQ(ierr);
236*129d8736Srezgarshakeri     ierr = ISGetSize(field_is, &num_fields); CHKERRQ(ierr);
237*129d8736Srezgarshakeri     for (PetscInt i = 0; i < num_fields; i++) {
238*129d8736Srezgarshakeri       if (dm_field == fields[i]) {
239*129d8736Srezgarshakeri         ds_field = i;
240*129d8736Srezgarshakeri         break;
241*129d8736Srezgarshakeri       }
242*129d8736Srezgarshakeri     }
243*129d8736Srezgarshakeri     ierr = ISRestoreIndices(field_is, &fields); CHKERRQ(ierr);
244*129d8736Srezgarshakeri   }
245*129d8736Srezgarshakeri   if (ds_field == -1) {
246*129d8736Srezgarshakeri     // LCOV_EXCL_START
247*129d8736Srezgarshakeri     SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP,
248*129d8736Srezgarshakeri              "Could not find dm_field %D in DS", dm_field);
249*129d8736Srezgarshakeri     // LCOV_EXCL_STOP
250*129d8736Srezgarshakeri   }
251*129d8736Srezgarshakeri 
252*129d8736Srezgarshakeri   // Get element information
253*129d8736Srezgarshakeri   {
254*129d8736Srezgarshakeri     PetscDualSpace dual_space;
255*129d8736Srezgarshakeri     PetscInt       num_dual_basis_vectors;
256*129d8736Srezgarshakeri 
257*129d8736Srezgarshakeri     ierr = PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe);
258*129d8736Srezgarshakeri     CHKERRQ(ierr);
259*129d8736Srezgarshakeri     ierr = PetscFEGetHeightSubspace(fe, height, &fe); CHKERRQ(ierr);
260*129d8736Srezgarshakeri     ierr = PetscFEGetSpatialDimension(fe, &dim); CHKERRQ(ierr);
261*129d8736Srezgarshakeri     ierr = PetscFEGetNumComponents(fe, &num_comp); CHKERRQ(ierr);
262*129d8736Srezgarshakeri     ierr = PetscFEGetDualSpace(fe, &dual_space); CHKERRQ(ierr);
263*129d8736Srezgarshakeri     ierr = PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors);
264*129d8736Srezgarshakeri     CHKERRQ(ierr);
265*129d8736Srezgarshakeri     P = num_dual_basis_vectors / num_comp;
266*129d8736Srezgarshakeri     ierr = PetscFEGetQuadrature(fe, &quadrature); CHKERRQ(ierr);
267*129d8736Srezgarshakeri     ierr = PetscQuadratureGetData(quadrature, NULL, NULL, &Q, NULL, NULL);
268*129d8736Srezgarshakeri     CHKERRQ(ierr);
269*129d8736Srezgarshakeri   }
270*129d8736Srezgarshakeri 
271*129d8736Srezgarshakeri   // Check if simplex or tensor-product mesh
272*129d8736Srezgarshakeri   ierr = DMPlexIsSimplex(dm, &is_simplex); CHKERRQ(ierr);
273*129d8736Srezgarshakeri   // Build libCEED basis
274*129d8736Srezgarshakeri   if (is_simplex) {
275*129d8736Srezgarshakeri     PetscInt          num_derivatives = 1, first_point;
276*129d8736Srezgarshakeri     PetscInt          ids[1] = {label_value};
277*129d8736Srezgarshakeri     PetscTabulation   basis_tabulation;
278*129d8736Srezgarshakeri     const PetscScalar *q_points, *q_weights;
279*129d8736Srezgarshakeri     DMLabel           depth_label;
280*129d8736Srezgarshakeri     DMPolytopeType    cell_type;
281*129d8736Srezgarshakeri     CeedElemTopology  elem_topo;
282*129d8736Srezgarshakeri     PetscScalar       *interp, *grad;
283*129d8736Srezgarshakeri 
284*129d8736Srezgarshakeri     // Use depth label if no domain label present
285*129d8736Srezgarshakeri     if (!domain_label) {
286*129d8736Srezgarshakeri       PetscInt depth;
287*129d8736Srezgarshakeri 
288*129d8736Srezgarshakeri       ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr);
289*129d8736Srezgarshakeri       ierr = DMPlexGetDepthLabel(dm, &depth_label); CHKERRQ(ierr);
290*129d8736Srezgarshakeri       ids[0] = depth - height;
291*129d8736Srezgarshakeri     }
292*129d8736Srezgarshakeri     // Get cell interp, grad, and quadrature data
293*129d8736Srezgarshakeri     ierr = PetscFEGetCellTabulation(fe, num_derivatives, &basis_tabulation);
294*129d8736Srezgarshakeri     CHKERRQ(ierr);
295*129d8736Srezgarshakeri     ierr = PetscQuadratureGetData(quadrature, NULL, NULL, NULL, &q_points,
296*129d8736Srezgarshakeri                                   &q_weights); CHKERRQ(ierr);
297*129d8736Srezgarshakeri     ierr = DMGetFirstLabeledPoint(dm, dm, domain_label ? domain_label : depth_label,
298*129d8736Srezgarshakeri                                   1, ids, height, &first_point, NULL);
299*129d8736Srezgarshakeri     CHKERRQ(ierr);
300*129d8736Srezgarshakeri     ierr = DMPlexGetCellType(dm, first_point, &cell_type); CHKERRQ(ierr);
301*129d8736Srezgarshakeri     elem_topo = ElemTopologyP2C(cell_type);
302*129d8736Srezgarshakeri     if (!elem_topo) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP,
303*129d8736Srezgarshakeri                               "DMPlex topology not supported");
304*129d8736Srezgarshakeri     // Convert to libCEED orientation
305*129d8736Srezgarshakeri     ierr = PetscCalloc(P * Q * sizeof(PetscScalar), &interp); CHKERRQ(ierr);
306*129d8736Srezgarshakeri     ierr = PetscCalloc(P * Q * dim * sizeof(PetscScalar), &grad); CHKERRQ(ierr);
307*129d8736Srezgarshakeri     const CeedInt c = 0;
308*129d8736Srezgarshakeri     for (CeedInt q = 0; q < Q; q++) {
309*129d8736Srezgarshakeri       for (CeedInt p = 0; p < P; p++) {
310*129d8736Srezgarshakeri         interp[q*P + p] = basis_tabulation->T[0][(q*P + p)*num_comp*num_comp + c];
311*129d8736Srezgarshakeri         for (CeedInt d = 0; d < dim; d++) {
312*129d8736Srezgarshakeri           grad[(d*Q + q)*P + p] = basis_tabulation->T[1][((q*P + p)*num_comp*num_comp + c)
313*129d8736Srezgarshakeri                                   *dim + d];
314*129d8736Srezgarshakeri         }
315*129d8736Srezgarshakeri       }
316*129d8736Srezgarshakeri     }
317*129d8736Srezgarshakeri     // Finaly, create libCEED basis
318*129d8736Srezgarshakeri     ierr = CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad,
319*129d8736Srezgarshakeri                              q_points, q_weights, basis);
320*129d8736Srezgarshakeri     CHKERRQ(ierr);
321*129d8736Srezgarshakeri     ierr = PetscFree(interp); CHKERRQ(ierr);
322*129d8736Srezgarshakeri     ierr = PetscFree(grad); CHKERRQ(ierr);
323*129d8736Srezgarshakeri   } else {
324*129d8736Srezgarshakeri     CeedInt P_1d = (CeedInt) round(pow(P, 1.0 / dim));
325*129d8736Srezgarshakeri     CeedInt Q_1d = (CeedInt) round(pow(Q, 1.0 / dim));
326*129d8736Srezgarshakeri 
327*129d8736Srezgarshakeri     ierr = CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp, P_1d, Q_1d,
328*129d8736Srezgarshakeri                                            CEED_GAUSS, basis);
329*129d8736Srezgarshakeri     CHKERRQ(ierr);
330*129d8736Srezgarshakeri   }
331*129d8736Srezgarshakeri 
332*129d8736Srezgarshakeri   PetscFunctionReturn(0);
333*129d8736Srezgarshakeri };
334*129d8736Srezgarshakeri 
335*129d8736Srezgarshakeri // -----------------------------------------------------------------------------
336