xref: /libCEED/examples/petsc/src/petscutils.c (revision 69f5adf13c38b25ccc4fa9a9f050d7775cff5bb5)
1e83e87a5Sjeremylt #include "../include/petscutils.h"
2*69f5adf1Sjeremylt #include "../include/petscmacros.h"
3e83e87a5Sjeremylt 
4e83e87a5Sjeremylt // -----------------------------------------------------------------------------
5e83e87a5Sjeremylt // Convert PETSc MemType to libCEED MemType
6e83e87a5Sjeremylt // -----------------------------------------------------------------------------
79b072555Sjeremylt CeedMemType MemTypeP2C(PetscMemType mem_type) {
89b072555Sjeremylt   return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST;
9e83e87a5Sjeremylt }
10e83e87a5Sjeremylt 
11e83e87a5Sjeremylt // -----------------------------------------------------------------------------
12e83e87a5Sjeremylt // Utility function taken from petsc/src/dm/impls/plex/examples/tutorials/ex7.c
13e83e87a5Sjeremylt // -----------------------------------------------------------------------------
14e83e87a5Sjeremylt PetscErrorCode ProjectToUnitSphere(DM dm) {
15e83e87a5Sjeremylt   Vec            coordinates;
16e83e87a5Sjeremylt   PetscScalar   *coords;
17e83e87a5Sjeremylt   PetscInt       Nv, v, dim, d;
18e83e87a5Sjeremylt   PetscErrorCode ierr;
19e83e87a5Sjeremylt 
20e83e87a5Sjeremylt   PetscFunctionBeginUser;
21e83e87a5Sjeremylt   ierr = DMGetCoordinatesLocal(dm, &coordinates); CHKERRQ(ierr);
22e83e87a5Sjeremylt   ierr = VecGetLocalSize(coordinates, &Nv); CHKERRQ(ierr);
23e83e87a5Sjeremylt   ierr = VecGetBlockSize(coordinates, &dim); CHKERRQ(ierr);
24e83e87a5Sjeremylt   Nv  /= dim;
25e83e87a5Sjeremylt   ierr = VecGetArray(coordinates, &coords); CHKERRQ(ierr);
26e83e87a5Sjeremylt   for (v = 0; v < Nv; ++v) {
27e83e87a5Sjeremylt     PetscReal r = 0.0;
28e83e87a5Sjeremylt 
29e83e87a5Sjeremylt     for (d = 0; d < dim; ++d) r += PetscSqr(PetscRealPart(coords[v*dim+d]));
30e83e87a5Sjeremylt     r = PetscSqrtReal(r);
31e83e87a5Sjeremylt     for (d = 0; d < dim; ++d) coords[v*dim+d] /= r;
32e83e87a5Sjeremylt   }
33e83e87a5Sjeremylt   ierr = VecRestoreArray(coordinates, &coords); CHKERRQ(ierr);
34e83e87a5Sjeremylt   PetscFunctionReturn(0);
35e83e87a5Sjeremylt };
36e83e87a5Sjeremylt 
37e83e87a5Sjeremylt // -----------------------------------------------------------------------------
382c58efb6Sjeremylt // Apply 3D Kershaw mesh transformation
392c58efb6Sjeremylt // -----------------------------------------------------------------------------
402c58efb6Sjeremylt // Transition from a value of "a" for x=0, to a value of "b" for x=1.  Optionally
412c58efb6Sjeremylt // smooth -- see the commented versions at the end.
422c58efb6Sjeremylt static double step(const double a, const double b, double x) {
432c58efb6Sjeremylt   if (x <= 0) return a;
442c58efb6Sjeremylt   if (x >= 1) return b;
452c58efb6Sjeremylt   return a + (b-a) * (x);
462c58efb6Sjeremylt }
472c58efb6Sjeremylt 
482c58efb6Sjeremylt // 1D transformation at the right boundary
492c58efb6Sjeremylt static double right(const double eps, const double x) {
502c58efb6Sjeremylt   return (x <= 0.5) ? (2-eps) * x : 1 + eps*(x-1);
512c58efb6Sjeremylt }
522c58efb6Sjeremylt 
532c58efb6Sjeremylt // 1D transformation at the left boundary
542c58efb6Sjeremylt static double left(const double eps, const double x) {
552c58efb6Sjeremylt   return 1-right(eps,1-x);
562c58efb6Sjeremylt }
572c58efb6Sjeremylt 
582c58efb6Sjeremylt // Apply 3D Kershaw mesh transformation
592c58efb6Sjeremylt // The eps parameters are in (0, 1]
602c58efb6Sjeremylt // Uniform mesh is recovered for eps=1
619b072555Sjeremylt PetscErrorCode Kershaw(DM dm_orig, PetscScalar eps) {
622c58efb6Sjeremylt   PetscErrorCode ierr;
632c58efb6Sjeremylt   Vec coord;
642c58efb6Sjeremylt   PetscInt ncoord;
652c58efb6Sjeremylt   PetscScalar *c;
662c58efb6Sjeremylt 
672c58efb6Sjeremylt   PetscFunctionBeginUser;
689b072555Sjeremylt   ierr = DMGetCoordinatesLocal(dm_orig, &coord); CHKERRQ(ierr);
692c58efb6Sjeremylt   ierr = VecGetLocalSize(coord, &ncoord); CHKERRQ(ierr);
702c58efb6Sjeremylt   ierr = VecGetArray(coord, &c); CHKERRQ(ierr);
712c58efb6Sjeremylt 
722c58efb6Sjeremylt   for (PetscInt i = 0; i < ncoord; i += 3) {
732c58efb6Sjeremylt     PetscScalar x = c[i], y = c[i+1], z = c[i+2];
742c58efb6Sjeremylt     PetscInt layer = x*6;
752c58efb6Sjeremylt     PetscScalar lambda = (x-layer/6.0)*6;
762c58efb6Sjeremylt     c[i] = x;
772c58efb6Sjeremylt 
782c58efb6Sjeremylt     switch (layer) {
792c58efb6Sjeremylt     case 0:
802c58efb6Sjeremylt       c[i+1] = left(eps, y);
812c58efb6Sjeremylt       c[i+2] = left(eps, z);
822c58efb6Sjeremylt       break;
832c58efb6Sjeremylt     case 1:
842c58efb6Sjeremylt     case 4:
852c58efb6Sjeremylt       c[i+1] = step(left(eps, y), right(eps, y), lambda);
862c58efb6Sjeremylt       c[i+2] = step(left(eps, z), right(eps, z), lambda);
872c58efb6Sjeremylt       break;
882c58efb6Sjeremylt     case 2:
892c58efb6Sjeremylt       c[i+1] = step(right(eps, y), left(eps, y), lambda/2);
902c58efb6Sjeremylt       c[i+2] = step(right(eps, z), left(eps, z), lambda/2);
912c58efb6Sjeremylt       break;
922c58efb6Sjeremylt     case 3:
932c58efb6Sjeremylt       c[i+1] = step(right(eps, y), left(eps, y), (1+lambda)/2);
942c58efb6Sjeremylt       c[i+2] = step(right(eps, z), left(eps, z), (1+lambda)/2);
952c58efb6Sjeremylt       break;
962c58efb6Sjeremylt     default:
972c58efb6Sjeremylt       c[i+1] = right(eps, y);
982c58efb6Sjeremylt       c[i+2] = right(eps, z);
992c58efb6Sjeremylt     }
1002c58efb6Sjeremylt   }
1012c58efb6Sjeremylt   ierr = VecRestoreArray(coord, &c); CHKERRQ(ierr);
1022c58efb6Sjeremylt   PetscFunctionReturn(0);
1032c58efb6Sjeremylt }
1042c58efb6Sjeremylt 
1052c58efb6Sjeremylt // -----------------------------------------------------------------------------
106e83e87a5Sjeremylt // PETSc FE Boilerplate
107e83e87a5Sjeremylt // -----------------------------------------------------------------------------
108e83e87a5Sjeremylt PetscErrorCode PetscFECreateByDegree(DM dm, PetscInt dim, PetscInt Nc,
109e83e87a5Sjeremylt                                      PetscBool isSimplex, const char prefix[],
110e83e87a5Sjeremylt                                      PetscInt order, PetscFE *fem) {
111e83e87a5Sjeremylt   PetscQuadrature q, fq;
112e83e87a5Sjeremylt   DM              K;
113e83e87a5Sjeremylt   PetscSpace      P;
114e83e87a5Sjeremylt   PetscDualSpace  Q;
115e83e87a5Sjeremylt   PetscInt        quadPointsPerEdge;
116e83e87a5Sjeremylt   PetscBool       tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE;
117e83e87a5Sjeremylt   PetscErrorCode  ierr;
118e83e87a5Sjeremylt 
119e83e87a5Sjeremylt   PetscFunctionBeginUser;
120e83e87a5Sjeremylt   /* Create space */
121e83e87a5Sjeremylt   ierr = PetscSpaceCreate(PetscObjectComm((PetscObject) dm), &P); CHKERRQ(ierr);
122e83e87a5Sjeremylt   ierr = PetscObjectSetOptionsPrefix((PetscObject) P, prefix); CHKERRQ(ierr);
123e83e87a5Sjeremylt   ierr = PetscSpacePolynomialSetTensor(P, tensor); CHKERRQ(ierr);
124e83e87a5Sjeremylt   ierr = PetscSpaceSetFromOptions(P); CHKERRQ(ierr);
125e83e87a5Sjeremylt   ierr = PetscSpaceSetNumComponents(P, Nc); CHKERRQ(ierr);
126e83e87a5Sjeremylt   ierr = PetscSpaceSetNumVariables(P, dim); CHKERRQ(ierr);
127e83e87a5Sjeremylt   ierr = PetscSpaceSetDegree(P, order, order); CHKERRQ(ierr);
128e83e87a5Sjeremylt   ierr = PetscSpaceSetUp(P); CHKERRQ(ierr);
129e83e87a5Sjeremylt   ierr = PetscSpacePolynomialGetTensor(P, &tensor); CHKERRQ(ierr);
130e83e87a5Sjeremylt   /* Create dual space */
131e83e87a5Sjeremylt   ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) dm), &Q);
132e83e87a5Sjeremylt   CHKERRQ(ierr);
133e83e87a5Sjeremylt   ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE); CHKERRQ(ierr);
134e83e87a5Sjeremylt   ierr = PetscObjectSetOptionsPrefix((PetscObject) Q, prefix); CHKERRQ(ierr);
135e83e87a5Sjeremylt   ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K); CHKERRQ(ierr);
136e83e87a5Sjeremylt   ierr = PetscDualSpaceSetDM(Q, K); CHKERRQ(ierr);
137e83e87a5Sjeremylt   ierr = DMDestroy(&K); CHKERRQ(ierr);
138e83e87a5Sjeremylt   ierr = PetscDualSpaceSetNumComponents(Q, Nc); CHKERRQ(ierr);
139e83e87a5Sjeremylt   ierr = PetscDualSpaceSetOrder(Q, order); CHKERRQ(ierr);
140e83e87a5Sjeremylt   ierr = PetscDualSpaceLagrangeSetTensor(Q, tensor); CHKERRQ(ierr);
141e83e87a5Sjeremylt   ierr = PetscDualSpaceSetFromOptions(Q); CHKERRQ(ierr);
142e83e87a5Sjeremylt   ierr = PetscDualSpaceSetUp(Q); CHKERRQ(ierr);
143e83e87a5Sjeremylt   /* Create element */
144e83e87a5Sjeremylt   ierr = PetscFECreate(PetscObjectComm((PetscObject) dm), fem); CHKERRQ(ierr);
145e83e87a5Sjeremylt   ierr = PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix); CHKERRQ(ierr);
146e83e87a5Sjeremylt   ierr = PetscFESetFromOptions(*fem); CHKERRQ(ierr);
147e83e87a5Sjeremylt   ierr = PetscFESetBasisSpace(*fem, P); CHKERRQ(ierr);
148e83e87a5Sjeremylt   ierr = PetscFESetDualSpace(*fem, Q); CHKERRQ(ierr);
149e83e87a5Sjeremylt   ierr = PetscFESetNumComponents(*fem, Nc); CHKERRQ(ierr);
150e83e87a5Sjeremylt   ierr = PetscFESetUp(*fem); CHKERRQ(ierr);
151e83e87a5Sjeremylt   ierr = PetscSpaceDestroy(&P); CHKERRQ(ierr);
152e83e87a5Sjeremylt   ierr = PetscDualSpaceDestroy(&Q); CHKERRQ(ierr);
153e83e87a5Sjeremylt   /* Create quadrature */
154e83e87a5Sjeremylt   quadPointsPerEdge = PetscMax(order + 1,1);
155e83e87a5Sjeremylt   if (isSimplex) {
156e83e87a5Sjeremylt     ierr = PetscDTStroudConicalQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0,
157e83e87a5Sjeremylt                                           &q); CHKERRQ(ierr);
158e83e87a5Sjeremylt     ierr = PetscDTStroudConicalQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0,
159e83e87a5Sjeremylt                                           &fq); CHKERRQ(ierr);
160e83e87a5Sjeremylt   } else {
161e83e87a5Sjeremylt     ierr = PetscDTGaussTensorQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0,
162e83e87a5Sjeremylt                                         &q); CHKERRQ(ierr);
163e83e87a5Sjeremylt     ierr = PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0,
164e83e87a5Sjeremylt                                         &fq); CHKERRQ(ierr);
165e83e87a5Sjeremylt   }
166e83e87a5Sjeremylt   ierr = PetscFESetQuadrature(*fem, q); CHKERRQ(ierr);
167e83e87a5Sjeremylt   ierr = PetscFESetFaceQuadrature(*fem, fq); CHKERRQ(ierr);
168e83e87a5Sjeremylt   ierr = PetscQuadratureDestroy(&q); CHKERRQ(ierr);
169e83e87a5Sjeremylt   ierr = PetscQuadratureDestroy(&fq); CHKERRQ(ierr);
170e83e87a5Sjeremylt 
171e83e87a5Sjeremylt   PetscFunctionReturn(0);
172e83e87a5Sjeremylt };
173e83e87a5Sjeremylt 
174e83e87a5Sjeremylt // -----------------------------------------------------------------------------
175e83e87a5Sjeremylt // Create BC label
176e83e87a5Sjeremylt // -----------------------------------------------------------------------------
177e83e87a5Sjeremylt static PetscErrorCode CreateBCLabel(DM dm, const char name[]) {
178e83e87a5Sjeremylt   int ierr;
179e83e87a5Sjeremylt   DMLabel label;
180e83e87a5Sjeremylt 
181e83e87a5Sjeremylt   PetscFunctionBeginUser;
182e83e87a5Sjeremylt 
183e83e87a5Sjeremylt   ierr = DMCreateLabel(dm, name); CHKERRQ(ierr);
184e83e87a5Sjeremylt   ierr = DMGetLabel(dm, name, &label); CHKERRQ(ierr);
185e83e87a5Sjeremylt   ierr = DMPlexMarkBoundaryFaces(dm, 1, label); CHKERRQ(ierr);
186e83e87a5Sjeremylt   ierr = DMPlexLabelComplete(dm, label); CHKERRQ(ierr);
187e83e87a5Sjeremylt 
188e83e87a5Sjeremylt   PetscFunctionReturn(0);
189e83e87a5Sjeremylt };
190e83e87a5Sjeremylt 
191e83e87a5Sjeremylt // -----------------------------------------------------------------------------
192e83e87a5Sjeremylt // This function sets up a DM for a given degree
193e83e87a5Sjeremylt // -----------------------------------------------------------------------------
1949b072555Sjeremylt PetscErrorCode SetupDMByDegree(DM dm, PetscInt degree, PetscInt num_comp_u,
1959b072555Sjeremylt                                PetscInt dim, bool enforce_bc, BCFunction bc_func) {
196e83e87a5Sjeremylt   PetscInt ierr, marker_ids[1] = {1};
197e83e87a5Sjeremylt   PetscFE fe;
198e83e87a5Sjeremylt 
199e83e87a5Sjeremylt   PetscFunctionBeginUser;
200e83e87a5Sjeremylt 
201e83e87a5Sjeremylt   // Setup FE
2029b072555Sjeremylt   ierr = PetscFECreateByDegree(dm, dim, num_comp_u, PETSC_FALSE, NULL, degree,
2039b072555Sjeremylt                                &fe);
204e83e87a5Sjeremylt   CHKERRQ(ierr);
205e83e87a5Sjeremylt   ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
206e83e87a5Sjeremylt   ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr);
207e83e87a5Sjeremylt 
208e83e87a5Sjeremylt   // Setup DM
209e83e87a5Sjeremylt   ierr = DMCreateDS(dm); CHKERRQ(ierr);
2109b072555Sjeremylt   if (enforce_bc) {
2119b072555Sjeremylt     PetscBool has_label;
2129b072555Sjeremylt     DMHasLabel(dm, "marker", &has_label);
2139b072555Sjeremylt     if (!has_label) {CreateBCLabel(dm, "marker");}
214*69f5adf1Sjeremylt     DMLabel label;
215*69f5adf1Sjeremylt     ierr = DMGetLabel(dm, "marker", &label); CHKERRQ(ierr);
216*69f5adf1Sjeremylt     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, "marker", 1,
217*69f5adf1Sjeremylt                          marker_ids, 0, 0, NULL,
218*69f5adf1Sjeremylt                          (void(*)(void))bc_func, NULL, NULL, NULL);
219e83e87a5Sjeremylt     CHKERRQ(ierr);
220e83e87a5Sjeremylt   }
221e83e87a5Sjeremylt   ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
222e83e87a5Sjeremylt   CHKERRQ(ierr);
223e83e87a5Sjeremylt   ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
224e83e87a5Sjeremylt 
225e83e87a5Sjeremylt   PetscFunctionReturn(0);
226e83e87a5Sjeremylt };
227e83e87a5Sjeremylt 
228e83e87a5Sjeremylt // -----------------------------------------------------------------------------
229e83e87a5Sjeremylt // Utility function - essential BC dofs are encoded in closure indices as -(i+1)
230e83e87a5Sjeremylt // -----------------------------------------------------------------------------
231e83e87a5Sjeremylt PetscInt Involute(PetscInt i) {
232e83e87a5Sjeremylt   return i >= 0 ? i : -(i + 1);
233e83e87a5Sjeremylt };
234e83e87a5Sjeremylt 
235e83e87a5Sjeremylt // -----------------------------------------------------------------------------
236e83e87a5Sjeremylt // Get CEED restriction data from DMPlex
237e83e87a5Sjeremylt // -----------------------------------------------------------------------------
238e83e87a5Sjeremylt PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P,
2399b072555Sjeremylt     CeedInt topo_dim, CeedInt height, DMLabel domain_label, CeedInt value,
2409b072555Sjeremylt     CeedElemRestriction *elem_restr) {
241e83e87a5Sjeremylt   PetscSection section;
2429b072555Sjeremylt   PetscInt p, num_elem, num_dof, *elem_restr_offsets, e_offset, num_fields, dim,
2439b072555Sjeremylt            depth;
2449b072555Sjeremylt   DMLabel depth_label;
2459b072555Sjeremylt   IS depth_is, iter_is;
2469b072555Sjeremylt   Vec U_loc;
2479b072555Sjeremylt   const PetscInt *iter_indices;
248e83e87a5Sjeremylt   PetscErrorCode ierr;
249e83e87a5Sjeremylt 
250e83e87a5Sjeremylt   PetscFunctionBeginUser;
251e83e87a5Sjeremylt 
252e83e87a5Sjeremylt   ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
253e83e87a5Sjeremylt   dim -= height;
254e83e87a5Sjeremylt   ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
2559b072555Sjeremylt   ierr = PetscSectionGetNumFields(section, &num_fields); CHKERRQ(ierr);
2569b072555Sjeremylt   PetscInt num_comp[num_fields], field_off[num_fields+1];
2579b072555Sjeremylt   field_off[0] = 0;
2589b072555Sjeremylt   for (PetscInt f = 0; f < num_fields; f++) {
2599b072555Sjeremylt     ierr = PetscSectionGetFieldComponents(section, f, &num_comp[f]); CHKERRQ(ierr);
2609b072555Sjeremylt     field_off[f+1] = field_off[f] + num_comp[f];
261e83e87a5Sjeremylt   }
262e83e87a5Sjeremylt 
263e83e87a5Sjeremylt   ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr);
2649b072555Sjeremylt   ierr = DMPlexGetDepthLabel(dm, &depth_label); CHKERRQ(ierr);
2659b072555Sjeremylt   ierr = DMLabelGetStratumIS(depth_label, depth - height, &depth_is);
2669b072555Sjeremylt   CHKERRQ(ierr);
2679b072555Sjeremylt   if (domain_label) {
2689b072555Sjeremylt     IS domain_is;
2699b072555Sjeremylt     ierr = DMLabelGetStratumIS(domain_label, value, &domain_is); CHKERRQ(ierr);
2709b072555Sjeremylt     if (domain_is) { // domain_is is non-empty
2719b072555Sjeremylt       ierr = ISIntersect(depth_is, domain_is, &iter_is); CHKERRQ(ierr);
2729b072555Sjeremylt       ierr = ISDestroy(&domain_is); CHKERRQ(ierr);
2739b072555Sjeremylt     } else { // domain_is is NULL (empty)
2749b072555Sjeremylt       iter_is = NULL;
275e83e87a5Sjeremylt     }
2769b072555Sjeremylt     ierr = ISDestroy(&depth_is); CHKERRQ(ierr);
277e83e87a5Sjeremylt   } else {
2789b072555Sjeremylt     iter_is = depth_is;
279e83e87a5Sjeremylt   }
2809b072555Sjeremylt   if (iter_is) {
2819b072555Sjeremylt     ierr = ISGetLocalSize(iter_is, &num_elem); CHKERRQ(ierr);
2829b072555Sjeremylt     ierr = ISGetIndices(iter_is, &iter_indices); CHKERRQ(ierr);
283e83e87a5Sjeremylt   } else {
2849b072555Sjeremylt     num_elem = 0;
2859b072555Sjeremylt     iter_indices = NULL;
286e83e87a5Sjeremylt   }
2879b072555Sjeremylt   ierr = PetscMalloc1(num_elem*PetscPowInt(P, topo_dim), &elem_restr_offsets);
2889b072555Sjeremylt   CHKERRQ(ierr);
2899b072555Sjeremylt   for (p = 0, e_offset = 0; p < num_elem; p++) {
2909b072555Sjeremylt     PetscInt c = iter_indices[p];
2919b072555Sjeremylt     PetscInt num_indices, *indices, num_nodes;
292e83e87a5Sjeremylt     ierr = DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE,
2939b072555Sjeremylt                                    &num_indices, &indices, NULL, NULL);
294e83e87a5Sjeremylt     CHKERRQ(ierr);
295e83e87a5Sjeremylt     bool flip = false;
296e83e87a5Sjeremylt     if (height > 0) {
2979b072555Sjeremylt       PetscInt num_cells, num_faces, start = -1;
298e83e87a5Sjeremylt       const PetscInt *orients, *faces, *cells;
299e83e87a5Sjeremylt       ierr = DMPlexGetSupport(dm, c, &cells); CHKERRQ(ierr);
3009b072555Sjeremylt       ierr = DMPlexGetSupportSize(dm, c, &num_cells); CHKERRQ(ierr);
3019b072555Sjeremylt       if (num_cells != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
302e83e87a5Sjeremylt                                      "Expected one cell in support of exterior face, but got %D cells",
3039b072555Sjeremylt                                      num_cells);
304e83e87a5Sjeremylt       ierr = DMPlexGetCone(dm, cells[0], &faces); CHKERRQ(ierr);
3059b072555Sjeremylt       ierr = DMPlexGetConeSize(dm, cells[0], &num_faces); CHKERRQ(ierr);
3069b072555Sjeremylt       for (PetscInt i=0; i<num_faces; i++) {if (faces[i] == c) start = i;}
307e83e87a5Sjeremylt       if (start < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT,
308e83e87a5Sjeremylt                                 "Could not find face %D in cone of its support",
309e83e87a5Sjeremylt                                 c);
310e83e87a5Sjeremylt       ierr = DMPlexGetConeOrientation(dm, cells[0], &orients); CHKERRQ(ierr);
311e83e87a5Sjeremylt       if (orients[start] < 0) flip = true;
312e83e87a5Sjeremylt     }
3139b072555Sjeremylt     if (num_indices % field_off[num_fields]) SETERRQ1(PETSC_COMM_SELF,
314e83e87a5Sjeremylt           PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D",
315e83e87a5Sjeremylt           c);
3169b072555Sjeremylt     num_nodes = num_indices / field_off[num_fields];
3179b072555Sjeremylt     for (PetscInt i = 0; i < num_nodes; i++) {
318e83e87a5Sjeremylt       PetscInt ii = i;
319e83e87a5Sjeremylt       if (flip) {
3209b072555Sjeremylt         if (P == num_nodes) ii = num_nodes - 1 - i;
3219b072555Sjeremylt         else if (P*P == num_nodes) {
322e83e87a5Sjeremylt           PetscInt row = i / P, col = i % P;
323e83e87a5Sjeremylt           ii = row + col * P;
324e83e87a5Sjeremylt         } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP,
325e83e87a5Sjeremylt                           "No support for flipping point with %D nodes != P (%D) or P^2",
3269b072555Sjeremylt                           num_nodes, P);
327e83e87a5Sjeremylt       }
328e83e87a5Sjeremylt       // Check that indices are blocked by node and thus can be coalesced as a single field with
3299b072555Sjeremylt       // field_off[num_fields] = sum(num_comp) components.
3309b072555Sjeremylt       for (PetscInt f = 0; f < num_fields; f++) {
3319b072555Sjeremylt         for (PetscInt j = 0; j < num_comp[f]; j++) {
3329b072555Sjeremylt           if (Involute(indices[field_off[f]*num_nodes + ii*num_comp[f] + j])
3339b072555Sjeremylt               != Involute(indices[ii*num_comp[0]]) + field_off[f] + j)
334e83e87a5Sjeremylt             SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
335e83e87a5Sjeremylt                      "Cell %D closure indices not interlaced for node %D field %D component %D",
336e83e87a5Sjeremylt                      c, ii, f, j);
337e83e87a5Sjeremylt         }
338e83e87a5Sjeremylt       }
339e83e87a5Sjeremylt       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
3409b072555Sjeremylt       PetscInt loc = Involute(indices[ii*num_comp[0]]);
3419b072555Sjeremylt       elem_restr_offsets[e_offset++] = loc;
342e83e87a5Sjeremylt     }
343e83e87a5Sjeremylt     ierr = DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE,
3449b072555Sjeremylt                                        &num_indices, &indices, NULL, NULL);
345e83e87a5Sjeremylt     CHKERRQ(ierr);
346e83e87a5Sjeremylt   }
3479b072555Sjeremylt   if (e_offset != num_elem*PetscPowInt(P, topo_dim))
348e83e87a5Sjeremylt     SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB,
3499b072555Sjeremylt              "ElemRestriction of size (%D,%D) initialized %D nodes", num_elem,
3509b072555Sjeremylt              PetscPowInt(P, topo_dim),e_offset);
3519b072555Sjeremylt   if (iter_is) {
3529b072555Sjeremylt     ierr = ISRestoreIndices(iter_is, &iter_indices); CHKERRQ(ierr);
353e83e87a5Sjeremylt   }
3549b072555Sjeremylt   ierr = ISDestroy(&iter_is); CHKERRQ(ierr);
355e83e87a5Sjeremylt 
3569b072555Sjeremylt   ierr = DMGetLocalVector(dm, &U_loc); CHKERRQ(ierr);
3579b072555Sjeremylt   ierr = VecGetLocalSize(U_loc, &num_dof); CHKERRQ(ierr);
3589b072555Sjeremylt   ierr = DMRestoreLocalVector(dm, &U_loc); CHKERRQ(ierr);
3599b072555Sjeremylt   CeedElemRestrictionCreate(ceed, num_elem, PetscPowInt(P, topo_dim),
3609b072555Sjeremylt                             field_off[num_fields], 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES,
3619b072555Sjeremylt                             elem_restr_offsets, elem_restr);
3629b072555Sjeremylt   ierr = PetscFree(elem_restr_offsets); CHKERRQ(ierr);
363e83e87a5Sjeremylt   PetscFunctionReturn(0);
364e83e87a5Sjeremylt };
365e83e87a5Sjeremylt 
366e83e87a5Sjeremylt // -----------------------------------------------------------------------------
367