xref: /petsc/src/dm/impls/plex/plexsubmesh.c (revision 4a189a86914a48eaa7599a21ab4dbba0dd62c2c5)
1 #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 #include <petscsf.h>
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "DMPlexMarkBoundaryFaces"
6 /*@
7   DMPlexMarkBoundaryFaces - Mark all faces on the boundary
8 
9   Not Collective
10 
11   Input Parameter:
12 . dm - The original DM
13 
14   Output Parameter:
15 . label - The DMLabel marking boundary faces with value 1
16 
17   Level: developer
18 
19 .seealso: DMLabelCreate(), DMPlexCreateLabel()
20 @*/
21 PetscErrorCode DMPlexMarkBoundaryFaces(DM dm, DMLabel label)
22 {
23   PetscInt       fStart, fEnd, f;
24   PetscErrorCode ierr;
25 
26   PetscFunctionBegin;
27   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
28   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
29   for (f = fStart; f < fEnd; ++f) {
30     PetscInt supportSize;
31 
32     ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr);
33     if (supportSize == 1) {
34       ierr = DMLabelSetValue(label, f, 1);CHKERRQ(ierr);
35     }
36   }
37   PetscFunctionReturn(0);
38 }
39 
40 #undef __FUNCT__
41 #define __FUNCT__ "DMPlexLabelComplete"
42 /*@
43   DMPlexLabelComplete - Starting with a label marking points on a surface, we add the transitive closure to the surface
44 
45   Input Parameters:
46 + dm - The DM
47 - label - A DMLabel marking the surface points
48 
49   Output Parameter:
50 . label - A DMLabel marking all surface points in the transitive closure
51 
52   Level: developer
53 
54 .seealso: DMPlexLabelCohesiveComplete()
55 @*/
56 PetscErrorCode DMPlexLabelComplete(DM dm, DMLabel label)
57 {
58   IS              valueIS;
59   const PetscInt *values;
60   PetscInt        numValues, v;
61   PetscErrorCode  ierr;
62 
63   PetscFunctionBegin;
64   ierr = DMLabelGetNumValues(label, &numValues);CHKERRQ(ierr);
65   ierr = DMLabelGetValueIS(label, &valueIS);CHKERRQ(ierr);
66   ierr = ISGetIndices(valueIS, &values);CHKERRQ(ierr);
67   for (v = 0; v < numValues; ++v) {
68     IS              pointIS;
69     const PetscInt *points;
70     PetscInt        numPoints, p;
71 
72     ierr = DMLabelGetStratumSize(label, values[v], &numPoints);CHKERRQ(ierr);
73     ierr = DMLabelGetStratumIS(label, values[v], &pointIS);CHKERRQ(ierr);
74     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
75     for (p = 0; p < numPoints; ++p) {
76       PetscInt *closure = NULL;
77       PetscInt  closureSize, c;
78 
79       ierr = DMPlexGetTransitiveClosure(dm, points[p], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
80       for (c = 0; c < closureSize*2; c += 2) {
81         ierr = DMLabelSetValue(label, closure[c], values[v]);CHKERRQ(ierr);
82       }
83       ierr = DMPlexRestoreTransitiveClosure(dm, points[p], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
84     }
85     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
86     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
87   }
88   ierr = ISRestoreIndices(valueIS, &values);CHKERRQ(ierr);
89   ierr = ISDestroy(&valueIS);CHKERRQ(ierr);
90   PetscFunctionReturn(0);
91 }
92 
93 #undef __FUNCT__
94 #define __FUNCT__ "DMPlexShiftPoint_Internal"
95 PETSC_STATIC_INLINE PetscInt DMPlexShiftPoint_Internal(PetscInt p, PetscInt depth, PetscInt depthEnd[], PetscInt depthShift[])
96 {
97   if (depth < 0) return p;
98   /* Cells    */ if (p < depthEnd[depth])   return p;
99   /* Vertices */ if (p < depthEnd[0])       return p + depthShift[depth];
100   /* Faces    */ if (p < depthEnd[depth-1]) return p + depthShift[depth] + depthShift[0];
101   /* Edges    */                            return p + depthShift[depth] + depthShift[0] + depthShift[depth-1];
102 }
103 
104 #undef __FUNCT__
105 #define __FUNCT__ "DMPlexShiftSizes_Internal"
106 static PetscErrorCode DMPlexShiftSizes_Internal(DM dm, PetscInt depthShift[], DM dmNew)
107 {
108   PetscInt      *depthEnd;
109   PetscInt       depth = 0, d, pStart, pEnd, p;
110   PetscErrorCode ierr;
111 
112   PetscFunctionBegin;
113   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
114   if (depth < 0) PetscFunctionReturn(0);
115   ierr = PetscMalloc((depth+1) * sizeof(PetscInt), &depthEnd);CHKERRQ(ierr);
116   /* Step 1: Expand chart */
117   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
118   for (d = 0; d <= depth; ++d) {
119     pEnd += depthShift[d];
120     ierr  = DMPlexGetDepthStratum(dm, d, NULL, &depthEnd[d]);CHKERRQ(ierr);
121   }
122   ierr = DMPlexSetChart(dmNew, pStart, pEnd);CHKERRQ(ierr);
123   /* Step 2: Set cone and support sizes */
124   for (d = 0; d <= depth; ++d) {
125     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
126     for (p = pStart; p < pEnd; ++p) {
127       PetscInt newp = DMPlexShiftPoint_Internal(p, depth, depthEnd, depthShift);
128       PetscInt size;
129 
130       ierr = DMPlexGetConeSize(dm, p, &size);CHKERRQ(ierr);
131       ierr = DMPlexSetConeSize(dmNew, newp, size);CHKERRQ(ierr);
132       ierr = DMPlexGetSupportSize(dm, p, &size);CHKERRQ(ierr);
133       ierr = DMPlexSetSupportSize(dmNew, newp, size);CHKERRQ(ierr);
134     }
135   }
136   ierr = PetscFree(depthEnd);CHKERRQ(ierr);
137   PetscFunctionReturn(0);
138 }
139 
140 #undef __FUNCT__
141 #define __FUNCT__ "DMPlexShiftPoints_Internal"
142 static PetscErrorCode DMPlexShiftPoints_Internal(DM dm, PetscInt depthShift[], DM dmNew)
143 {
144   PetscInt      *depthEnd, *newpoints;
145   PetscInt       depth = 0, d, maxConeSize, maxSupportSize, pStart, pEnd, p;
146   PetscErrorCode ierr;
147 
148   PetscFunctionBegin;
149   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
150   if (depth < 0) PetscFunctionReturn(0);
151   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
152   ierr = PetscMalloc2(depth+1,PetscInt,&depthEnd,PetscMax(maxConeSize, maxSupportSize),PetscInt,&newpoints);CHKERRQ(ierr);
153   for (d = 0; d <= depth; ++d) {
154     ierr = DMPlexGetDepthStratum(dm, d, NULL, &depthEnd[d]);CHKERRQ(ierr);
155   }
156   /* Step 5: Set cones and supports */
157   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
158   for (p = pStart; p < pEnd; ++p) {
159     const PetscInt *points = NULL, *orientations = NULL;
160     PetscInt        size, i, newp = DMPlexShiftPoint_Internal(p, depth, depthEnd, depthShift);
161 
162     ierr = DMPlexGetConeSize(dm, p, &size);CHKERRQ(ierr);
163     ierr = DMPlexGetCone(dm, p, &points);CHKERRQ(ierr);
164     ierr = DMPlexGetConeOrientation(dm, p, &orientations);CHKERRQ(ierr);
165     for (i = 0; i < size; ++i) {
166       newpoints[i] = DMPlexShiftPoint_Internal(points[i], depth, depthEnd, depthShift);
167     }
168     ierr = DMPlexSetCone(dmNew, newp, newpoints);CHKERRQ(ierr);
169     ierr = DMPlexSetConeOrientation(dmNew, newp, orientations);CHKERRQ(ierr);
170     ierr = DMPlexGetSupportSize(dm, p, &size);CHKERRQ(ierr);
171     ierr = DMPlexGetSupport(dm, p, &points);CHKERRQ(ierr);
172     for (i = 0; i < size; ++i) {
173       newpoints[i] = DMPlexShiftPoint_Internal(points[i], depth, depthEnd, depthShift);
174     }
175     ierr = DMPlexSetSupport(dmNew, newp, newpoints);CHKERRQ(ierr);
176   }
177   ierr = PetscFree2(depthEnd,newpoints);CHKERRQ(ierr);
178   PetscFunctionReturn(0);
179 }
180 
181 #undef __FUNCT__
182 #define __FUNCT__ "DMPlexShiftCoordinates_Internal"
183 static PetscErrorCode DMPlexShiftCoordinates_Internal(DM dm, PetscInt depthShift[], DM dmNew)
184 {
185   PetscSection   coordSection, newCoordSection;
186   Vec            coordinates, newCoordinates;
187   PetscScalar   *coords, *newCoords;
188   PetscInt      *depthEnd, coordSize;
189   PetscInt       dim, depth = 0, d, vStart, vEnd, vStartNew, vEndNew, v;
190   PetscErrorCode ierr;
191 
192   PetscFunctionBegin;
193   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
194   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
195   ierr = PetscMalloc((depth+1) * sizeof(PetscInt), &depthEnd);CHKERRQ(ierr);
196   for (d = 0; d <= depth; ++d) {
197     ierr = DMPlexGetDepthStratum(dm, d, NULL, &depthEnd[d]);CHKERRQ(ierr);
198   }
199   /* Step 8: Convert coordinates */
200   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
201   ierr = DMPlexGetDepthStratum(dmNew, 0, &vStartNew, &vEndNew);CHKERRQ(ierr);
202   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
203   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &newCoordSection);CHKERRQ(ierr);
204   ierr = PetscSectionSetNumFields(newCoordSection, 1);CHKERRQ(ierr);
205   ierr = PetscSectionSetFieldComponents(newCoordSection, 0, dim);CHKERRQ(ierr);
206   ierr = PetscSectionSetChart(newCoordSection, vStartNew, vEndNew);CHKERRQ(ierr);
207   for (v = vStartNew; v < vEndNew; ++v) {
208     ierr = PetscSectionSetDof(newCoordSection, v, dim);CHKERRQ(ierr);
209     ierr = PetscSectionSetFieldDof(newCoordSection, v, 0, dim);CHKERRQ(ierr);
210   }
211   ierr = PetscSectionSetUp(newCoordSection);CHKERRQ(ierr);
212   ierr = DMPlexSetCoordinateSection(dmNew, newCoordSection);CHKERRQ(ierr);
213   ierr = PetscSectionGetStorageSize(newCoordSection, &coordSize);CHKERRQ(ierr);
214   ierr = VecCreate(PetscObjectComm((PetscObject)dm), &newCoordinates);CHKERRQ(ierr);
215   ierr = PetscObjectSetName((PetscObject) newCoordinates, "coordinates");CHKERRQ(ierr);
216   ierr = VecSetSizes(newCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
217   ierr = VecSetFromOptions(newCoordinates);CHKERRQ(ierr);
218   ierr = DMSetCoordinatesLocal(dmNew, newCoordinates);CHKERRQ(ierr);
219   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
220   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
221   ierr = VecGetArray(newCoordinates, &newCoords);CHKERRQ(ierr);
222   for (v = vStart; v < vEnd; ++v) {
223     PetscInt dof, off, noff, d;
224 
225     ierr = PetscSectionGetDof(coordSection, v, &dof);CHKERRQ(ierr);
226     ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
227     ierr = PetscSectionGetOffset(newCoordSection, DMPlexShiftPoint_Internal(v, depth, depthEnd, depthShift), &noff);CHKERRQ(ierr);
228     for (d = 0; d < dof; ++d) {
229       newCoords[noff+d] = coords[off+d];
230     }
231   }
232   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
233   ierr = VecRestoreArray(newCoordinates, &newCoords);CHKERRQ(ierr);
234   ierr = VecDestroy(&newCoordinates);CHKERRQ(ierr);
235   ierr = PetscSectionDestroy(&newCoordSection);CHKERRQ(ierr);
236   ierr = PetscFree(depthEnd);CHKERRQ(ierr);
237   PetscFunctionReturn(0);
238 }
239 
240 #undef __FUNCT__
241 #define __FUNCT__ "DMPlexShiftSF_Internal"
242 static PetscErrorCode DMPlexShiftSF_Internal(DM dm, PetscInt depthShift[], DM dmNew)
243 {
244   PetscInt          *depthEnd;
245   PetscInt           depth = 0, d;
246   PetscSF            sfPoint, sfPointNew;
247   const PetscSFNode *remotePoints;
248   PetscSFNode       *gremotePoints;
249   const PetscInt    *localPoints;
250   PetscInt          *glocalPoints, *newLocation, *newRemoteLocation;
251   PetscInt           numRoots, numLeaves, l, pStart, pEnd, totShift = 0;
252   PetscMPIInt        numProcs;
253   PetscErrorCode     ierr;
254 
255   PetscFunctionBegin;
256   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
257   ierr = PetscMalloc((depth+1) * sizeof(PetscInt), &depthEnd);CHKERRQ(ierr);
258   for (d = 0; d <= depth; ++d) {
259     totShift += depthShift[d];
260     ierr      = DMPlexGetDepthStratum(dm, d, NULL, &depthEnd[d]);CHKERRQ(ierr);
261   }
262   /* Step 9: Convert pointSF */
263   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr);
264   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
265   ierr = DMGetPointSF(dmNew, &sfPointNew);CHKERRQ(ierr);
266   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
267   ierr = PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr);
268   if (numRoots >= 0) {
269     ierr = PetscMalloc2(numRoots,PetscInt,&newLocation,pEnd-pStart,PetscInt,&newRemoteLocation);CHKERRQ(ierr);
270     for (l=0; l<numRoots; l++) newLocation[l] = DMPlexShiftPoint_Internal(l, depth, depthEnd, depthShift);
271     ierr = PetscSFBcastBegin(sfPoint, MPIU_INT, newLocation, newRemoteLocation);CHKERRQ(ierr);
272     ierr = PetscSFBcastEnd(sfPoint, MPIU_INT, newLocation, newRemoteLocation);CHKERRQ(ierr);
273     ierr = PetscMalloc(numLeaves * sizeof(PetscInt),    &glocalPoints);CHKERRQ(ierr);
274     ierr = PetscMalloc(numLeaves * sizeof(PetscSFNode), &gremotePoints);CHKERRQ(ierr);
275     for (l = 0; l < numLeaves; ++l) {
276       glocalPoints[l]        = DMPlexShiftPoint_Internal(localPoints[l], depth, depthEnd, depthShift);
277       gremotePoints[l].rank  = remotePoints[l].rank;
278       gremotePoints[l].index = newRemoteLocation[localPoints[l]];
279     }
280     ierr = PetscFree2(newLocation,newRemoteLocation);CHKERRQ(ierr);
281     ierr = PetscSFSetGraph(sfPointNew, numRoots + totShift, numLeaves, glocalPoints, PETSC_OWN_POINTER, gremotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
282   }
283   ierr = PetscFree(depthEnd);CHKERRQ(ierr);
284   PetscFunctionReturn(0);
285 }
286 
287 #undef __FUNCT__
288 #define __FUNCT__ "DMPlexShiftLabels_Internal"
289 static PetscErrorCode DMPlexShiftLabels_Internal(DM dm, PetscInt depthShift[], DM dmNew)
290 {
291   PetscSF            sfPoint;
292   DMLabel            vtkLabel, ghostLabel;
293   PetscInt          *depthEnd;
294   const PetscSFNode *leafRemote;
295   const PetscInt    *leafLocal;
296   PetscInt           depth = 0, d, numLeaves, numLabels, l, cStart, cEnd, c, fStart, fEnd, f;
297   PetscMPIInt        rank;
298   PetscErrorCode     ierr;
299 
300   PetscFunctionBegin;
301   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
302   ierr = PetscMalloc((depth+1) * sizeof(PetscInt), &depthEnd);CHKERRQ(ierr);
303   for (d = 0; d <= depth; ++d) {
304     ierr = DMPlexGetDepthStratum(dm, d, NULL, &depthEnd[d]);CHKERRQ(ierr);
305   }
306   /* Step 10: Convert labels */
307   ierr = DMPlexGetNumLabels(dm, &numLabels);CHKERRQ(ierr);
308   for (l = 0; l < numLabels; ++l) {
309     DMLabel         label, newlabel;
310     const char     *lname;
311     PetscBool       isDepth;
312     IS              valueIS;
313     const PetscInt *values;
314     PetscInt        numValues, val;
315 
316     ierr = DMPlexGetLabelName(dm, l, &lname);CHKERRQ(ierr);
317     ierr = PetscStrcmp(lname, "depth", &isDepth);CHKERRQ(ierr);
318     if (isDepth) continue;
319     ierr = DMPlexCreateLabel(dmNew, lname);CHKERRQ(ierr);
320     ierr = DMPlexGetLabel(dm, lname, &label);CHKERRQ(ierr);
321     ierr = DMPlexGetLabel(dmNew, lname, &newlabel);CHKERRQ(ierr);
322     ierr = DMLabelGetValueIS(label, &valueIS);CHKERRQ(ierr);
323     ierr = ISGetLocalSize(valueIS, &numValues);CHKERRQ(ierr);
324     ierr = ISGetIndices(valueIS, &values);CHKERRQ(ierr);
325     for (val = 0; val < numValues; ++val) {
326       IS              pointIS;
327       const PetscInt *points;
328       PetscInt        numPoints, p;
329 
330       ierr = DMLabelGetStratumIS(label, values[val], &pointIS);CHKERRQ(ierr);
331       ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr);
332       ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
333       for (p = 0; p < numPoints; ++p) {
334         const PetscInt newpoint = DMPlexShiftPoint_Internal(points[p], depth, depthEnd, depthShift);
335 
336         ierr = DMLabelSetValue(newlabel, newpoint, values[val]);CHKERRQ(ierr);
337       }
338       ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
339       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
340     }
341     ierr = ISRestoreIndices(valueIS, &values);CHKERRQ(ierr);
342     ierr = ISDestroy(&valueIS);CHKERRQ(ierr);
343   }
344   ierr = PetscFree(depthEnd);CHKERRQ(ierr);
345   /* Step 11: Make label for output (vtk) and to mark ghost points (ghost) */
346   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
347   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
348   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
349   ierr = PetscSFGetGraph(sfPoint, NULL, &numLeaves, &leafLocal, &leafRemote);CHKERRQ(ierr);
350   ierr = DMPlexCreateLabel(dmNew, "vtk");CHKERRQ(ierr);
351   ierr = DMPlexCreateLabel(dmNew, "ghost");CHKERRQ(ierr);
352   ierr = DMPlexGetLabel(dmNew, "vtk", &vtkLabel);CHKERRQ(ierr);
353   ierr = DMPlexGetLabel(dmNew, "ghost", &ghostLabel);CHKERRQ(ierr);
354   for (l = 0, c = cStart; l < numLeaves && c < cEnd; ++l, ++c) {
355     for (; c < leafLocal[l] && c < cEnd; ++c) {
356       ierr = DMLabelSetValue(vtkLabel, c, 1);CHKERRQ(ierr);
357     }
358     if (leafLocal[l] >= cEnd) break;
359     if (leafRemote[l].rank == rank) {
360       ierr = DMLabelSetValue(vtkLabel, c, 1);CHKERRQ(ierr);
361     } else {
362       ierr = DMLabelSetValue(ghostLabel, c, 2);CHKERRQ(ierr);
363     }
364   }
365   for (; c < cEnd; ++c) {
366     ierr = DMLabelSetValue(vtkLabel, c, 1);CHKERRQ(ierr);
367   }
368   if (0) {
369     ierr = PetscViewerASCIISynchronizedAllow(PETSC_VIEWER_STDOUT_WORLD, PETSC_TRUE);CHKERRQ(ierr);
370     ierr = DMLabelView(vtkLabel, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
371     ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
372   }
373   ierr = DMPlexGetHeightStratum(dmNew, 1, &fStart, &fEnd);CHKERRQ(ierr);
374   for (f = fStart; f < fEnd; ++f) {
375     PetscInt numCells;
376 
377     ierr = DMPlexGetSupportSize(dmNew, f, &numCells);CHKERRQ(ierr);
378     if (numCells < 2) {
379       ierr = DMLabelSetValue(ghostLabel, f, 1);CHKERRQ(ierr);
380     } else {
381       const PetscInt *cells = NULL;
382       PetscInt        vA, vB;
383 
384       ierr = DMPlexGetSupport(dmNew, f, &cells);CHKERRQ(ierr);
385       ierr = DMLabelGetValue(vtkLabel, cells[0], &vA);CHKERRQ(ierr);
386       ierr = DMLabelGetValue(vtkLabel, cells[1], &vB);CHKERRQ(ierr);
387       if (!vA && !vB) {ierr = DMLabelSetValue(ghostLabel, f, 1);CHKERRQ(ierr);}
388     }
389   }
390   if (0) {
391     ierr = PetscViewerASCIISynchronizedAllow(PETSC_VIEWER_STDOUT_WORLD, PETSC_TRUE);CHKERRQ(ierr);
392     ierr = DMLabelView(ghostLabel, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
393     ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
394   }
395   PetscFunctionReturn(0);
396 }
397 
398 #undef __FUNCT__
399 #define __FUNCT__ "DMPlexConstructGhostCells_Internal"
400 static PetscErrorCode DMPlexConstructGhostCells_Internal(DM dm, DMLabel label, PetscInt *numGhostCells, DM gdm)
401 {
402   IS              valueIS;
403   const PetscInt *values;
404   PetscInt       *depthShift;
405   PetscInt        depth = 0, numFS, fs, ghostCell, cEnd, c;
406   PetscErrorCode  ierr;
407 
408   PetscFunctionBegin;
409   /* Count ghost cells */
410   ierr = DMLabelGetValueIS(label, &valueIS);CHKERRQ(ierr);
411   ierr = ISGetLocalSize(valueIS, &numFS);CHKERRQ(ierr);
412   ierr = ISGetIndices(valueIS, &values);CHKERRQ(ierr);
413 
414   *numGhostCells = 0;
415   for (fs = 0; fs < numFS; ++fs) {
416     PetscInt numBdFaces;
417 
418     ierr = DMLabelGetStratumSize(label, values[fs], &numBdFaces);CHKERRQ(ierr);
419 
420     *numGhostCells += numBdFaces;
421   }
422   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
423   ierr = PetscMalloc((depth+1) * sizeof(PetscInt), &depthShift);CHKERRQ(ierr);
424   ierr = PetscMemzero(depthShift, (depth+1) * sizeof(PetscInt));CHKERRQ(ierr);
425   if (depth >= 0) depthShift[depth] = *numGhostCells;
426   ierr = DMPlexShiftSizes_Internal(dm, depthShift, gdm);CHKERRQ(ierr);
427   /* Step 3: Set cone/support sizes for new points */
428   ierr = DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);CHKERRQ(ierr);
429   for (c = cEnd; c < cEnd + *numGhostCells; ++c) {
430     ierr = DMPlexSetConeSize(gdm, c, 1);CHKERRQ(ierr);
431   }
432   for (fs = 0; fs < numFS; ++fs) {
433     IS              faceIS;
434     const PetscInt *faces;
435     PetscInt        numFaces, f;
436 
437     ierr = DMLabelGetStratumIS(label, values[fs], &faceIS);CHKERRQ(ierr);
438     ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
439     ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
440     for (f = 0; f < numFaces; ++f) {
441       PetscInt size;
442 
443       ierr = DMPlexGetSupportSize(dm, faces[f], &size);CHKERRQ(ierr);
444       if (size != 1) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM has boundary face %d with %d support cells", faces[f], size);
445       ierr = DMPlexSetSupportSize(gdm, faces[f] + *numGhostCells, 2);CHKERRQ(ierr);
446     }
447     ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
448     ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
449   }
450   /* Step 4: Setup ghosted DM */
451   ierr = DMSetUp(gdm);CHKERRQ(ierr);
452   ierr = DMPlexShiftPoints_Internal(dm, depthShift, gdm);CHKERRQ(ierr);
453   /* Step 6: Set cones and supports for new points */
454   ghostCell = cEnd;
455   for (fs = 0; fs < numFS; ++fs) {
456     IS              faceIS;
457     const PetscInt *faces;
458     PetscInt        numFaces, f;
459 
460     ierr = DMLabelGetStratumIS(label, values[fs], &faceIS);CHKERRQ(ierr);
461     ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
462     ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
463     for (f = 0; f < numFaces; ++f, ++ghostCell) {
464       PetscInt newFace = faces[f] + *numGhostCells;
465 
466       ierr = DMPlexSetCone(gdm, ghostCell, &newFace);CHKERRQ(ierr);
467       ierr = DMPlexInsertSupport(gdm, newFace, 1, ghostCell);CHKERRQ(ierr);
468     }
469     ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
470     ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
471   }
472   ierr = ISRestoreIndices(valueIS, &values);CHKERRQ(ierr);
473   ierr = ISDestroy(&valueIS);CHKERRQ(ierr);
474   /* Step 7: Stratify */
475   ierr = DMPlexStratify(gdm);CHKERRQ(ierr);
476   ierr = DMPlexShiftCoordinates_Internal(dm, depthShift, gdm);CHKERRQ(ierr);
477   ierr = DMPlexShiftSF_Internal(dm, depthShift, gdm);CHKERRQ(ierr);
478   ierr = DMPlexShiftLabels_Internal(dm, depthShift, gdm);CHKERRQ(ierr);
479   ierr = PetscFree(depthShift);CHKERRQ(ierr);
480   PetscFunctionReturn(0);
481 }
482 
483 #undef __FUNCT__
484 #define __FUNCT__ "DMPlexConstructGhostCells"
485 /*@C
486   DMPlexConstructGhostCells - Construct ghost cells which connect to every boundary face
487 
488   Collective on dm
489 
490   Input Parameters:
491 + dm - The original DM
492 - labelName - The label specifying the boundary faces, or "Face Sets" if this is NULL
493 
494   Output Parameters:
495 + numGhostCells - The number of ghost cells added to the DM
496 - dmGhosted - The new DM
497 
498   Note: If no label exists of that name, one will be created marking all boundary faces
499 
500   Level: developer
501 
502 .seealso: DMCreate()
503 */
504 PetscErrorCode DMPlexConstructGhostCells(DM dm, const char labelName[], PetscInt *numGhostCells, DM *dmGhosted)
505 {
506   DM             gdm;
507   DMLabel        label;
508   const char    *name = labelName ? labelName : "Face Sets";
509   PetscInt       dim;
510   PetscErrorCode ierr;
511 
512   PetscFunctionBegin;
513   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
514   PetscValidPointer(numGhostCells, 3);
515   PetscValidPointer(dmGhosted, 4);
516   ierr = DMCreate(PetscObjectComm((PetscObject)dm), &gdm);CHKERRQ(ierr);
517   ierr = DMSetType(gdm, DMPLEX);CHKERRQ(ierr);
518   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
519   ierr = DMPlexSetDimension(gdm, dim);CHKERRQ(ierr);
520   ierr = DMPlexGetLabel(dm, name, &label);CHKERRQ(ierr);
521   if (!label) {
522     /* Get label for boundary faces */
523     ierr = DMPlexCreateLabel(dm, name);CHKERRQ(ierr);
524     ierr = DMPlexGetLabel(dm, name, &label);CHKERRQ(ierr);
525     ierr = DMPlexMarkBoundaryFaces(dm, label);CHKERRQ(ierr);
526   }
527   ierr = DMPlexConstructGhostCells_Internal(dm, label, numGhostCells, gdm);CHKERRQ(ierr);
528   ierr = DMSetFromOptions(gdm);CHKERRQ(ierr);
529   *dmGhosted = gdm;
530   PetscFunctionReturn(0);
531 }
532 
533 #undef __FUNCT__
534 #define __FUNCT__ "DMPlexConstructCohesiveCells_Internal"
535 static PetscErrorCode DMPlexConstructCohesiveCells_Internal(DM dm, DMLabel label, DM sdm)
536 {
537   MPI_Comm        comm;
538   IS              valueIS, *pointIS;
539   const PetscInt *values, **splitPoints;
540   PetscSection    coordSection;
541   Vec             coordinates;
542   PetscScalar    *coords;
543   PetscInt       *depthShift, *depthOffset, *pMaxNew, *numSplitPoints, *coneNew, *coneONew, *supportNew;
544   PetscInt        shift = 100, depth = 0, dep, dim, d, numSP = 0, sp, maxConeSize, maxSupportSize, numLabels, vStart, vEnd, pEnd, p, v;
545   PetscErrorCode  ierr;
546 
547   PetscFunctionBegin;
548   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
549   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
550   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
551   /* Count split points and add cohesive cells */
552   if (label) {
553     ierr = DMLabelGetValueIS(label, &valueIS);CHKERRQ(ierr);
554     ierr = ISGetLocalSize(valueIS, &numSP);CHKERRQ(ierr);
555     ierr = ISGetIndices(valueIS, &values);CHKERRQ(ierr);
556   }
557   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
558   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
559   ierr = PetscMalloc6(depth+1,PetscInt,&depthShift,depth+1,PetscInt,&depthOffset,depth+1,PetscInt,&pMaxNew,maxConeSize*3,PetscInt,&coneNew,maxConeSize*3,PetscInt,&coneONew,maxSupportSize,PetscInt,&supportNew);CHKERRQ(ierr);
560   ierr = PetscMalloc3(depth+1,IS,&pointIS,depth+1,PetscInt,&numSplitPoints,depth+1,const PetscInt*,&splitPoints);CHKERRQ(ierr);
561   ierr = PetscMemzero(depthShift, (depth+1) * sizeof(PetscInt));CHKERRQ(ierr);
562   for (d = 0; d <= depth; ++d) {
563     ierr              = DMPlexGetDepthStratum(dm, d, NULL, &pMaxNew[d]);CHKERRQ(ierr);
564     numSplitPoints[d] = 0;
565     splitPoints[d]    = NULL;
566     pointIS[d]        = NULL;
567   }
568   for (sp = 0; sp < numSP; ++sp) {
569     const PetscInt dep = values[sp];
570 
571     if ((dep < 0) || (dep > depth)) continue;
572     ierr = DMLabelGetStratumSize(label, dep, &depthShift[dep]);CHKERRQ(ierr);
573     ierr = DMLabelGetStratumIS(label, dep, &pointIS[dep]);CHKERRQ(ierr);
574     if (pointIS[dep]) {
575       ierr = ISGetLocalSize(pointIS[dep], &numSplitPoints[dep]);CHKERRQ(ierr);
576       ierr = ISGetIndices(pointIS[dep], &splitPoints[dep]);CHKERRQ(ierr);
577     }
578   }
579   if (depth >= 0) {
580     /* Calculate number of additional points */
581     depthShift[depth] = depthShift[depth-1]; /* There is a cohesive cell for every split face   */
582     depthShift[1]    += depthShift[0];       /* There is a cohesive edge for every split vertex */
583     /* Calculate hybrid bound for each dimension */
584     pMaxNew[0] += depthShift[depth];
585     if (depth > 1) pMaxNew[dim-1] += depthShift[depth] + depthShift[0];
586     if (depth > 2) pMaxNew[1]     += depthShift[depth] + depthShift[0] + depthShift[dim-1];
587 
588     /* Calculate point offset for each dimension */
589     depthOffset[depth] = 0;
590     depthOffset[0]     = depthOffset[depth] + depthShift[depth];
591     if (depth > 1) depthOffset[dim-1] = depthOffset[0]     + depthShift[0];
592     if (depth > 2) depthOffset[1]     = depthOffset[dim-1] + depthShift[dim-1];
593   }
594   ierr = DMPlexShiftSizes_Internal(dm, depthShift, sdm);CHKERRQ(ierr);
595   /* Step 3: Set cone/support sizes for new points */
596   for (dep = 0; dep <= depth; ++dep) {
597     for (p = 0; p < numSplitPoints[dep]; ++p) {
598       const PetscInt  oldp   = splitPoints[dep][p];
599       const PetscInt  newp   = depthOffset[dep] + oldp;
600       const PetscInt  splitp = pMaxNew[dep] + p;
601       const PetscInt *support;
602       PetscInt        coneSize, supportSize, q, e;
603 
604       ierr = DMPlexGetConeSize(dm, oldp, &coneSize);CHKERRQ(ierr);
605       ierr = DMPlexSetConeSize(sdm, splitp, coneSize);CHKERRQ(ierr);
606       ierr = DMPlexGetSupportSize(dm, oldp, &supportSize);CHKERRQ(ierr);
607       ierr = DMPlexSetSupportSize(sdm, splitp, supportSize);CHKERRQ(ierr);
608       if (dep == depth-1) {
609         const PetscInt ccell = pMaxNew[depth] + p;
610         /* Add cohesive cells, they are prisms */
611         ierr = DMPlexSetConeSize(sdm, ccell, 2 + coneSize);CHKERRQ(ierr);
612       } else if (dep == 0) {
613         const PetscInt cedge = pMaxNew[1] + (depthShift[1] - depthShift[0]) + p;
614 
615         ierr = DMPlexGetSupport(dm, oldp, &support);CHKERRQ(ierr);
616         /* Split old vertex: Edges in old split faces and new cohesive edge */
617         for (e = 0, q = 0; e < supportSize; ++e) {
618           PetscInt val;
619 
620           ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr);
621           if ((val == 1) || (val == (shift + 1))) ++q;
622         }
623         ierr = DMPlexSetSupportSize(sdm, newp, q+1);CHKERRQ(ierr);
624         /* Split new vertex: Edges in new split faces and new cohesive edge */
625         for (e = 0, q = 0; e < supportSize; ++e) {
626           PetscInt val;
627 
628           ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr);
629           if ((val == 1) || (val == -(shift + 1))) ++q;
630         }
631         ierr = DMPlexSetSupportSize(sdm, splitp, q+1);CHKERRQ(ierr);
632         /* Add cohesive edges */
633         ierr = DMPlexSetConeSize(sdm, cedge, 2);CHKERRQ(ierr);
634         /* Punt for now on support, you loop over closure, extract faces, check which ones are in the label */
635       } else if (dep == dim-2) {
636         ierr = DMPlexGetSupport(dm, oldp, &support);CHKERRQ(ierr);
637         /* Split old edge: Faces in positive side cells and old split faces */
638         for (e = 0, q = 0; e < supportSize; ++e) {
639           PetscInt val;
640 
641           ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr);
642           if ((val == dim-1) || (val == (shift + dim-1))) ++q;
643         }
644         ierr = DMPlexSetSupportSize(sdm, newp, q);CHKERRQ(ierr);
645         /* Split new edge: Faces in negative side cells and new split faces */
646         for (e = 0, q = 0; e < supportSize; ++e) {
647           PetscInt val;
648 
649           ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr);
650           if ((val == dim-1) || (val == -(shift + dim-1))) ++q;
651         }
652         ierr = DMPlexSetSupportSize(sdm, splitp, q);CHKERRQ(ierr);
653       }
654     }
655   }
656   /* Step 4: Setup split DM */
657   ierr = DMSetUp(sdm);CHKERRQ(ierr);
658   ierr = DMPlexShiftPoints_Internal(dm, depthShift, sdm);CHKERRQ(ierr);
659   /* Step 6: Set cones and supports for new points */
660   for (dep = 0; dep <= depth; ++dep) {
661     for (p = 0; p < numSplitPoints[dep]; ++p) {
662       const PetscInt  oldp   = splitPoints[dep][p];
663       const PetscInt  newp   = depthOffset[dep] + oldp;
664       const PetscInt  splitp = pMaxNew[dep] + p;
665       const PetscInt *cone, *support, *ornt;
666       PetscInt        coneSize, supportSize, q, v, e, s;
667 
668       ierr = DMPlexGetConeSize(dm, oldp, &coneSize);CHKERRQ(ierr);
669       ierr = DMPlexGetCone(dm, oldp, &cone);CHKERRQ(ierr);
670       ierr = DMPlexGetConeOrientation(dm, oldp, &ornt);CHKERRQ(ierr);
671       ierr = DMPlexGetSupportSize(dm, oldp, &supportSize);CHKERRQ(ierr);
672       ierr = DMPlexGetSupport(dm, oldp, &support);CHKERRQ(ierr);
673       if (dep == depth-1) {
674         const PetscInt  ccell = pMaxNew[depth] + p;
675         const PetscInt *supportF;
676 
677         /* Split face:       copy in old face to new face to start */
678         ierr = DMPlexGetSupport(sdm, newp,  &supportF);CHKERRQ(ierr);
679         ierr = DMPlexSetSupport(sdm, splitp, supportF);CHKERRQ(ierr);
680         /* Split old face:   old vertices/edges in cone so no change */
681         /* Split new face:   new vertices/edges in cone */
682         for (q = 0; q < coneSize; ++q) {
683           ierr = PetscFindInt(cone[q], numSplitPoints[dim-2], splitPoints[dim-2], &v);CHKERRQ(ierr);
684 
685           coneNew[2+q] = pMaxNew[dim-2] + v;
686         }
687         ierr = DMPlexSetCone(sdm, splitp, &coneNew[2]);CHKERRQ(ierr);
688         ierr = DMPlexSetConeOrientation(sdm, splitp, ornt);CHKERRQ(ierr);
689         /* Face support */
690         for (s = 0; s < supportSize; ++s) {
691           PetscInt val;
692 
693           ierr = DMLabelGetValue(label, support[s], &val);CHKERRQ(ierr);
694           if (val < 0) {
695             /* Split old face:   Replace negative side cell with cohesive cell */
696              ierr = DMPlexInsertSupport(sdm, newp, s, ccell);CHKERRQ(ierr);
697           } else {
698             /* Split new face:   Replace positive side cell with cohesive cell */
699             ierr = DMPlexInsertSupport(sdm, splitp, s, ccell);CHKERRQ(ierr);
700             /* Get orientation for cohesive face */
701             {
702               const PetscInt *ncone, *nconeO;
703               PetscInt        nconeSize, nc;
704 
705               ierr = DMPlexGetConeSize(dm, support[s], &nconeSize);CHKERRQ(ierr);
706               ierr = DMPlexGetCone(dm, support[s], &ncone);CHKERRQ(ierr);
707               ierr = DMPlexGetConeOrientation(dm, support[s], &nconeO);CHKERRQ(ierr);
708               for (nc = 0; nc < nconeSize; ++nc) {
709                 if (ncone[nc] == oldp) {
710                   coneONew[0] = nconeO[nc];
711                   break;
712                 }
713               }
714               if (nc >= nconeSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not locate face %d in neighboring cell %d", oldp, support[s]);
715             }
716           }
717         }
718         /* Cohesive cell:    Old and new split face, then new cohesive edges */
719         coneNew[0] = newp;   /* Extracted negative side orientation above */
720         coneNew[1] = splitp; coneONew[1] = coneONew[0];
721         if (dim > 2) {
722           PetscInt *closure = NULL, closureSize, cl;
723 
724           ierr = DMPlexGetTransitiveClosure(dm, oldp, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
725           for (cl = 0, q = 0; cl < closureSize*2; cl += 2) {
726             const PetscInt clp = closure[cl];
727 
728             if ((clp >= vStart) && (clp < vEnd)) {
729               ierr = PetscFindInt(clp, numSplitPoints[0], splitPoints[0], &v);CHKERRQ(ierr);
730               coneNew[2+q]  = pMaxNew[1] + (depthShift[1] - depthShift[0]) + v;
731               coneONew[2+q] = 0;
732               ++q;
733             }
734           }
735           ierr = DMPlexRestoreTransitiveClosure(dm, oldp, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
736           if (q != coneSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of split face vertices %d should be %d", q, coneSize);
737         } else {
738           for (q = 0; q < coneSize; ++q) {
739             coneNew[2+q]  = (pMaxNew[1] - pMaxNew[dim-2]) + (depthShift[1] - depthShift[0]) + coneNew[2+q];
740             coneONew[2+q] = 0;
741           }
742         }
743         ierr = DMPlexSetCone(sdm, ccell, coneNew);CHKERRQ(ierr);
744         ierr = DMPlexSetConeOrientation(sdm, ccell, coneONew);CHKERRQ(ierr);
745       } else if (dep == 0) {
746         const PetscInt cedge = pMaxNew[1] + (depthShift[1] - depthShift[0]) + p;
747 
748         /* Split old vertex: Edges in old split faces and new cohesive edge */
749         for (e = 0, q = 0; e < supportSize; ++e) {
750           PetscInt val;
751 
752           ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr);
753           if ((val == 1) || (val == (shift + 1))) {
754             supportNew[q++] = depthOffset[1] + support[e];
755           }
756         }
757         supportNew[q] = cedge;
758 
759         ierr = DMPlexSetSupport(sdm, newp, supportNew);CHKERRQ(ierr);
760         /* Split new vertex: Edges in new split faces and new cohesive edge */
761         for (e = 0, q = 0; e < supportSize; ++e) {
762           PetscInt val, edge;
763 
764           ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr);
765           if (val == 1) {
766             ierr = PetscFindInt(support[e], numSplitPoints[1], splitPoints[1], &edge);CHKERRQ(ierr);
767             if (edge < 0) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Edge %d is not a split edge", support[e]);
768             supportNew[q++] = pMaxNew[1] + edge;
769           } else if (val == -(shift + 1)) {
770             supportNew[q++] = depthOffset[1] + support[e];
771           }
772         }
773         supportNew[q] = cedge;
774         ierr          = DMPlexSetSupport(sdm, splitp, supportNew);CHKERRQ(ierr);
775         /* Cohesive edge:    Old and new split vertex, punting on support */
776         coneNew[0] = newp;
777         coneNew[1] = splitp;
778         ierr       = DMPlexSetCone(sdm, cedge, coneNew);CHKERRQ(ierr);
779       } else if (dep == dim-2) {
780         /* Split old edge:   old vertices in cone so no change */
781         /* Split new edge:   new vertices in cone */
782         for (q = 0; q < coneSize; ++q) {
783           ierr = PetscFindInt(cone[q], numSplitPoints[dim-3], splitPoints[dim-3], &v);CHKERRQ(ierr);
784 
785           coneNew[q] = pMaxNew[dim-3] + v;
786         }
787         ierr = DMPlexSetCone(sdm, splitp, coneNew);CHKERRQ(ierr);
788         /* Split old edge: Faces in positive side cells and old split faces */
789         for (e = 0, q = 0; e < supportSize; ++e) {
790           PetscInt val;
791 
792           ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr);
793           if ((val == dim-1) || (val == (shift + dim-1))) {
794             supportNew[q++] = depthOffset[dim-1] + support[e];
795           }
796         }
797         ierr = DMPlexSetSupport(sdm, newp, supportNew);CHKERRQ(ierr);
798         /* Split new edge: Faces in negative side cells and new split faces */
799         for (e = 0, q = 0; e < supportSize; ++e) {
800           PetscInt val, face;
801 
802           ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr);
803           if (val == dim-1) {
804             ierr = PetscFindInt(support[e], numSplitPoints[dim-1], splitPoints[dim-1], &face);CHKERRQ(ierr);
805             if (face < 0) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Face %d is not a split face", support[e]);
806             supportNew[q++] = pMaxNew[dim-1] + face;
807           } else if (val == -(shift + dim-1)) {
808             supportNew[q++] = depthOffset[dim-1] + support[e];
809           }
810         }
811         ierr = DMPlexSetSupport(sdm, splitp, supportNew);CHKERRQ(ierr);
812       }
813     }
814   }
815   /* Step 6b: Replace split points in negative side cones */
816   for (sp = 0; sp < numSP; ++sp) {
817     PetscInt        dep = values[sp];
818     IS              pIS;
819     PetscInt        numPoints;
820     const PetscInt *points;
821 
822     if (dep >= 0) continue;
823     ierr = DMLabelGetStratumIS(label, dep, &pIS);CHKERRQ(ierr);
824     if (!pIS) continue;
825     dep  = -dep - shift;
826     ierr = ISGetLocalSize(pIS, &numPoints);CHKERRQ(ierr);
827     ierr = ISGetIndices(pIS, &points);CHKERRQ(ierr);
828     for (p = 0; p < numPoints; ++p) {
829       const PetscInt  oldp = points[p];
830       const PetscInt  newp = depthOffset[dep] + oldp;
831       const PetscInt *cone;
832       PetscInt        coneSize, c;
833       PetscBool       replaced = PETSC_FALSE;
834 
835       /* Negative edge: replace split vertex */
836       /* Negative cell: replace split face */
837       ierr = DMPlexGetConeSize(sdm, newp, &coneSize);CHKERRQ(ierr);
838       ierr = DMPlexGetCone(sdm, newp, &cone);CHKERRQ(ierr);
839       for (c = 0; c < coneSize; ++c) {
840         const PetscInt coldp = cone[c] - depthOffset[dep-1];
841         PetscInt       csplitp, cp, val;
842 
843         ierr = DMLabelGetValue(label, coldp, &val);CHKERRQ(ierr);
844         if (val == dep-1) {
845           ierr = PetscFindInt(coldp, numSplitPoints[dep-1], splitPoints[dep-1], &cp);CHKERRQ(ierr);
846           if (cp < 0) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "Point %d is not a split point of dimension %d", oldp, dep-1);
847           csplitp  = pMaxNew[dep-1] + cp;
848           ierr     = DMPlexInsertCone(sdm, newp, c, csplitp);CHKERRQ(ierr);
849           replaced = PETSC_TRUE;
850         }
851       }
852       /* Cells with only a vertex or edge on the submesh have no replacement */
853       /* if (!replaced) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "The cone of point %d does not contain split points", oldp); */
854     }
855     ierr = ISRestoreIndices(pIS, &points);CHKERRQ(ierr);
856     ierr = ISDestroy(&pIS);CHKERRQ(ierr);
857   }
858   /* Step 7: Stratify */
859   ierr = DMPlexStratify(sdm);CHKERRQ(ierr);
860   /* Step 8: Coordinates */
861   ierr = DMPlexShiftCoordinates_Internal(dm, depthShift, sdm);CHKERRQ(ierr);
862   ierr = DMPlexGetCoordinateSection(sdm, &coordSection);CHKERRQ(ierr);
863   ierr = DMGetCoordinatesLocal(sdm, &coordinates);CHKERRQ(ierr);
864   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
865   for (v = 0; v < (numSplitPoints ? numSplitPoints[0] : 0); ++v) {
866     const PetscInt newp   = depthOffset[0] + splitPoints[0][v];
867     const PetscInt splitp = pMaxNew[0] + v;
868     PetscInt       dof, off, soff, d;
869 
870     ierr = PetscSectionGetDof(coordSection, newp, &dof);CHKERRQ(ierr);
871     ierr = PetscSectionGetOffset(coordSection, newp, &off);CHKERRQ(ierr);
872     ierr = PetscSectionGetOffset(coordSection, splitp, &soff);CHKERRQ(ierr);
873     for (d = 0; d < dof; ++d) coords[soff+d] = coords[off+d];
874   }
875   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
876   /* Step 9: SF, if I can figure this out we can split the mesh in parallel */
877   ierr = DMPlexShiftSF_Internal(dm, depthShift, sdm);CHKERRQ(ierr);
878   /* Step 10: Labels */
879   ierr = DMPlexShiftLabels_Internal(dm, depthShift, sdm);CHKERRQ(ierr);
880   ierr = DMPlexGetNumLabels(sdm, &numLabels);CHKERRQ(ierr);
881   for (dep = 0; dep <= depth; ++dep) {
882     for (p = 0; p < numSplitPoints[dep]; ++p) {
883       const PetscInt newp   = depthOffset[dep] + splitPoints[dep][p];
884       const PetscInt splitp = pMaxNew[dep] + p;
885       PetscInt       l;
886 
887       for (l = 0; l < numLabels; ++l) {
888         DMLabel     mlabel;
889         const char *lname;
890         PetscInt    val;
891 
892         ierr = DMPlexGetLabelName(sdm, l, &lname);CHKERRQ(ierr);
893         ierr = DMPlexGetLabel(sdm, lname, &mlabel);CHKERRQ(ierr);
894         ierr = DMLabelGetValue(mlabel, newp, &val);CHKERRQ(ierr);
895         if (val >= 0) {
896           ierr = DMLabelSetValue(mlabel, splitp, val);CHKERRQ(ierr);
897           if (dep == 0) {
898             const PetscInt cedge = pMaxNew[1] + (depthShift[1] - depthShift[0]) + p;
899             ierr = DMLabelSetValue(mlabel, cedge, val);CHKERRQ(ierr);
900           }
901         }
902       }
903     }
904   }
905   for (sp = 0; sp < numSP; ++sp) {
906     const PetscInt dep = values[sp];
907 
908     if ((dep < 0) || (dep > depth)) continue;
909     if (pointIS[dep]) {ierr = ISRestoreIndices(pointIS[dep], &splitPoints[dep]);CHKERRQ(ierr);}
910     ierr = ISDestroy(&pointIS[dep]);CHKERRQ(ierr);
911   }
912   if (label) {
913     ierr = ISRestoreIndices(valueIS, &values);CHKERRQ(ierr);
914     ierr = ISDestroy(&valueIS);CHKERRQ(ierr);
915   }
916   ierr = DMPlexGetChart(sdm, NULL, &pEnd);CHKERRQ(ierr);
917   if (depth > 0) pMaxNew[0] += depthShift[0];        /* Account for shadow vertices */
918   if (depth > 1) pMaxNew[1]  = pEnd - depthShift[0]; /* There is a hybrid edge for every shadow vertex */
919   if (depth > 2) pMaxNew[2]  = -1;                   /* There are no hybrid faces */
920   ierr = DMPlexSetHybridBounds(sdm, depth >= 0 ? pMaxNew[depth] : PETSC_DETERMINE, depth>1 ? pMaxNew[depth-1] : PETSC_DETERMINE, depth>2 ? pMaxNew[1] : PETSC_DETERMINE, depth >= 0 ? pMaxNew[0] : PETSC_DETERMINE);CHKERRQ(ierr);
921   ierr = PetscFree6(depthShift, depthOffset, pMaxNew, coneNew, coneONew, supportNew);CHKERRQ(ierr);
922   ierr = PetscFree3(pointIS, numSplitPoints, splitPoints);CHKERRQ(ierr);
923   PetscFunctionReturn(0);
924 }
925 
926 #undef __FUNCT__
927 #define __FUNCT__ "DMPlexConstructCohesiveCells"
928 /*@C
929   DMPlexConstructCohesiveCells - Construct cohesive cells which split the face along an internal interface
930 
931   Collective on dm
932 
933   Input Parameters:
934 + dm - The original DM
935 - labelName - The label specifying the boundary faces (this could be auto-generated)
936 
937   Output Parameters:
938 - dmSplit - The new DM
939 
940   Level: developer
941 
942 .seealso: DMCreate(), DMPlexLabelCohesiveComplete()
943 @*/
944 PetscErrorCode DMPlexConstructCohesiveCells(DM dm, DMLabel label, DM *dmSplit)
945 {
946   DM             sdm;
947   PetscInt       dim;
948   PetscErrorCode ierr;
949 
950   PetscFunctionBegin;
951   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
952   PetscValidPointer(dmSplit, 4);
953   ierr = DMCreate(PetscObjectComm((PetscObject)dm), &sdm);CHKERRQ(ierr);
954   ierr = DMSetType(sdm, DMPLEX);CHKERRQ(ierr);
955   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
956   ierr = DMPlexSetDimension(sdm, dim);CHKERRQ(ierr);
957   switch (dim) {
958   case 2:
959   case 3:
960     ierr = DMPlexConstructCohesiveCells_Internal(dm, label, sdm);CHKERRQ(ierr);
961     break;
962   default:
963     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct cohesive cells for dimension %d", dim);
964   }
965   *dmSplit = sdm;
966   PetscFunctionReturn(0);
967 }
968 
969 #undef __FUNCT__
970 #define __FUNCT__ "DMPlexLabelCohesiveComplete"
971 /*@
972   DMPlexLabelCohesiveComplete - Starting with a label marking vertices on an internal surface, we add all other mesh pieces
973   to complete the surface
974 
975   Input Parameters:
976 + dm - The DM
977 . label - A DMLabel marking the surface vertices
978 . flip  - Flag to flip the submesh normal and replace points on the other side
979 - subdm - The subDM associated with the label, or NULL
980 
981   Output Parameter:
982 . label - A DMLabel marking all surface points
983 
984   Level: developer
985 
986 .seealso: DMPlexConstructCohesiveCells(), DMPlexLabelComplete()
987 @*/
988 PetscErrorCode DMPlexLabelCohesiveComplete(DM dm, DMLabel label, PetscBool flip, DM subdm)
989 {
990   IS              dimIS, subpointIS;
991   const PetscInt *points, *subpoints;
992   const PetscInt  rev   = flip ? -1 : 1;
993   PetscInt        shift = 100, dim, dep, cStart, cEnd, numPoints, numSubpoints, p, val;
994   PetscErrorCode  ierr;
995 
996   PetscFunctionBegin;
997   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
998   if (subdm) {
999     ierr = DMPlexCreateSubpointIS(subdm, &subpointIS);CHKERRQ(ierr);
1000     if (subpointIS) {
1001       ierr = ISGetLocalSize(subpointIS, &numSubpoints);CHKERRQ(ierr);
1002       ierr = ISGetIndices(subpointIS, &subpoints);CHKERRQ(ierr);
1003     }
1004   }
1005   /* Cell orientation for face gives the side of the fault */
1006   ierr = DMLabelGetStratumIS(label, dim-1, &dimIS);CHKERRQ(ierr);
1007   if (!dimIS) PetscFunctionReturn(0);
1008   ierr = ISGetLocalSize(dimIS, &numPoints);CHKERRQ(ierr);
1009   ierr = ISGetIndices(dimIS, &points);CHKERRQ(ierr);
1010   for (p = 0; p < numPoints; ++p) {
1011     const PetscInt *support;
1012     PetscInt        supportSize, s;
1013 
1014     ierr = DMPlexGetSupportSize(dm, points[p], &supportSize);CHKERRQ(ierr);
1015     if (supportSize != 2) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Split face %d has %d != 2 supports", points[p], supportSize);
1016     ierr = DMPlexGetSupport(dm, points[p], &support);CHKERRQ(ierr);
1017     for (s = 0; s < supportSize; ++s) {
1018       const PetscInt *cone, *ornt;
1019       PetscInt        coneSize, c;
1020       PetscBool       pos = PETSC_TRUE;
1021 
1022       ierr = DMPlexGetConeSize(dm, support[s], &coneSize);CHKERRQ(ierr);
1023       ierr = DMPlexGetCone(dm, support[s], &cone);CHKERRQ(ierr);
1024       ierr = DMPlexGetConeOrientation(dm, support[s], &ornt);CHKERRQ(ierr);
1025       for (c = 0; c < coneSize; ++c) {
1026         if (cone[c] == points[p]) {
1027           PetscInt o = ornt[c];
1028 
1029           if (subdm) {
1030             const PetscInt *subcone, *subornt;
1031             PetscInt        subpoint, subface, subconeSize, sc;
1032 
1033             ierr = PetscFindInt(support[s], numSubpoints, subpoints, &subpoint);CHKERRQ(ierr);
1034             ierr = PetscFindInt(points[p],  numSubpoints, subpoints, &subface);CHKERRQ(ierr);
1035             ierr = DMPlexGetConeSize(subdm, subpoint, &subconeSize);CHKERRQ(ierr);
1036             ierr = DMPlexGetCone(subdm, subpoint, &subcone);CHKERRQ(ierr);
1037             ierr = DMPlexGetConeOrientation(subdm, subpoint, &subornt);CHKERRQ(ierr);
1038             for (sc = 0; sc < subconeSize; ++sc) {
1039               if (subcone[sc] == subface) {
1040                 o = subornt[0];
1041                 break;
1042               }
1043             }
1044             if (sc >= subconeSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find point %d in cone for subpoint %d", points[p], subpoint);
1045           }
1046           if (o >= 0) {
1047             ierr = DMLabelSetValue(label, support[s],  rev*(shift+dim));CHKERRQ(ierr);
1048           } else {
1049             ierr = DMLabelSetValue(label, support[s], -rev*(shift+dim));CHKERRQ(ierr);
1050             pos  = PETSC_FALSE;
1051           }
1052           break;
1053         }
1054       }
1055       if (c == coneSize) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Cell split face %d support does not have it in the cone", points[p]);
1056       /* Put faces touching the fault in the label */
1057       for (c = 0; c < coneSize; ++c) {
1058         const PetscInt point = cone[c];
1059 
1060         ierr = DMLabelGetValue(label, point, &val);CHKERRQ(ierr);
1061         if (val == -1) {
1062           PetscInt *closure = NULL;
1063           PetscInt  closureSize, cl;
1064 
1065           ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1066           for (cl = 0; cl < closureSize*2; cl += 2) {
1067             const PetscInt clp = closure[cl];
1068 
1069             ierr = DMLabelGetValue(label, clp, &val);CHKERRQ(ierr);
1070             if ((val >= 0) && (val < dim-1)) {
1071               ierr = DMLabelSetValue(label, point, pos == PETSC_TRUE ? shift+dim-1 : -(shift+dim-1));CHKERRQ(ierr);
1072               break;
1073             }
1074           }
1075           ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1076         }
1077       }
1078     }
1079   }
1080   if (subdm) {
1081     if (subpointIS) {ierr = ISRestoreIndices(subpointIS, &subpoints);CHKERRQ(ierr);}
1082     ierr = ISDestroy(&subpointIS);CHKERRQ(ierr);
1083   }
1084   ierr = ISRestoreIndices(dimIS, &points);CHKERRQ(ierr);
1085   ierr = ISDestroy(&dimIS);CHKERRQ(ierr);
1086   /* Search for other cells/faces/edges connected to the fault by a vertex */
1087   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1088   ierr = DMLabelGetStratumIS(label, 0, &dimIS);CHKERRQ(ierr);
1089   if (!dimIS) PetscFunctionReturn(0);
1090   ierr = ISGetLocalSize(dimIS, &numPoints);CHKERRQ(ierr);
1091   ierr = ISGetIndices(dimIS, &points);CHKERRQ(ierr);
1092   for (p = 0; p < numPoints; ++p) {
1093     PetscInt *star = NULL;
1094     PetscInt  starSize, s;
1095     PetscInt  again = 1;  /* 0: Finished 1: Keep iterating after a change 2: No change */
1096 
1097     /* First mark cells connected to the fault */
1098     ierr = DMPlexGetTransitiveClosure(dm, points[p], PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
1099     while (again) {
1100       if (again > 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Could not classify all cells connected to the fault");
1101       again = 0;
1102       for (s = 0; s < starSize*2; s += 2) {
1103         const PetscInt  point = star[s];
1104         const PetscInt *cone;
1105         PetscInt        coneSize, c;
1106 
1107         if ((point < cStart) || (point >= cEnd)) continue;
1108         ierr = DMLabelGetValue(label, point, &val);CHKERRQ(ierr);
1109         if (val != -1) continue;
1110         again = 2;
1111         ierr  = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
1112         ierr  = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
1113         for (c = 0; c < coneSize; ++c) {
1114           ierr = DMLabelGetValue(label, cone[c], &val);CHKERRQ(ierr);
1115           if (val != -1) {
1116             if (abs(val) < shift) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Face %d on cell %d has an invalid label %d", cone[c], point, val);
1117             if (val > 0) {
1118               ierr = DMLabelSetValue(label, point,   shift+dim);CHKERRQ(ierr);
1119             } else {
1120               ierr = DMLabelSetValue(label, point, -(shift+dim));CHKERRQ(ierr);
1121             }
1122             again = 1;
1123             break;
1124           }
1125         }
1126       }
1127     }
1128     /* Classify the rest by cell membership */
1129     for (s = 0; s < starSize*2; s += 2) {
1130       const PetscInt point = star[s];
1131 
1132       ierr = DMLabelGetValue(label, point, &val);CHKERRQ(ierr);
1133       if (val == -1) {
1134         PetscInt  *sstar = NULL;
1135         PetscInt   sstarSize, ss;
1136         PetscBool  marked = PETSC_FALSE;
1137 
1138         ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &sstarSize, &sstar);CHKERRQ(ierr);
1139         for (ss = 0; ss < sstarSize*2; ss += 2) {
1140           const PetscInt spoint = sstar[ss];
1141 
1142           if ((spoint < cStart) || (spoint >= cEnd)) continue;
1143           ierr = DMLabelGetValue(label, spoint, &val);CHKERRQ(ierr);
1144           if (val == -1) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Cell %d in star of %d does not have a valid label", spoint, point);
1145           ierr = DMPlexGetLabelValue(dm, "depth", point, &dep);CHKERRQ(ierr);
1146           if (val > 0) {
1147             ierr = DMLabelSetValue(label, point,   shift+dep);CHKERRQ(ierr);
1148           } else {
1149             ierr = DMLabelSetValue(label, point, -(shift+dep));CHKERRQ(ierr);
1150           }
1151           marked = PETSC_TRUE;
1152           break;
1153         }
1154         ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, &sstarSize, &sstar);CHKERRQ(ierr);
1155         if (!marked) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Point %d could not be classified", point);
1156       }
1157     }
1158     ierr = DMPlexRestoreTransitiveClosure(dm, points[p], PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
1159   }
1160   ierr = ISRestoreIndices(dimIS, &points);CHKERRQ(ierr);
1161   ierr = ISDestroy(&dimIS);CHKERRQ(ierr);
1162   PetscFunctionReturn(0);
1163 }
1164 
1165 #undef __FUNCT__
1166 #define __FUNCT__ "DMPlexMarkSubmesh_Uninterpolated"
1167 /* Here we need the explicit assumption that:
1168 
1169      For any marked cell, the marked vertices constitute a single face
1170 */
1171 static PetscErrorCode DMPlexMarkSubmesh_Uninterpolated(DM dm, DMLabel vertexLabel, PetscInt value, DMLabel subpointMap, PetscInt *numFaces, PetscInt *nFV, DM subdm)
1172 {
1173   IS               subvertexIS;
1174   const PetscInt  *subvertices;
1175   PetscInt        *pStart, *pEnd, *pMax, pSize;
1176   PetscInt         depth, dim, d, numSubVerticesInitial = 0, v;
1177   PetscErrorCode   ierr;
1178 
1179   PetscFunctionBegin;
1180   *numFaces = 0;
1181   *nFV      = 0;
1182   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1183   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1184   pSize = PetscMax(depth, dim) + 1;
1185   ierr = PetscMalloc3(pSize,PetscInt,&pStart,pSize,PetscInt,&pEnd,pSize,PetscInt,&pMax);CHKERRQ(ierr);
1186   ierr = DMPlexGetHybridBounds(dm, &pMax[depth], depth>1 ? &pMax[depth-1] : NULL, depth > 2 ? &pMax[1] : NULL, &pMax[0]);CHKERRQ(ierr);
1187   for (d = 0; d <= depth; ++d) {
1188     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr);
1189     if (pMax[d] >= 0) pEnd[d] = PetscMin(pEnd[d], pMax[d]);
1190   }
1191   /* Loop over initial vertices and mark all faces in the collective star() */
1192   ierr = DMLabelGetStratumIS(vertexLabel, value, &subvertexIS);CHKERRQ(ierr);
1193   if (subvertexIS) {
1194     ierr = ISGetSize(subvertexIS, &numSubVerticesInitial);CHKERRQ(ierr);
1195     ierr = ISGetIndices(subvertexIS, &subvertices);CHKERRQ(ierr);
1196   }
1197   for (v = 0; v < numSubVerticesInitial; ++v) {
1198     const PetscInt vertex = subvertices[v];
1199     PetscInt      *star   = NULL;
1200     PetscInt       starSize, s, numCells = 0, c;
1201 
1202     ierr = DMPlexGetTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
1203     for (s = 0; s < starSize*2; s += 2) {
1204       const PetscInt point = star[s];
1205       if ((point >= pStart[depth]) && (point < pEnd[depth])) star[numCells++] = point;
1206     }
1207     for (c = 0; c < numCells; ++c) {
1208       const PetscInt cell    = star[c];
1209       PetscInt      *closure = NULL;
1210       PetscInt       closureSize, cl;
1211       PetscInt       cellLoc, numCorners = 0, faceSize = 0;
1212 
1213       ierr = DMLabelGetValue(subpointMap, cell, &cellLoc);CHKERRQ(ierr);
1214       if (cellLoc == 2) continue;
1215       if (cellLoc >= 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Cell %d has dimension %d in the surface label", cell, cellLoc);
1216       ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1217       for (cl = 0; cl < closureSize*2; cl += 2) {
1218         const PetscInt point = closure[cl];
1219         PetscInt       vertexLoc;
1220 
1221         if ((point >= pStart[0]) && (point < pEnd[0])) {
1222           ++numCorners;
1223           ierr = DMLabelGetValue(vertexLabel, point, &vertexLoc);CHKERRQ(ierr);
1224           if (vertexLoc == value) closure[faceSize++] = point;
1225         }
1226       }
1227       if (!(*nFV)) {ierr = DMPlexGetNumFaceVertices(dm, dim, numCorners, nFV);CHKERRQ(ierr);}
1228       if (faceSize > *nFV) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid submesh: Too many vertices %d of an element on the surface", faceSize);
1229       if (faceSize == *nFV) {
1230         const PetscInt *cells = NULL;
1231         PetscInt        numCells, nc;
1232 
1233         ++(*numFaces);
1234         for (cl = 0; cl < faceSize; ++cl) {
1235           ierr = DMLabelSetValue(subpointMap, closure[cl], 0);CHKERRQ(ierr);
1236         }
1237         ierr = DMPlexGetJoin(dm, faceSize, closure, &numCells, &cells);CHKERRQ(ierr);
1238         for (nc = 0; nc < numCells; ++nc) {
1239           ierr = DMLabelSetValue(subpointMap, cells[nc], 2);CHKERRQ(ierr);
1240         }
1241         ierr = DMPlexRestoreJoin(dm, faceSize, closure, &numCells, &cells);CHKERRQ(ierr);
1242       }
1243       ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1244     }
1245     ierr = DMPlexRestoreTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
1246   }
1247   if (subvertexIS) {
1248     ierr = ISRestoreIndices(subvertexIS, &subvertices);CHKERRQ(ierr);
1249   }
1250   ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr);
1251   ierr = PetscFree3(pStart,pEnd,pMax);CHKERRQ(ierr);
1252   PetscFunctionReturn(0);
1253 }
1254 
1255 #undef __FUNCT__
1256 #define __FUNCT__ "DMPlexMarkSubmesh_Interpolated"
1257 static PetscErrorCode DMPlexMarkSubmesh_Interpolated(DM dm, DMLabel vertexLabel, PetscInt value, DMLabel subpointMap, DM subdm)
1258 {
1259   IS               subvertexIS;
1260   const PetscInt  *subvertices;
1261   PetscInt        *pStart, *pEnd, *pMax;
1262   PetscInt         dim, d, numSubVerticesInitial = 0, v;
1263   PetscErrorCode   ierr;
1264 
1265   PetscFunctionBegin;
1266   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1267   ierr = PetscMalloc3(dim+1,PetscInt,&pStart,dim+1,PetscInt,&pEnd,dim+1,PetscInt,&pMax);CHKERRQ(ierr);
1268   ierr = DMPlexGetHybridBounds(dm, &pMax[dim], dim>1 ? &pMax[dim-1] : NULL, dim > 2 ? &pMax[1] : NULL, &pMax[0]);CHKERRQ(ierr);
1269   for (d = 0; d <= dim; ++d) {
1270     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr);
1271     if (pMax[d] >= 0) pEnd[d] = PetscMin(pEnd[d], pMax[d]);
1272   }
1273   /* Loop over initial vertices and mark all faces in the collective star() */
1274   ierr = DMLabelGetStratumIS(vertexLabel, value, &subvertexIS);CHKERRQ(ierr);
1275   if (subvertexIS) {
1276     ierr = ISGetSize(subvertexIS, &numSubVerticesInitial);CHKERRQ(ierr);
1277     ierr = ISGetIndices(subvertexIS, &subvertices);CHKERRQ(ierr);
1278   }
1279   for (v = 0; v < numSubVerticesInitial; ++v) {
1280     const PetscInt vertex = subvertices[v];
1281     PetscInt      *star   = NULL;
1282     PetscInt       starSize, s, numFaces = 0, f;
1283 
1284     ierr = DMPlexGetTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
1285     for (s = 0; s < starSize*2; s += 2) {
1286       const PetscInt point = star[s];
1287       if ((point >= pStart[dim-1]) && (point < pEnd[dim-1])) star[numFaces++] = point;
1288     }
1289     for (f = 0; f < numFaces; ++f) {
1290       const PetscInt face    = star[f];
1291       PetscInt      *closure = NULL;
1292       PetscInt       closureSize, c;
1293       PetscInt       faceLoc;
1294 
1295       ierr = DMLabelGetValue(subpointMap, face, &faceLoc);CHKERRQ(ierr);
1296       if (faceLoc == dim-1) continue;
1297       if (faceLoc >= 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Face %d has dimension %d in the surface label", face, faceLoc);
1298       ierr = DMPlexGetTransitiveClosure(dm, face, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1299       for (c = 0; c < closureSize*2; c += 2) {
1300         const PetscInt point = closure[c];
1301         PetscInt       vertexLoc;
1302 
1303         if ((point >= pStart[0]) && (point < pEnd[0])) {
1304           ierr = DMLabelGetValue(vertexLabel, point, &vertexLoc);CHKERRQ(ierr);
1305           if (vertexLoc != value) break;
1306         }
1307       }
1308       if (c == closureSize*2) {
1309         const PetscInt *support;
1310         PetscInt        supportSize, s;
1311 
1312         for (c = 0; c < closureSize*2; c += 2) {
1313           const PetscInt point = closure[c];
1314 
1315           for (d = 0; d < dim; ++d) {
1316             if ((point >= pStart[d]) && (point < pEnd[d])) {
1317               ierr = DMLabelSetValue(subpointMap, point, d);CHKERRQ(ierr);
1318               break;
1319             }
1320           }
1321         }
1322         ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
1323         ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr);
1324         for (s = 0; s < supportSize; ++s) {
1325           ierr = DMLabelSetValue(subpointMap, support[s], dim);CHKERRQ(ierr);
1326         }
1327       }
1328       ierr = DMPlexRestoreTransitiveClosure(dm, face, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1329     }
1330     ierr = DMPlexRestoreTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
1331   }
1332   if (subvertexIS) {
1333     ierr = ISRestoreIndices(subvertexIS, &subvertices);CHKERRQ(ierr);
1334   }
1335   ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr);
1336   ierr = PetscFree3(pStart,pEnd,pMax);CHKERRQ(ierr);
1337   PetscFunctionReturn(0);
1338 }
1339 
1340 #undef __FUNCT__
1341 #define __FUNCT__ "DMPlexMarkCohesiveSubmesh_Uninterpolated"
1342 static PetscErrorCode DMPlexMarkCohesiveSubmesh_Uninterpolated(DM dm, PetscBool hasLagrange, DMLabel subpointMap, PetscInt *numFaces, PetscInt *nFV, PetscInt *subCells[], DM subdm)
1343 {
1344   const PetscInt *cone;
1345   PetscInt        dim, cMax, cEnd, c, p, coneSize;
1346   PetscErrorCode  ierr;
1347 
1348   PetscFunctionBegin;
1349   *numFaces = 0;
1350   *nFV = 0;
1351   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1352   ierr = DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);CHKERRQ(ierr);
1353   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
1354   if (cMax < 0) PetscFunctionReturn(0);
1355   ierr = DMPlexGetConeSize(dm, cMax, &coneSize);CHKERRQ(ierr);
1356   *numFaces = cEnd - cMax;
1357   *nFV      = hasLagrange ? coneSize/3 : coneSize/2;
1358   ierr = PetscMalloc(*numFaces *2 * sizeof(PetscInt), subCells);CHKERRQ(ierr);
1359   for (c = cMax; c < cEnd; ++c) {
1360     const PetscInt *cells;
1361     PetscInt        numCells;
1362 
1363     ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr);
1364     for (p = 0; p < *nFV; ++p) {
1365       ierr = DMLabelSetValue(subpointMap, cone[p], 0);CHKERRQ(ierr);
1366     }
1367     /* Negative face */
1368     ierr = DMPlexGetJoin(dm, *nFV, cone, &numCells, &cells);CHKERRQ(ierr);
1369     if (numCells != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells");
1370     for (p = 0; p < numCells; ++p) {
1371       ierr = DMLabelSetValue(subpointMap, cells[p], 2);CHKERRQ(ierr);
1372       (*subCells)[(c-cMax)*2+p] = cells[p];
1373     }
1374     ierr = DMPlexRestoreJoin(dm, *nFV, cone, &numCells, &cells);CHKERRQ(ierr);
1375     /* Positive face is not included */
1376   }
1377   PetscFunctionReturn(0);
1378 }
1379 
1380 #undef __FUNCT__
1381 #define __FUNCT__ "DMPlexMarkCohesiveSubmesh_Interpolated"
1382 static PetscErrorCode DMPlexMarkCohesiveSubmesh_Interpolated(DM dm, PetscBool hasLagrange, DMLabel subpointMap, DM subdm)
1383 {
1384   PetscInt      *pStart, *pEnd;
1385   PetscInt       dim, cMax, cEnd, c, d;
1386   PetscErrorCode ierr;
1387 
1388   PetscFunctionBegin;
1389   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1390   ierr = DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);CHKERRQ(ierr);
1391   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
1392   if (cMax < 0) PetscFunctionReturn(0);
1393   ierr = PetscMalloc2(dim+1,PetscInt,&pStart,dim+1,PetscInt,&pEnd);CHKERRQ(ierr);
1394   for (d = 0; d <= dim; ++d) {
1395     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr);
1396   }
1397   for (c = cMax; c < cEnd; ++c) {
1398     const PetscInt *cone;
1399     PetscInt       *closure = NULL;
1400     PetscInt        coneSize, closureSize, cl;
1401 
1402     ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
1403     if (hasLagrange) {
1404       if (coneSize != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells");
1405     } else {
1406       if (coneSize != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells");
1407     }
1408     /* Negative face */
1409     ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr);
1410     ierr = DMPlexGetTransitiveClosure(dm, cone[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1411     for (cl = 0; cl < closureSize*2; cl += 2) {
1412       const PetscInt point = closure[cl];
1413 
1414       for (d = 0; d <= dim; ++d) {
1415         if ((point >= pStart[d]) && (point < pEnd[d])) {
1416           ierr = DMLabelSetValue(subpointMap, point, d);CHKERRQ(ierr);
1417           break;
1418         }
1419       }
1420     }
1421     ierr = DMPlexRestoreTransitiveClosure(dm, cone[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1422     /* Cells -- positive face is not included */
1423     for (cl = 0; cl < 1; ++cl) {
1424       const PetscInt *support;
1425       PetscInt        supportSize, s;
1426 
1427       ierr = DMPlexGetSupportSize(dm, cone[cl], &supportSize);CHKERRQ(ierr);
1428       if (supportSize != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive faces should separate two cells");
1429       ierr = DMPlexGetSupport(dm, cone[cl], &support);CHKERRQ(ierr);
1430       for (s = 0; s < supportSize; ++s) {
1431         ierr = DMLabelSetValue(subpointMap, support[s], dim);CHKERRQ(ierr);
1432       }
1433     }
1434   }
1435   ierr = PetscFree2(pStart, pEnd);CHKERRQ(ierr);
1436   PetscFunctionReturn(0);
1437 }
1438 
1439 #undef __FUNCT__
1440 #define __FUNCT__ "DMPlexGetFaceOrientation"
1441 PetscErrorCode DMPlexGetFaceOrientation(DM dm, PetscInt cell, PetscInt numCorners, PetscInt indices[], PetscInt oppositeVertex, PetscInt origVertices[], PetscInt faceVertices[], PetscBool *posOriented)
1442 {
1443   MPI_Comm       comm;
1444   PetscBool      posOrient = PETSC_FALSE;
1445   const PetscInt debug     = 0;
1446   PetscInt       cellDim, faceSize, f;
1447   PetscErrorCode ierr;
1448 
1449   PetscFunctionBegin;
1450   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1451   ierr = DMPlexGetDimension(dm, &cellDim);CHKERRQ(ierr);
1452   if (debug) {PetscPrintf(comm, "cellDim: %d numCorners: %d\n", cellDim, numCorners);CHKERRQ(ierr);}
1453 
1454   if (cellDim == numCorners-1) {
1455     /* Simplices */
1456     faceSize  = numCorners-1;
1457     posOrient = !(oppositeVertex%2) ? PETSC_TRUE : PETSC_FALSE;
1458   } else if (cellDim == 1 && numCorners == 3) {
1459     /* Quadratic line */
1460     faceSize  = 1;
1461     posOrient = PETSC_TRUE;
1462   } else if (cellDim == 2 && numCorners == 4) {
1463     /* Quads */
1464     faceSize = 2;
1465     if ((indices[1] > indices[0]) && (indices[1] - indices[0] == 1)) {
1466       posOrient = PETSC_TRUE;
1467     } else if ((indices[0] == 3) && (indices[1] == 0)) {
1468       posOrient = PETSC_TRUE;
1469     } else {
1470       if (((indices[0] > indices[1]) && (indices[0] - indices[1] == 1)) || ((indices[0] == 0) && (indices[1] == 3))) {
1471         posOrient = PETSC_FALSE;
1472       } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid quad crossedge");
1473     }
1474   } else if (cellDim == 2 && numCorners == 6) {
1475     /* Quadratic triangle (I hate this) */
1476     /* Edges are determined by the first 2 vertices (corners of edges) */
1477     const PetscInt faceSizeTri = 3;
1478     PetscInt       sortedIndices[3], i, iFace;
1479     PetscBool      found                    = PETSC_FALSE;
1480     PetscInt       faceVerticesTriSorted[9] = {
1481       0, 3,  4, /* bottom */
1482       1, 4,  5, /* right */
1483       2, 3,  5, /* left */
1484     };
1485     PetscInt       faceVerticesTri[9] = {
1486       0, 3,  4, /* bottom */
1487       1, 4,  5, /* right */
1488       2, 5,  3, /* left */
1489     };
1490 
1491     faceSize = faceSizeTri;
1492     for (i = 0; i < faceSizeTri; ++i) sortedIndices[i] = indices[i];
1493     ierr = PetscSortInt(faceSizeTri, sortedIndices);CHKERRQ(ierr);
1494     for (iFace = 0; iFace < 3; ++iFace) {
1495       const PetscInt ii = iFace*faceSizeTri;
1496       PetscInt       fVertex, cVertex;
1497 
1498       if ((sortedIndices[0] == faceVerticesTriSorted[ii+0]) &&
1499           (sortedIndices[1] == faceVerticesTriSorted[ii+1])) {
1500         for (fVertex = 0; fVertex < faceSizeTri; ++fVertex) {
1501           for (cVertex = 0; cVertex < faceSizeTri; ++cVertex) {
1502             if (indices[cVertex] == faceVerticesTri[ii+fVertex]) {
1503               faceVertices[fVertex] = origVertices[cVertex];
1504               break;
1505             }
1506           }
1507         }
1508         found = PETSC_TRUE;
1509         break;
1510       }
1511     }
1512     if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid tri crossface");
1513     if (posOriented) *posOriented = PETSC_TRUE;
1514     PetscFunctionReturn(0);
1515   } else if (cellDim == 2 && numCorners == 9) {
1516     /* Quadratic quad (I hate this) */
1517     /* Edges are determined by the first 2 vertices (corners of edges) */
1518     const PetscInt faceSizeQuad = 3;
1519     PetscInt       sortedIndices[3], i, iFace;
1520     PetscBool      found                      = PETSC_FALSE;
1521     PetscInt       faceVerticesQuadSorted[12] = {
1522       0, 1,  4, /* bottom */
1523       1, 2,  5, /* right */
1524       2, 3,  6, /* top */
1525       0, 3,  7, /* left */
1526     };
1527     PetscInt       faceVerticesQuad[12] = {
1528       0, 1,  4, /* bottom */
1529       1, 2,  5, /* right */
1530       2, 3,  6, /* top */
1531       3, 0,  7, /* left */
1532     };
1533 
1534     faceSize = faceSizeQuad;
1535     for (i = 0; i < faceSizeQuad; ++i) sortedIndices[i] = indices[i];
1536     ierr = PetscSortInt(faceSizeQuad, sortedIndices);CHKERRQ(ierr);
1537     for (iFace = 0; iFace < 4; ++iFace) {
1538       const PetscInt ii = iFace*faceSizeQuad;
1539       PetscInt       fVertex, cVertex;
1540 
1541       if ((sortedIndices[0] == faceVerticesQuadSorted[ii+0]) &&
1542           (sortedIndices[1] == faceVerticesQuadSorted[ii+1])) {
1543         for (fVertex = 0; fVertex < faceSizeQuad; ++fVertex) {
1544           for (cVertex = 0; cVertex < faceSizeQuad; ++cVertex) {
1545             if (indices[cVertex] == faceVerticesQuad[ii+fVertex]) {
1546               faceVertices[fVertex] = origVertices[cVertex];
1547               break;
1548             }
1549           }
1550         }
1551         found = PETSC_TRUE;
1552         break;
1553       }
1554     }
1555     if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid quad crossface");
1556     if (posOriented) *posOriented = PETSC_TRUE;
1557     PetscFunctionReturn(0);
1558   } else if (cellDim == 3 && numCorners == 8) {
1559     /* Hexes
1560        A hex is two oriented quads with the normal of the first
1561        pointing up at the second.
1562 
1563           7---6
1564          /|  /|
1565         4---5 |
1566         | 3-|-2
1567         |/  |/
1568         0---1
1569 
1570         Faces are determined by the first 4 vertices (corners of faces) */
1571     const PetscInt faceSizeHex = 4;
1572     PetscInt       sortedIndices[4], i, iFace;
1573     PetscBool      found                     = PETSC_FALSE;
1574     PetscInt       faceVerticesHexSorted[24] = {
1575       0, 1, 2, 3,  /* bottom */
1576       4, 5, 6, 7,  /* top */
1577       0, 1, 4, 5,  /* front */
1578       1, 2, 5, 6,  /* right */
1579       2, 3, 6, 7,  /* back */
1580       0, 3, 4, 7,  /* left */
1581     };
1582     PetscInt       faceVerticesHex[24] = {
1583       3, 2, 1, 0,  /* bottom */
1584       4, 5, 6, 7,  /* top */
1585       0, 1, 5, 4,  /* front */
1586       1, 2, 6, 5,  /* right */
1587       2, 3, 7, 6,  /* back */
1588       3, 0, 4, 7,  /* left */
1589     };
1590 
1591     faceSize = faceSizeHex;
1592     for (i = 0; i < faceSizeHex; ++i) sortedIndices[i] = indices[i];
1593     ierr = PetscSortInt(faceSizeHex, sortedIndices);CHKERRQ(ierr);
1594     for (iFace = 0; iFace < 6; ++iFace) {
1595       const PetscInt ii = iFace*faceSizeHex;
1596       PetscInt       fVertex, cVertex;
1597 
1598       if ((sortedIndices[0] == faceVerticesHexSorted[ii+0]) &&
1599           (sortedIndices[1] == faceVerticesHexSorted[ii+1]) &&
1600           (sortedIndices[2] == faceVerticesHexSorted[ii+2]) &&
1601           (sortedIndices[3] == faceVerticesHexSorted[ii+3])) {
1602         for (fVertex = 0; fVertex < faceSizeHex; ++fVertex) {
1603           for (cVertex = 0; cVertex < faceSizeHex; ++cVertex) {
1604             if (indices[cVertex] == faceVerticesHex[ii+fVertex]) {
1605               faceVertices[fVertex] = origVertices[cVertex];
1606               break;
1607             }
1608           }
1609         }
1610         found = PETSC_TRUE;
1611         break;
1612       }
1613     }
1614     if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid hex crossface");
1615     if (posOriented) *posOriented = PETSC_TRUE;
1616     PetscFunctionReturn(0);
1617   } else if (cellDim == 3 && numCorners == 10) {
1618     /* Quadratic tet */
1619     /* Faces are determined by the first 3 vertices (corners of faces) */
1620     const PetscInt faceSizeTet = 6;
1621     PetscInt       sortedIndices[6], i, iFace;
1622     PetscBool      found                     = PETSC_FALSE;
1623     PetscInt       faceVerticesTetSorted[24] = {
1624       0, 1, 2,  6, 7, 8, /* bottom */
1625       0, 3, 4,  6, 7, 9,  /* front */
1626       1, 4, 5,  7, 8, 9,  /* right */
1627       2, 3, 5,  6, 8, 9,  /* left */
1628     };
1629     PetscInt       faceVerticesTet[24] = {
1630       0, 1, 2,  6, 7, 8, /* bottom */
1631       0, 4, 3,  6, 7, 9,  /* front */
1632       1, 5, 4,  7, 8, 9,  /* right */
1633       2, 3, 5,  8, 6, 9,  /* left */
1634     };
1635 
1636     faceSize = faceSizeTet;
1637     for (i = 0; i < faceSizeTet; ++i) sortedIndices[i] = indices[i];
1638     ierr = PetscSortInt(faceSizeTet, sortedIndices);CHKERRQ(ierr);
1639     for (iFace=0; iFace < 4; ++iFace) {
1640       const PetscInt ii = iFace*faceSizeTet;
1641       PetscInt       fVertex, cVertex;
1642 
1643       if ((sortedIndices[0] == faceVerticesTetSorted[ii+0]) &&
1644           (sortedIndices[1] == faceVerticesTetSorted[ii+1]) &&
1645           (sortedIndices[2] == faceVerticesTetSorted[ii+2]) &&
1646           (sortedIndices[3] == faceVerticesTetSorted[ii+3])) {
1647         for (fVertex = 0; fVertex < faceSizeTet; ++fVertex) {
1648           for (cVertex = 0; cVertex < faceSizeTet; ++cVertex) {
1649             if (indices[cVertex] == faceVerticesTet[ii+fVertex]) {
1650               faceVertices[fVertex] = origVertices[cVertex];
1651               break;
1652             }
1653           }
1654         }
1655         found = PETSC_TRUE;
1656         break;
1657       }
1658     }
1659     if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid tet crossface");
1660     if (posOriented) *posOriented = PETSC_TRUE;
1661     PetscFunctionReturn(0);
1662   } else if (cellDim == 3 && numCorners == 27) {
1663     /* Quadratic hexes (I hate this)
1664        A hex is two oriented quads with the normal of the first
1665        pointing up at the second.
1666 
1667          7---6
1668         /|  /|
1669        4---5 |
1670        | 3-|-2
1671        |/  |/
1672        0---1
1673 
1674        Faces are determined by the first 4 vertices (corners of faces) */
1675     const PetscInt faceSizeQuadHex = 9;
1676     PetscInt       sortedIndices[9], i, iFace;
1677     PetscBool      found                         = PETSC_FALSE;
1678     PetscInt       faceVerticesQuadHexSorted[54] = {
1679       0, 1, 2, 3,  8, 9, 10, 11,  24, /* bottom */
1680       4, 5, 6, 7,  12, 13, 14, 15,  25, /* top */
1681       0, 1, 4, 5,  8, 12, 16, 17,  22, /* front */
1682       1, 2, 5, 6,  9, 13, 17, 18,  21, /* right */
1683       2, 3, 6, 7,  10, 14, 18, 19,  23, /* back */
1684       0, 3, 4, 7,  11, 15, 16, 19,  20, /* left */
1685     };
1686     PetscInt       faceVerticesQuadHex[54] = {
1687       3, 2, 1, 0,  10, 9, 8, 11,  24, /* bottom */
1688       4, 5, 6, 7,  12, 13, 14, 15,  25, /* top */
1689       0, 1, 5, 4,  8, 17, 12, 16,  22, /* front */
1690       1, 2, 6, 5,  9, 18, 13, 17,  21, /* right */
1691       2, 3, 7, 6,  10, 19, 14, 18,  23, /* back */
1692       3, 0, 4, 7,  11, 16, 15, 19,  20 /* left */
1693     };
1694 
1695     faceSize = faceSizeQuadHex;
1696     for (i = 0; i < faceSizeQuadHex; ++i) sortedIndices[i] = indices[i];
1697     ierr = PetscSortInt(faceSizeQuadHex, sortedIndices);CHKERRQ(ierr);
1698     for (iFace = 0; iFace < 6; ++iFace) {
1699       const PetscInt ii = iFace*faceSizeQuadHex;
1700       PetscInt       fVertex, cVertex;
1701 
1702       if ((sortedIndices[0] == faceVerticesQuadHexSorted[ii+0]) &&
1703           (sortedIndices[1] == faceVerticesQuadHexSorted[ii+1]) &&
1704           (sortedIndices[2] == faceVerticesQuadHexSorted[ii+2]) &&
1705           (sortedIndices[3] == faceVerticesQuadHexSorted[ii+3])) {
1706         for (fVertex = 0; fVertex < faceSizeQuadHex; ++fVertex) {
1707           for (cVertex = 0; cVertex < faceSizeQuadHex; ++cVertex) {
1708             if (indices[cVertex] == faceVerticesQuadHex[ii+fVertex]) {
1709               faceVertices[fVertex] = origVertices[cVertex];
1710               break;
1711             }
1712           }
1713         }
1714         found = PETSC_TRUE;
1715         break;
1716       }
1717     }
1718     if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid hex crossface");
1719     if (posOriented) *posOriented = PETSC_TRUE;
1720     PetscFunctionReturn(0);
1721   } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Unknown cell type for faceOrientation().");
1722   if (!posOrient) {
1723     if (debug) {ierr = PetscPrintf(comm, "  Reversing initial face orientation\n");CHKERRQ(ierr);}
1724     for (f = 0; f < faceSize; ++f) faceVertices[f] = origVertices[faceSize-1 - f];
1725   } else {
1726     if (debug) {ierr = PetscPrintf(comm, "  Keeping initial face orientation\n");CHKERRQ(ierr);}
1727     for (f = 0; f < faceSize; ++f) faceVertices[f] = origVertices[f];
1728   }
1729   if (posOriented) *posOriented = posOrient;
1730   PetscFunctionReturn(0);
1731 }
1732 
1733 #undef __FUNCT__
1734 #define __FUNCT__ "DMPlexGetOrientedFace"
1735 /*
1736     Given a cell and a face, as a set of vertices,
1737       return the oriented face, as a set of vertices, in faceVertices
1738     The orientation is such that the face normal points out of the cell
1739 */
1740 PetscErrorCode DMPlexGetOrientedFace(DM dm, PetscInt cell, PetscInt faceSize, const PetscInt face[], PetscInt numCorners, PetscInt indices[], PetscInt origVertices[], PetscInt faceVertices[], PetscBool *posOriented)
1741 {
1742   const PetscInt *cone = NULL;
1743   PetscInt        coneSize, v, f, v2;
1744   PetscInt        oppositeVertex = -1;
1745   PetscErrorCode  ierr;
1746 
1747   PetscFunctionBegin;
1748   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
1749   ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
1750   for (v = 0, v2 = 0; v < coneSize; ++v) {
1751     PetscBool found = PETSC_FALSE;
1752 
1753     for (f = 0; f < faceSize; ++f) {
1754       if (face[f] == cone[v]) {
1755         found = PETSC_TRUE; break;
1756       }
1757     }
1758     if (found) {
1759       indices[v2]      = v;
1760       origVertices[v2] = cone[v];
1761       ++v2;
1762     } else {
1763       oppositeVertex = v;
1764     }
1765   }
1766   ierr = DMPlexGetFaceOrientation(dm, cell, numCorners, indices, oppositeVertex, origVertices, faceVertices, posOriented);CHKERRQ(ierr);
1767   PetscFunctionReturn(0);
1768 }
1769 
1770 #undef __FUNCT__
1771 #define __FUNCT__ "DMPlexInsertFace_Internal"
1772 /*
1773   DMPlexInsertFace_Internal - Puts a face into the mesh
1774 
1775   Not collective
1776 
1777   Input Parameters:
1778   + dm              - The DMPlex
1779   . numFaceVertex   - The number of vertices in the face
1780   . faceVertices    - The vertices in the face for dm
1781   . subfaceVertices - The vertices in the face for subdm
1782   . numCorners      - The number of vertices in the cell
1783   . cell            - A cell in dm containing the face
1784   . subcell         - A cell in subdm containing the face
1785   . firstFace       - First face in the mesh
1786   - newFacePoint    - Next face in the mesh
1787 
1788   Output Parameters:
1789   . newFacePoint - Contains next face point number on input, updated on output
1790 
1791   Level: developer
1792 */
1793 static PetscErrorCode DMPlexInsertFace_Internal(DM dm, DM subdm, PetscInt numFaceVertices, const PetscInt faceVertices[], const PetscInt subfaceVertices[], PetscInt numCorners, PetscInt cell, PetscInt subcell, PetscInt firstFace, PetscInt *newFacePoint)
1794 {
1795   MPI_Comm        comm;
1796   DM_Plex        *submesh = (DM_Plex*) subdm->data;
1797   const PetscInt *faces;
1798   PetscInt        numFaces, coneSize;
1799   PetscErrorCode  ierr;
1800 
1801   PetscFunctionBegin;
1802   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1803   ierr = DMPlexGetConeSize(subdm, subcell, &coneSize);CHKERRQ(ierr);
1804   if (coneSize != 1) SETERRQ2(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cone size of cell %d is %d != 1", cell, coneSize);
1805 #if 0
1806   /* Cannot use this because support() has not been constructed yet */
1807   ierr = DMPlexGetJoin(subdm, numFaceVertices, subfaceVertices, &numFaces, &faces);CHKERRQ(ierr);
1808 #else
1809   {
1810     PetscInt f;
1811 
1812     numFaces = 0;
1813     ierr     = DMGetWorkArray(subdm, 1, PETSC_INT, (void **) &faces);CHKERRQ(ierr);
1814     for (f = firstFace; f < *newFacePoint; ++f) {
1815       PetscInt dof, off, d;
1816 
1817       ierr = PetscSectionGetDof(submesh->coneSection, f, &dof);CHKERRQ(ierr);
1818       ierr = PetscSectionGetOffset(submesh->coneSection, f, &off);CHKERRQ(ierr);
1819       /* Yes, I know this is quadratic, but I expect the sizes to be <5 */
1820       for (d = 0; d < dof; ++d) {
1821         const PetscInt p = submesh->cones[off+d];
1822         PetscInt       v;
1823 
1824         for (v = 0; v < numFaceVertices; ++v) {
1825           if (subfaceVertices[v] == p) break;
1826         }
1827         if (v == numFaceVertices) break;
1828       }
1829       if (d == dof) {
1830         numFaces               = 1;
1831         ((PetscInt*) faces)[0] = f;
1832       }
1833     }
1834   }
1835 #endif
1836   if (numFaces > 1) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Vertex set had %d faces, not one", numFaces);
1837   else if (numFaces == 1) {
1838     /* Add the other cell neighbor for this face */
1839     ierr = DMPlexSetCone(subdm, subcell, faces);CHKERRQ(ierr);
1840   } else {
1841     PetscInt *indices, *origVertices, *orientedVertices, *orientedSubVertices, v, ov;
1842     PetscBool posOriented;
1843 
1844     ierr                = DMGetWorkArray(subdm, 4*numFaceVertices * sizeof(PetscInt), PETSC_INT, &orientedVertices);CHKERRQ(ierr);
1845     origVertices        = &orientedVertices[numFaceVertices];
1846     indices             = &orientedVertices[numFaceVertices*2];
1847     orientedSubVertices = &orientedVertices[numFaceVertices*3];
1848     ierr                = DMPlexGetOrientedFace(dm, cell, numFaceVertices, faceVertices, numCorners, indices, origVertices, orientedVertices, &posOriented);CHKERRQ(ierr);
1849     /* TODO: I know that routine should return a permutation, not the indices */
1850     for (v = 0; v < numFaceVertices; ++v) {
1851       const PetscInt vertex = faceVertices[v], subvertex = subfaceVertices[v];
1852       for (ov = 0; ov < numFaceVertices; ++ov) {
1853         if (orientedVertices[ov] == vertex) {
1854           orientedSubVertices[ov] = subvertex;
1855           break;
1856         }
1857       }
1858       if (ov == numFaceVertices) SETERRQ1(comm, PETSC_ERR_PLIB, "Could not find face vertex %d in orientated set", vertex);
1859     }
1860     ierr = DMPlexSetCone(subdm, *newFacePoint, orientedSubVertices);CHKERRQ(ierr);
1861     ierr = DMPlexSetCone(subdm, subcell, newFacePoint);CHKERRQ(ierr);
1862     ierr = DMRestoreWorkArray(subdm, 4*numFaceVertices * sizeof(PetscInt), PETSC_INT, &orientedVertices);CHKERRQ(ierr);
1863     ++(*newFacePoint);
1864   }
1865 #if 0
1866   ierr = DMPlexRestoreJoin(subdm, numFaceVertices, subfaceVertices, &numFaces, &faces);CHKERRQ(ierr);
1867 #else
1868   ierr = DMRestoreWorkArray(subdm, 1, PETSC_INT, (void **) &faces);CHKERRQ(ierr);
1869 #endif
1870   PetscFunctionReturn(0);
1871 }
1872 
1873 #undef __FUNCT__
1874 #define __FUNCT__ "DMPlexCreateSubmesh_Uninterpolated"
1875 static PetscErrorCode DMPlexCreateSubmesh_Uninterpolated(DM dm, const char vertexLabelName[], PetscInt value, DM subdm)
1876 {
1877   MPI_Comm        comm;
1878   DMLabel         vertexLabel, subpointMap;
1879   IS              subvertexIS,  subcellIS;
1880   const PetscInt *subVertices, *subCells;
1881   PetscInt        numSubVertices, firstSubVertex, numSubCells;
1882   PetscInt       *subface, maxConeSize, numSubFaces, firstSubFace, newFacePoint, nFV;
1883   PetscInt        vStart, vEnd, c, f;
1884   PetscErrorCode  ierr;
1885 
1886   PetscFunctionBegin;
1887   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1888   /* Create subpointMap which marks the submesh */
1889   ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr);
1890   ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr);
1891   ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr);
1892   ierr = DMPlexGetLabel(dm, vertexLabelName, &vertexLabel);CHKERRQ(ierr);
1893   ierr = DMPlexMarkSubmesh_Uninterpolated(dm, vertexLabel, value, subpointMap, &numSubFaces, &nFV, subdm);CHKERRQ(ierr);
1894   /* Setup chart */
1895   ierr = DMLabelGetStratumSize(subpointMap, 0, &numSubVertices);CHKERRQ(ierr);
1896   ierr = DMLabelGetStratumSize(subpointMap, 2, &numSubCells);CHKERRQ(ierr);
1897   ierr = DMPlexSetChart(subdm, 0, numSubCells+numSubFaces+numSubVertices);CHKERRQ(ierr);
1898   ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr);
1899   /* Set cone sizes */
1900   firstSubVertex = numSubCells;
1901   firstSubFace   = numSubCells+numSubVertices;
1902   newFacePoint   = firstSubFace;
1903   ierr = DMLabelGetStratumIS(subpointMap, 0, &subvertexIS);CHKERRQ(ierr);
1904   if (subvertexIS) {ierr = ISGetIndices(subvertexIS, &subVertices);CHKERRQ(ierr);}
1905   ierr = DMLabelGetStratumIS(subpointMap, 2, &subcellIS);CHKERRQ(ierr);
1906   if (subcellIS) {ierr = ISGetIndices(subcellIS, &subCells);CHKERRQ(ierr);}
1907   for (c = 0; c < numSubCells; ++c) {
1908     ierr = DMPlexSetConeSize(subdm, c, 1);CHKERRQ(ierr);
1909   }
1910   for (f = firstSubFace; f < firstSubFace+numSubFaces; ++f) {
1911     ierr = DMPlexSetConeSize(subdm, f, nFV);CHKERRQ(ierr);
1912   }
1913   ierr = DMSetUp(subdm);CHKERRQ(ierr);
1914   /* Create face cones */
1915   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
1916   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
1917   ierr = DMGetWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr);
1918   for (c = 0; c < numSubCells; ++c) {
1919     const PetscInt cell    = subCells[c];
1920     const PetscInt subcell = c;
1921     PetscInt      *closure = NULL;
1922     PetscInt       closureSize, cl, numCorners = 0, faceSize = 0;
1923 
1924     ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1925     for (cl = 0; cl < closureSize*2; cl += 2) {
1926       const PetscInt point = closure[cl];
1927       PetscInt       subVertex;
1928 
1929       if ((point >= vStart) && (point < vEnd)) {
1930         ++numCorners;
1931         ierr = PetscFindInt(point, numSubVertices, subVertices, &subVertex);CHKERRQ(ierr);
1932         if (subVertex >= 0) {
1933           closure[faceSize] = point;
1934           subface[faceSize] = firstSubVertex+subVertex;
1935           ++faceSize;
1936         }
1937       }
1938     }
1939     if (faceSize > nFV) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Invalid submesh: Too many vertices %d of an element on the surface", faceSize);
1940     if (faceSize == nFV) {
1941       ierr = DMPlexInsertFace_Internal(dm, subdm, faceSize, closure, subface, numCorners, cell, subcell, firstSubFace, &newFacePoint);CHKERRQ(ierr);
1942     }
1943     ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1944   }
1945   ierr = DMRestoreWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr);
1946   ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr);
1947   ierr = DMPlexStratify(subdm);CHKERRQ(ierr);
1948   /* Build coordinates */
1949   {
1950     PetscSection coordSection, subCoordSection;
1951     Vec          coordinates, subCoordinates;
1952     PetscScalar *coords, *subCoords;
1953     PetscInt     numComp, coordSize, v;
1954 
1955     ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1956     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1957     ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr);
1958     ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr);
1959     ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr);
1960     ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr);
1961     ierr = PetscSectionSetChart(subCoordSection, firstSubVertex, firstSubVertex+numSubVertices);CHKERRQ(ierr);
1962     for (v = 0; v < numSubVertices; ++v) {
1963       const PetscInt vertex    = subVertices[v];
1964       const PetscInt subvertex = firstSubVertex+v;
1965       PetscInt       dof;
1966 
1967       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
1968       ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr);
1969       ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr);
1970     }
1971     ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr);
1972     ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr);
1973     ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr);
1974     ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
1975     ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr);
1976     if (coordSize) {
1977       ierr = VecGetArray(coordinates,    &coords);CHKERRQ(ierr);
1978       ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr);
1979       for (v = 0; v < numSubVertices; ++v) {
1980         const PetscInt vertex    = subVertices[v];
1981         const PetscInt subvertex = firstSubVertex+v;
1982         PetscInt       dof, off, sdof, soff, d;
1983 
1984         ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
1985         ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr);
1986         ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr);
1987         ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr);
1988         if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof);
1989         for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d];
1990       }
1991       ierr = VecRestoreArray(coordinates,    &coords);CHKERRQ(ierr);
1992       ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr);
1993     }
1994     ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr);
1995     ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr);
1996   }
1997   /* Cleanup */
1998   if (subvertexIS) {ierr = ISRestoreIndices(subvertexIS, &subVertices);CHKERRQ(ierr);}
1999   ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr);
2000   if (subcellIS) {ierr = ISRestoreIndices(subcellIS, &subCells);CHKERRQ(ierr);}
2001   ierr = ISDestroy(&subcellIS);CHKERRQ(ierr);
2002   PetscFunctionReturn(0);
2003 }
2004 
2005 #undef __FUNCT__
2006 #define __FUNCT__ "DMPlexCreateSubmesh_Interpolated"
2007 static PetscErrorCode DMPlexCreateSubmesh_Interpolated(DM dm, const char vertexLabelName[], PetscInt value, DM subdm)
2008 {
2009   MPI_Comm         comm;
2010   DMLabel          subpointMap, vertexLabel;
2011   IS              *subpointIS;
2012   const PetscInt **subpoints;
2013   PetscInt        *numSubPoints, *firstSubPoint, *coneNew, *coneONew;
2014   PetscInt         totSubPoints = 0, maxConeSize, dim, p, d, v;
2015   PetscErrorCode   ierr;
2016 
2017   PetscFunctionBegin;
2018   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
2019   /* Create subpointMap which marks the submesh */
2020   ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr);
2021   ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr);
2022   ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr);
2023   ierr = DMPlexGetLabel(dm, vertexLabelName, &vertexLabel);CHKERRQ(ierr);
2024   ierr = DMPlexMarkSubmesh_Interpolated(dm, vertexLabel, value, subpointMap, subdm);CHKERRQ(ierr);
2025   /* Setup chart */
2026   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
2027   ierr = PetscMalloc4(dim+1,PetscInt,&numSubPoints,dim+1,PetscInt,&firstSubPoint,dim+1,IS,&subpointIS,dim+1,const PetscInt *,&subpoints);CHKERRQ(ierr);
2028   for (d = 0; d <= dim; ++d) {
2029     ierr = DMLabelGetStratumSize(subpointMap, d, &numSubPoints[d]);CHKERRQ(ierr);
2030     totSubPoints += numSubPoints[d];
2031   }
2032   ierr = DMPlexSetChart(subdm, 0, totSubPoints);CHKERRQ(ierr);
2033   ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr);
2034   /* Set cone sizes */
2035   firstSubPoint[dim] = 0;
2036   firstSubPoint[0]   = firstSubPoint[dim] + numSubPoints[dim];
2037   if (dim > 1) {firstSubPoint[dim-1] = firstSubPoint[0]     + numSubPoints[0];}
2038   if (dim > 2) {firstSubPoint[dim-2] = firstSubPoint[dim-1] + numSubPoints[dim-1];}
2039   for (d = 0; d <= dim; ++d) {
2040     ierr = DMLabelGetStratumIS(subpointMap, d, &subpointIS[d]);CHKERRQ(ierr);
2041     ierr = ISGetIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr);
2042   }
2043   for (d = 0; d <= dim; ++d) {
2044     for (p = 0; p < numSubPoints[d]; ++p) {
2045       const PetscInt  point    = subpoints[d][p];
2046       const PetscInt  subpoint = firstSubPoint[d] + p;
2047       const PetscInt *cone;
2048       PetscInt        coneSize, coneSizeNew, c, val;
2049 
2050       ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
2051       ierr = DMPlexSetConeSize(subdm, subpoint, coneSize);CHKERRQ(ierr);
2052       if (d == dim) {
2053         ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
2054         for (c = 0, coneSizeNew = 0; c < coneSize; ++c) {
2055           ierr = DMLabelGetValue(subpointMap, cone[c], &val);CHKERRQ(ierr);
2056           if (val >= 0) coneSizeNew++;
2057         }
2058         ierr = DMPlexSetConeSize(subdm, subpoint, coneSizeNew);CHKERRQ(ierr);
2059       }
2060     }
2061   }
2062   ierr = DMSetUp(subdm);CHKERRQ(ierr);
2063   /* Set cones */
2064   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
2065   ierr = PetscMalloc2(maxConeSize,PetscInt,&coneNew,maxConeSize,PetscInt,&coneONew);CHKERRQ(ierr);
2066   for (d = 0; d <= dim; ++d) {
2067     for (p = 0; p < numSubPoints[d]; ++p) {
2068       const PetscInt  point    = subpoints[d][p];
2069       const PetscInt  subpoint = firstSubPoint[d] + p;
2070       const PetscInt *cone, *ornt;
2071       PetscInt        coneSize, subconeSize, coneSizeNew, c, subc;
2072 
2073       ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
2074       ierr = DMPlexGetConeSize(subdm, subpoint, &subconeSize);CHKERRQ(ierr);
2075       ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
2076       ierr = DMPlexGetConeOrientation(dm, point, &ornt);CHKERRQ(ierr);
2077       for (c = 0, coneSizeNew = 0; c < coneSize; ++c) {
2078         ierr = PetscFindInt(cone[c], numSubPoints[d-1], subpoints[d-1], &subc);CHKERRQ(ierr);
2079         if (subc >= 0) {
2080           coneNew[coneSizeNew]  = firstSubPoint[d-1] + subc;
2081           coneONew[coneSizeNew] = ornt[c];
2082           ++coneSizeNew;
2083         }
2084       }
2085       if (coneSizeNew != subconeSize) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of cone points located %d does not match subcone size %d", coneSizeNew, subconeSize);
2086       ierr = DMPlexSetCone(subdm, subpoint, coneNew);CHKERRQ(ierr);
2087       ierr = DMPlexSetConeOrientation(subdm, subpoint, coneONew);CHKERRQ(ierr);
2088     }
2089   }
2090   ierr = PetscFree2(coneNew,coneONew);CHKERRQ(ierr);
2091   ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr);
2092   ierr = DMPlexStratify(subdm);CHKERRQ(ierr);
2093   /* Build coordinates */
2094   {
2095     PetscSection coordSection, subCoordSection;
2096     Vec          coordinates, subCoordinates;
2097     PetscScalar *coords, *subCoords;
2098     PetscInt     numComp, coordSize;
2099 
2100     ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
2101     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
2102     ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr);
2103     ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr);
2104     ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr);
2105     ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr);
2106     ierr = PetscSectionSetChart(subCoordSection, firstSubPoint[0], firstSubPoint[0]+numSubPoints[0]);CHKERRQ(ierr);
2107     for (v = 0; v < numSubPoints[0]; ++v) {
2108       const PetscInt vertex    = subpoints[0][v];
2109       const PetscInt subvertex = firstSubPoint[0]+v;
2110       PetscInt       dof;
2111 
2112       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
2113       ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr);
2114       ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr);
2115     }
2116     ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr);
2117     ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr);
2118     ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr);
2119     ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
2120     ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr);
2121     ierr = VecGetArray(coordinates,    &coords);CHKERRQ(ierr);
2122     ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr);
2123     for (v = 0; v < numSubPoints[0]; ++v) {
2124       const PetscInt vertex    = subpoints[0][v];
2125       const PetscInt subvertex = firstSubPoint[0]+v;
2126       PetscInt dof, off, sdof, soff, d;
2127 
2128       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
2129       ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr);
2130       ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr);
2131       ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr);
2132       if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof);
2133       for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d];
2134     }
2135     ierr = VecRestoreArray(coordinates,    &coords);CHKERRQ(ierr);
2136     ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr);
2137     ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr);
2138     ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr);
2139   }
2140   /* Cleanup */
2141   for (d = 0; d <= dim; ++d) {
2142     ierr = ISRestoreIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr);
2143     ierr = ISDestroy(&subpointIS[d]);CHKERRQ(ierr);
2144   }
2145   ierr = PetscFree4(numSubPoints,firstSubPoint,subpointIS,subpoints);CHKERRQ(ierr);
2146   PetscFunctionReturn(0);
2147 }
2148 
2149 #undef __FUNCT__
2150 #define __FUNCT__ "DMPlexCreateSubmesh"
2151 /*@C
2152   DMPlexCreateSubmesh - Extract a hypersurface from the mesh using vertices defined by a label
2153 
2154   Input Parameters:
2155 + dm           - The original mesh
2156 . vertexLabel  - The DMLabel marking vertices contained in the surface
2157 - value        - The label value to use
2158 
2159   Output Parameter:
2160 . subdm - The surface mesh
2161 
2162   Note: This function produces a DMLabel mapping original points in the submesh to their depth. This can be obtained using DMPlexGetSubpointMap().
2163 
2164   Level: developer
2165 
2166 .seealso: DMPlexGetSubpointMap(), DMPlexGetLabel(), DMLabelSetValue()
2167 @*/
2168 PetscErrorCode DMPlexCreateSubmesh(DM dm, const char vertexLabel[], PetscInt value, DM *subdm)
2169 {
2170   PetscInt       dim, depth;
2171   PetscErrorCode ierr;
2172 
2173   PetscFunctionBegin;
2174   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2175   PetscValidCharPointer(vertexLabel, 2);
2176   PetscValidPointer(subdm, 3);
2177   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
2178   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
2179   ierr = DMCreate(PetscObjectComm((PetscObject)dm), subdm);CHKERRQ(ierr);
2180   ierr = DMSetType(*subdm, DMPLEX);CHKERRQ(ierr);
2181   ierr = DMPlexSetDimension(*subdm, dim-1);CHKERRQ(ierr);
2182   if (depth == dim) {
2183     ierr = DMPlexCreateSubmesh_Interpolated(dm, vertexLabel, value, *subdm);CHKERRQ(ierr);
2184   } else {
2185     ierr = DMPlexCreateSubmesh_Uninterpolated(dm, vertexLabel, value, *subdm);CHKERRQ(ierr);
2186   }
2187   PetscFunctionReturn(0);
2188 }
2189 
2190 #undef __FUNCT__
2191 #define __FUNCT__ "DMPlexCreateCohesiveSubmesh_Uninterpolated"
2192 static PetscErrorCode DMPlexCreateCohesiveSubmesh_Uninterpolated(DM dm, PetscBool hasLagrange, DM subdm)
2193 {
2194   MPI_Comm        comm;
2195   DMLabel         subpointMap;
2196   IS              subvertexIS;
2197   const PetscInt *subVertices;
2198   PetscInt        numSubVertices, firstSubVertex, numSubCells, *subCells;
2199   PetscInt       *subface, maxConeSize, numSubFaces, firstSubFace, newFacePoint, nFV;
2200   PetscInt        cMax, c, f;
2201   PetscErrorCode  ierr;
2202 
2203   PetscFunctionBegin;
2204   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
2205   /* Create subpointMap which marks the submesh */
2206   ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr);
2207   ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr);
2208   ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr);
2209   ierr = DMPlexMarkCohesiveSubmesh_Uninterpolated(dm, hasLagrange, subpointMap, &numSubFaces, &nFV, &subCells, subdm);CHKERRQ(ierr);
2210   /* Setup chart */
2211   ierr = DMLabelGetStratumSize(subpointMap, 0, &numSubVertices);CHKERRQ(ierr);
2212   ierr = DMLabelGetStratumSize(subpointMap, 2, &numSubCells);CHKERRQ(ierr);
2213   ierr = DMPlexSetChart(subdm, 0, numSubCells+numSubFaces+numSubVertices);CHKERRQ(ierr);
2214   ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr);
2215   /* Set cone sizes */
2216   firstSubVertex = numSubCells;
2217   firstSubFace   = numSubCells+numSubVertices;
2218   newFacePoint   = firstSubFace;
2219   ierr = DMLabelGetStratumIS(subpointMap, 0, &subvertexIS);CHKERRQ(ierr);
2220   if (subvertexIS) {ierr = ISGetIndices(subvertexIS, &subVertices);CHKERRQ(ierr);}
2221   for (c = 0; c < numSubCells; ++c) {
2222     ierr = DMPlexSetConeSize(subdm, c, 1);CHKERRQ(ierr);
2223   }
2224   for (f = firstSubFace; f < firstSubFace+numSubFaces; ++f) {
2225     ierr = DMPlexSetConeSize(subdm, f, nFV);CHKERRQ(ierr);
2226   }
2227   ierr = DMSetUp(subdm);CHKERRQ(ierr);
2228   /* Create face cones */
2229   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
2230   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
2231   ierr = DMGetWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr);
2232   for (c = 0; c < numSubCells; ++c) {
2233     const PetscInt  cell    = subCells[c];
2234     const PetscInt  subcell = c;
2235     const PetscInt *cone, *cells;
2236     PetscInt        numCells, subVertex, p, v;
2237 
2238     if (cell < cMax) continue;
2239     ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
2240     for (v = 0; v < nFV; ++v) {
2241       ierr = PetscFindInt(cone[v], numSubVertices, subVertices, &subVertex);CHKERRQ(ierr);
2242       subface[v] = firstSubVertex+subVertex;
2243     }
2244     ierr = DMPlexSetCone(subdm, newFacePoint, subface);CHKERRQ(ierr);
2245     ierr = DMPlexSetCone(subdm, subcell, &newFacePoint);CHKERRQ(ierr);
2246     ierr = DMPlexGetJoin(dm, nFV, cone, &numCells, &cells);CHKERRQ(ierr);
2247     if (numCells != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells");
2248     for (p = 0; p < numCells; ++p) {
2249       PetscInt negsubcell;
2250 
2251       if (cells[p] >= cMax) continue;
2252       /* I know this is a crap search */
2253       for (negsubcell = 0; negsubcell < numSubCells; ++negsubcell) {
2254         if (subCells[negsubcell] == cells[p]) break;
2255       }
2256       if (negsubcell == numSubCells) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find negative face neighbor for cohesive cell %d", cell);
2257       ierr = DMPlexSetCone(subdm, negsubcell, &newFacePoint);CHKERRQ(ierr);
2258     }
2259     ierr = DMPlexRestoreJoin(dm, nFV, cone, &numCells, &cells);CHKERRQ(ierr);
2260     ++newFacePoint;
2261   }
2262   ierr = DMRestoreWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr);
2263   ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr);
2264   ierr = DMPlexStratify(subdm);CHKERRQ(ierr);
2265   /* Build coordinates */
2266   {
2267     PetscSection coordSection, subCoordSection;
2268     Vec          coordinates, subCoordinates;
2269     PetscScalar *coords, *subCoords;
2270     PetscInt     numComp, coordSize, v;
2271 
2272     ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
2273     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
2274     ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr);
2275     ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr);
2276     ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr);
2277     ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr);
2278     ierr = PetscSectionSetChart(subCoordSection, firstSubVertex, firstSubVertex+numSubVertices);CHKERRQ(ierr);
2279     for (v = 0; v < numSubVertices; ++v) {
2280       const PetscInt vertex    = subVertices[v];
2281       const PetscInt subvertex = firstSubVertex+v;
2282       PetscInt       dof;
2283 
2284       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
2285       ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr);
2286       ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr);
2287     }
2288     ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr);
2289     ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr);
2290     ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr);
2291     if (coordSize) {
2292       ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
2293       ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr);
2294       ierr = VecGetArray(coordinates,    &coords);CHKERRQ(ierr);
2295       ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr);
2296       for (v = 0; v < numSubVertices; ++v) {
2297         const PetscInt vertex    = subVertices[v];
2298         const PetscInt subvertex = firstSubVertex+v;
2299         PetscInt       dof, off, sdof, soff, d;
2300 
2301         ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
2302         ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr);
2303         ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr);
2304         ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr);
2305         if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof);
2306         for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d];
2307       }
2308       ierr = VecRestoreArray(coordinates,    &coords);CHKERRQ(ierr);
2309       ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr);
2310     }
2311     ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr);
2312     ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr);
2313   }
2314   /* Cleanup */
2315   if (subvertexIS) {ierr = ISRestoreIndices(subvertexIS, &subVertices);CHKERRQ(ierr);}
2316   ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr);
2317   ierr = PetscFree(subCells);CHKERRQ(ierr);
2318   PetscFunctionReturn(0);
2319 }
2320 
2321 #undef __FUNCT__
2322 #define __FUNCT__ "DMPlexCreateCohesiveSubmesh_Interpolated"
2323 static PetscErrorCode DMPlexCreateCohesiveSubmesh_Interpolated(DM dm, PetscBool hasLagrange, DM subdm)
2324 {
2325   MPI_Comm         comm;
2326   DMLabel          subpointMap;
2327   IS              *subpointIS;
2328   const PetscInt **subpoints;
2329   PetscInt        *numSubPoints, *firstSubPoint, *coneNew;
2330   PetscInt         totSubPoints = 0, maxConeSize, dim, p, d, v;
2331   PetscErrorCode   ierr;
2332 
2333   PetscFunctionBegin;
2334   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
2335   /* Create subpointMap which marks the submesh */
2336   ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr);
2337   ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr);
2338   ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr);
2339   ierr = DMPlexMarkCohesiveSubmesh_Interpolated(dm, hasLagrange, subpointMap, subdm);CHKERRQ(ierr);
2340   /* Setup chart */
2341   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
2342   ierr = PetscMalloc4(dim+1,PetscInt,&numSubPoints,dim+1,PetscInt,&firstSubPoint,dim+1,IS,&subpointIS,dim+1,const PetscInt *,&subpoints);CHKERRQ(ierr);
2343   for (d = 0; d <= dim; ++d) {
2344     ierr = DMLabelGetStratumSize(subpointMap, d, &numSubPoints[d]);CHKERRQ(ierr);
2345     totSubPoints += numSubPoints[d];
2346   }
2347   ierr = DMPlexSetChart(subdm, 0, totSubPoints);CHKERRQ(ierr);
2348   ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr);
2349   /* Set cone sizes */
2350   firstSubPoint[dim] = 0;
2351   firstSubPoint[0]   = firstSubPoint[dim] + numSubPoints[dim];
2352   if (dim > 1) {firstSubPoint[dim-1] = firstSubPoint[0]     + numSubPoints[0];}
2353   if (dim > 2) {firstSubPoint[dim-2] = firstSubPoint[dim-1] + numSubPoints[dim-1];}
2354   for (d = 0; d <= dim; ++d) {
2355     ierr = DMLabelGetStratumIS(subpointMap, d, &subpointIS[d]);CHKERRQ(ierr);
2356     ierr = ISGetIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr);
2357   }
2358   for (d = 0; d <= dim; ++d) {
2359     for (p = 0; p < numSubPoints[d]; ++p) {
2360       const PetscInt  point    = subpoints[d][p];
2361       const PetscInt  subpoint = firstSubPoint[d] + p;
2362       const PetscInt *cone;
2363       PetscInt        coneSize, coneSizeNew, c, val;
2364 
2365       ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
2366       ierr = DMPlexSetConeSize(subdm, subpoint, coneSize);CHKERRQ(ierr);
2367       if (d == dim) {
2368         ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
2369         for (c = 0, coneSizeNew = 0; c < coneSize; ++c) {
2370           ierr = DMLabelGetValue(subpointMap, cone[c], &val);CHKERRQ(ierr);
2371           if (val >= 0) coneSizeNew++;
2372         }
2373         ierr = DMPlexSetConeSize(subdm, subpoint, coneSizeNew);CHKERRQ(ierr);
2374       }
2375     }
2376   }
2377   ierr = DMSetUp(subdm);CHKERRQ(ierr);
2378   /* Set cones */
2379   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
2380   ierr = PetscMalloc(maxConeSize * sizeof(PetscInt), &coneNew);CHKERRQ(ierr);
2381   for (d = 0; d <= dim; ++d) {
2382     for (p = 0; p < numSubPoints[d]; ++p) {
2383       const PetscInt  point    = subpoints[d][p];
2384       const PetscInt  subpoint = firstSubPoint[d] + p;
2385       const PetscInt *cone;
2386       PetscInt        coneSize, subconeSize, coneSizeNew, c, subc;
2387 
2388       ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
2389       ierr = DMPlexGetConeSize(subdm, subpoint, &subconeSize);CHKERRQ(ierr);
2390       ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
2391       for (c = 0, coneSizeNew = 0; c < coneSize; ++c) {
2392         ierr = PetscFindInt(cone[c], numSubPoints[d-1], subpoints[d-1], &subc);CHKERRQ(ierr);
2393         if (subc >= 0) coneNew[coneSizeNew++] = firstSubPoint[d-1] + subc;
2394       }
2395       if (coneSizeNew != subconeSize) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of cone points located %d does not match subcone size %d", coneSizeNew, subconeSize);
2396       ierr = DMPlexSetCone(subdm, subpoint, coneNew);CHKERRQ(ierr);
2397     }
2398   }
2399   ierr = PetscFree(coneNew);CHKERRQ(ierr);
2400   ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr);
2401   ierr = DMPlexStratify(subdm);CHKERRQ(ierr);
2402   /* Build coordinates */
2403   {
2404     PetscSection coordSection, subCoordSection;
2405     Vec          coordinates, subCoordinates;
2406     PetscScalar *coords, *subCoords;
2407     PetscInt     numComp, coordSize;
2408 
2409     ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
2410     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
2411     ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr);
2412     ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr);
2413     ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr);
2414     ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr);
2415     ierr = PetscSectionSetChart(subCoordSection, firstSubPoint[0], firstSubPoint[0]+numSubPoints[0]);CHKERRQ(ierr);
2416     for (v = 0; v < numSubPoints[0]; ++v) {
2417       const PetscInt vertex    = subpoints[0][v];
2418       const PetscInt subvertex = firstSubPoint[0]+v;
2419       PetscInt       dof;
2420 
2421       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
2422       ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr);
2423       ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr);
2424     }
2425     ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr);
2426     ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr);
2427     ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr);
2428     ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
2429     ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr);
2430     ierr = VecGetArray(coordinates,    &coords);CHKERRQ(ierr);
2431     ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr);
2432     for (v = 0; v < numSubPoints[0]; ++v) {
2433       const PetscInt vertex    = subpoints[0][v];
2434       const PetscInt subvertex = firstSubPoint[0]+v;
2435       PetscInt dof, off, sdof, soff, d;
2436 
2437       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
2438       ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr);
2439       ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr);
2440       ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr);
2441       if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof);
2442       for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d];
2443     }
2444     ierr = VecRestoreArray(coordinates,    &coords);CHKERRQ(ierr);
2445     ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr);
2446     ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr);
2447     ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr);
2448   }
2449   /* Cleanup */
2450   for (d = 0; d <= dim; ++d) {
2451     ierr = ISRestoreIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr);
2452     ierr = ISDestroy(&subpointIS[d]);CHKERRQ(ierr);
2453   }
2454   ierr = PetscFree4(numSubPoints,firstSubPoint,subpointIS,subpoints);CHKERRQ(ierr);
2455   PetscFunctionReturn(0);
2456 }
2457 
2458 #undef __FUNCT__
2459 #define __FUNCT__ "DMPlexCreateCohesiveSubmesh"
2460 /*
2461   DMPlexCreateCohesiveSubmesh - Extract from a mesh with cohesive cells the hypersurface defined by one face of the cells.
2462 
2463   Input Parameters:
2464 + dm          - The original mesh
2465 - hasLagrange - The mesh has Lagrange unknowns in the cohesive cells
2466 
2467   Output Parameter:
2468 . subdm - The surface mesh
2469 
2470   Note: This function produces a DMLabel mapping original points in the submesh to their depth. This can be obtained using DMPlexGetSubpointMap().
2471 
2472   Level: developer
2473 
2474 .seealso: DMPlexGetSubpointMap(), DMPlexCreateSubmesh()
2475 */
2476 PetscErrorCode DMPlexCreateCohesiveSubmesh(DM dm, PetscBool hasLagrange, DM *subdm)
2477 {
2478   PetscInt       dim, depth;
2479   PetscErrorCode ierr;
2480 
2481   PetscFunctionBegin;
2482   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2483   PetscValidPointer(subdm, 3);
2484   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
2485   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
2486   ierr = DMCreate(PetscObjectComm((PetscObject)dm), subdm);CHKERRQ(ierr);
2487   ierr = DMSetType(*subdm, DMPLEX);CHKERRQ(ierr);
2488   ierr = DMPlexSetDimension(*subdm, dim-1);CHKERRQ(ierr);
2489   if (depth == dim) {
2490     ierr = DMPlexCreateCohesiveSubmesh_Interpolated(dm, hasLagrange, *subdm);CHKERRQ(ierr);
2491   } else {
2492     ierr = DMPlexCreateCohesiveSubmesh_Uninterpolated(dm, hasLagrange, *subdm);CHKERRQ(ierr);
2493   }
2494   PetscFunctionReturn(0);
2495 }
2496 
2497 #undef __FUNCT__
2498 #define __FUNCT__ "DMPlexGetSubpointMap"
2499 /*@
2500   DMPlexGetSubpointMap - Returns a DMLabel with point dimension as values
2501 
2502   Input Parameter:
2503 . dm - The submesh DM
2504 
2505   Output Parameter:
2506 . subpointMap - The DMLabel of all the points from the original mesh in this submesh, or NULL if this is not a submesh
2507 
2508   Level: developer
2509 
2510 .seealso: DMPlexCreateSubmesh(), DMPlexCreateSubpointIS()
2511 @*/
2512 PetscErrorCode DMPlexGetSubpointMap(DM dm, DMLabel *subpointMap)
2513 {
2514   DM_Plex *mesh = (DM_Plex*) dm->data;
2515 
2516   PetscFunctionBegin;
2517   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2518   PetscValidPointer(subpointMap, 2);
2519   *subpointMap = mesh->subpointMap;
2520   PetscFunctionReturn(0);
2521 }
2522 
2523 #undef __FUNCT__
2524 #define __FUNCT__ "DMPlexSetSubpointMap"
2525 /* Note: Should normally not be called by the user, since it is set in DMPlexCreateSubmesh() */
2526 PetscErrorCode DMPlexSetSubpointMap(DM dm, DMLabel subpointMap)
2527 {
2528   DM_Plex       *mesh = (DM_Plex *) dm->data;
2529   DMLabel        tmp;
2530   PetscErrorCode ierr;
2531 
2532   PetscFunctionBegin;
2533   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2534   tmp  = mesh->subpointMap;
2535   mesh->subpointMap = subpointMap;
2536   ++mesh->subpointMap->refct;
2537   ierr = DMLabelDestroy(&tmp);CHKERRQ(ierr);
2538   PetscFunctionReturn(0);
2539 }
2540 
2541 #undef __FUNCT__
2542 #define __FUNCT__ "DMPlexCreateSubpointIS"
2543 /*@
2544   DMPlexCreateSubpointIS - Creates an IS covering the entire subdm chart with the original points as data
2545 
2546   Input Parameter:
2547 . dm - The submesh DM
2548 
2549   Output Parameter:
2550 . subpointIS - The IS of all the points from the original mesh in this submesh, or NULL if this is not a submesh
2551 
2552   Note: This is IS is guaranteed to be sorted by the construction of the submesh
2553 
2554   Level: developer
2555 
2556 .seealso: DMPlexCreateSubmesh(), DMPlexGetSubpointMap()
2557 @*/
2558 PetscErrorCode DMPlexCreateSubpointIS(DM dm, IS *subpointIS)
2559 {
2560   MPI_Comm        comm;
2561   DMLabel         subpointMap;
2562   IS              is;
2563   const PetscInt *opoints;
2564   PetscInt       *points, *depths;
2565   PetscInt        depth, depStart, depEnd, d, pStart, pEnd, p, n, off;
2566   PetscErrorCode  ierr;
2567 
2568   PetscFunctionBegin;
2569   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2570   PetscValidPointer(subpointIS, 2);
2571   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
2572   *subpointIS = NULL;
2573   ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr);
2574   if (subpointMap) {
2575     ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
2576     ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
2577     if (pStart) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Submeshes must start the point numbering at 0, not %d", pStart);
2578     ierr = DMGetWorkArray(dm, depth+1, PETSC_INT, &depths);CHKERRQ(ierr);
2579     depths[0] = depth;
2580     depths[1] = 0;
2581     for(d = 2; d <= depth; ++d) {depths[d] = depth+1 - d;}
2582     ierr = PetscMalloc(pEnd * sizeof(PetscInt), &points);CHKERRQ(ierr);
2583     for(d = 0, off = 0; d <= depth; ++d) {
2584       const PetscInt dep = depths[d];
2585 
2586       ierr = DMPlexGetDepthStratum(dm, dep, &depStart, &depEnd);CHKERRQ(ierr);
2587       ierr = DMLabelGetStratumSize(subpointMap, dep, &n);CHKERRQ(ierr);
2588       if (((d < 2) && (depth > 1)) || (d == 1)) { /* Only check vertices and cells for now since the map is broken for others */
2589         if (n != depEnd-depStart) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "The number of mapped submesh points %d at depth %d should be %d", n, dep, depEnd-depStart);
2590       } else {
2591         if (!n) {
2592           if (d == 0) {
2593             /* Missing cells */
2594             for(p = 0; p < depEnd-depStart; ++p, ++off) points[off] = -1;
2595           } else {
2596             /* Missing faces */
2597             for(p = 0; p < depEnd-depStart; ++p, ++off) points[off] = PETSC_MAX_INT;
2598           }
2599         }
2600       }
2601       if (n) {
2602         ierr = DMLabelGetStratumIS(subpointMap, dep, &is);CHKERRQ(ierr);
2603         ierr = ISGetIndices(is, &opoints);CHKERRQ(ierr);
2604         for(p = 0; p < n; ++p, ++off) points[off] = opoints[p];
2605         ierr = ISRestoreIndices(is, &opoints);CHKERRQ(ierr);
2606         ierr = ISDestroy(&is);CHKERRQ(ierr);
2607       }
2608     }
2609     ierr = DMRestoreWorkArray(dm, depth+1, PETSC_INT, &depths);CHKERRQ(ierr);
2610     if (off != pEnd) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "The number of mapped submesh points %d should be %d", off, pEnd);
2611     ierr = ISCreateGeneral(comm, pEnd, points, PETSC_OWN_POINTER, subpointIS);CHKERRQ(ierr);
2612   }
2613   PetscFunctionReturn(0);
2614 }
2615