xref: /petsc/src/dm/impls/plex/plexsubmesh.c (revision f17f3550ad786ae781113eda50b18aadf6f31384)
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, *supportNew;
544   PetscInt        shift = 100, depth = 0, dep, dim, d, numSP = 0, sp, maxConeSize, maxSupportSize, numLabels, 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   /* Count split points and add cohesive cells */
551   if (label) {
552     ierr = DMLabelGetValueIS(label, &valueIS);CHKERRQ(ierr);
553     ierr = ISGetLocalSize(valueIS, &numSP);CHKERRQ(ierr);
554     ierr = ISGetIndices(valueIS, &values);CHKERRQ(ierr);
555   }
556   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
557   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
558   ierr = PetscMalloc5(depth+1,PetscInt,&depthShift,depth+1,PetscInt,&depthOffset,depth+1,PetscInt,&pMaxNew,maxConeSize*3,PetscInt,&coneNew,maxSupportSize,PetscInt,&supportNew);CHKERRQ(ierr);
559   ierr = PetscMalloc3(depth+1,IS,&pointIS,depth+1,PetscInt,&numSplitPoints,depth+1,const PetscInt*,&splitPoints);CHKERRQ(ierr);
560   ierr = PetscMemzero(depthShift, (depth+1) * sizeof(PetscInt));CHKERRQ(ierr);
561   for (d = 0; d <= depth; ++d) {
562     ierr              = DMPlexGetDepthStratum(dm, d, NULL, &pMaxNew[d]);CHKERRQ(ierr);
563     numSplitPoints[d] = 0;
564     splitPoints[d]    = NULL;
565     pointIS[d]        = NULL;
566   }
567   for (sp = 0; sp < numSP; ++sp) {
568     const PetscInt dep = values[sp];
569 
570     if ((dep < 0) || (dep > depth)) continue;
571     ierr = DMLabelGetStratumSize(label, dep, &depthShift[dep]);CHKERRQ(ierr);
572     ierr = DMLabelGetStratumIS(label, dep, &pointIS[dep]);CHKERRQ(ierr);
573     if (pointIS[dep]) {
574       ierr = ISGetLocalSize(pointIS[dep], &numSplitPoints[dep]);CHKERRQ(ierr);
575       ierr = ISGetIndices(pointIS[dep], &splitPoints[dep]);CHKERRQ(ierr);
576     }
577   }
578   if (depth >= 0) {
579     /* Calculate number of additional points */
580     depthShift[depth] = depthShift[depth-1]; /* There is a cohesive cell for every split face   */
581     depthShift[1]    += depthShift[0];       /* There is a cohesive edge for every split vertex */
582     /* Calculate hybrid bound for each dimension */
583     pMaxNew[0] += depthShift[depth];
584     if (depth > 1) pMaxNew[dim-1] += depthShift[depth] + depthShift[0];
585     if (depth > 2) pMaxNew[1]     += depthShift[depth] + depthShift[0] + depthShift[dim-1];
586 
587     /* Calculate point offset for each dimension */
588     depthOffset[depth] = 0;
589     depthOffset[0]     = depthOffset[depth] + depthShift[depth];
590     if (depth > 1) depthOffset[dim-1] = depthOffset[0]     + depthShift[0];
591     if (depth > 2) depthOffset[1]     = depthOffset[dim-1] + depthShift[dim-1];
592   }
593   ierr = DMPlexShiftSizes_Internal(dm, depthShift, sdm);CHKERRQ(ierr);
594   /* Step 3: Set cone/support sizes for new points */
595   for (dep = 0; dep <= depth; ++dep) {
596     for (p = 0; p < numSplitPoints[dep]; ++p) {
597       const PetscInt  oldp   = splitPoints[dep][p];
598       const PetscInt  newp   = depthOffset[dep] + oldp;
599       const PetscInt  splitp = pMaxNew[dep] + p;
600       const PetscInt *support;
601       PetscInt        coneSize, supportSize, q, e;
602 
603       ierr = DMPlexGetConeSize(dm, oldp, &coneSize);CHKERRQ(ierr);
604       ierr = DMPlexSetConeSize(sdm, splitp, coneSize);CHKERRQ(ierr);
605       ierr = DMPlexGetSupportSize(dm, oldp, &supportSize);CHKERRQ(ierr);
606       ierr = DMPlexSetSupportSize(sdm, splitp, supportSize);CHKERRQ(ierr);
607       if (dep == depth-1) {
608         const PetscInt ccell = pMaxNew[depth] + p;
609         /* Add cohesive cells, they are prisms */
610         ierr = DMPlexSetConeSize(sdm, ccell, 2 + coneSize);CHKERRQ(ierr);
611       } else if (dep == 0) {
612         const PetscInt cedge = pMaxNew[1] + (depthShift[1] - depthShift[0]) + p;
613 
614         ierr = DMPlexGetSupport(dm, oldp, &support);CHKERRQ(ierr);
615         /* Split old vertex: Edges in old split faces and new cohesive edge */
616         for (e = 0, q = 0; e < supportSize; ++e) {
617           PetscInt val;
618 
619           ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr);
620           if ((val == 1) || (val == (shift + 1))) ++q;
621         }
622         ierr = DMPlexSetSupportSize(sdm, newp, q+1);CHKERRQ(ierr);
623         /* Split new vertex: Edges in new split faces and new cohesive edge */
624         for (e = 0, q = 0; e < supportSize; ++e) {
625           PetscInt val;
626 
627           ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr);
628           if ((val == 1) || (val == -(shift + 1))) ++q;
629         }
630         ierr = DMPlexSetSupportSize(sdm, splitp, q+1);CHKERRQ(ierr);
631         /* Add cohesive edges */
632         ierr = DMPlexSetConeSize(sdm, cedge, 2);CHKERRQ(ierr);
633         /* Punt for now on support, you loop over closure, extract faces, check which ones are in the label */
634       } else if (dep == dim-2) {
635         ierr = DMPlexGetSupport(dm, oldp, &support);CHKERRQ(ierr);
636         /* Split old edge: Faces in positive side cells and old split faces */
637         for (e = 0, q = 0; e < supportSize; ++e) {
638           PetscInt val;
639 
640           ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr);
641           if ((val == dim-1) || (val == (shift + dim-1))) ++q;
642         }
643         ierr = DMPlexSetSupportSize(sdm, newp, q);CHKERRQ(ierr);
644         /* Split new edge: Faces in negative side cells and new split faces */
645         for (e = 0, q = 0; e < supportSize; ++e) {
646           PetscInt val;
647 
648           ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr);
649           if ((val == dim-1) || (val == -(shift + dim-1))) ++q;
650         }
651         ierr = DMPlexSetSupportSize(sdm, splitp, q);CHKERRQ(ierr);
652       }
653     }
654   }
655   /* Step 4: Setup split DM */
656   ierr = DMSetUp(sdm);CHKERRQ(ierr);
657   ierr = DMPlexShiftPoints_Internal(dm, depthShift, sdm);CHKERRQ(ierr);
658   /* Step 6: Set cones and supports for new points */
659   for (dep = 0; dep <= depth; ++dep) {
660     for (p = 0; p < numSplitPoints[dep]; ++p) {
661       const PetscInt  oldp   = splitPoints[dep][p];
662       const PetscInt  newp   = depthOffset[dep] + oldp;
663       const PetscInt  splitp = pMaxNew[dep] + p;
664       const PetscInt *cone, *support, *ornt;
665       PetscInt        coneSize, supportSize, q, v, e, s;
666 
667       ierr = DMPlexGetConeSize(dm, oldp, &coneSize);CHKERRQ(ierr);
668       ierr = DMPlexGetCone(dm, oldp, &cone);CHKERRQ(ierr);
669       ierr = DMPlexGetConeOrientation(dm, oldp, &ornt);CHKERRQ(ierr);
670       ierr = DMPlexGetSupportSize(dm, oldp, &supportSize);CHKERRQ(ierr);
671       ierr = DMPlexGetSupport(dm, oldp, &support);CHKERRQ(ierr);
672       if (dep == depth-1) {
673         const PetscInt  ccell = pMaxNew[depth] + p;
674         const PetscInt *supportF;
675 
676         /* Split face:       copy in old face to new face to start */
677         ierr = DMPlexGetSupport(sdm, newp,  &supportF);CHKERRQ(ierr);
678         ierr = DMPlexSetSupport(sdm, splitp, supportF);CHKERRQ(ierr);
679         /* Split old face:   old vertices/edges in cone so no change */
680         /* Split new face:   new vertices/edges in cone */
681         for (q = 0; q < coneSize; ++q) {
682           ierr = PetscFindInt(cone[q], numSplitPoints[dim-2], splitPoints[dim-2], &v);CHKERRQ(ierr);
683 
684           coneNew[2+q] = pMaxNew[dim-2] + v;
685         }
686         ierr = DMPlexSetCone(sdm, splitp, &coneNew[2]);CHKERRQ(ierr);
687         ierr = DMPlexSetConeOrientation(sdm, splitp, ornt);CHKERRQ(ierr);
688         /* Cohesive cell:    Old and new split face, then new cohesive edges */
689         coneNew[0] = newp;
690         coneNew[1] = splitp;
691         for (q = 0; q < coneSize; ++q) {
692           coneNew[2+q] = (pMaxNew[1] - pMaxNew[dim-2]) + (depthShift[1] - depthShift[0]) + coneNew[2+q];
693         }
694         ierr = DMPlexSetCone(sdm, ccell, coneNew);CHKERRQ(ierr);
695 
696 
697         for (s = 0; s < supportSize; ++s) {
698           PetscInt val;
699 
700           ierr = DMLabelGetValue(label, support[s], &val);CHKERRQ(ierr);
701           if (val < 0) {
702             /* Split old face:   Replace negative side cell with cohesive cell */
703             ierr = DMPlexInsertSupport(sdm, newp, s, ccell);CHKERRQ(ierr);
704           } else {
705             /* Split new face:   Replace positive side cell with cohesive cell */
706             ierr = DMPlexInsertSupport(sdm, splitp, s, ccell);CHKERRQ(ierr);
707           }
708         }
709       } else if (dep == 0) {
710         const PetscInt cedge = pMaxNew[1] + (depthShift[1] - depthShift[0]) + p;
711 
712         /* Split old vertex: Edges in old split faces and new cohesive edge */
713         for (e = 0, q = 0; e < supportSize; ++e) {
714           PetscInt val;
715 
716           ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr);
717           if ((val == 1) || (val == (shift + 1))) {
718             supportNew[q++] = depthOffset[1] + support[e];
719           }
720         }
721         supportNew[q] = cedge;
722 
723         ierr = DMPlexSetSupport(sdm, newp, supportNew);CHKERRQ(ierr);
724         /* Split new vertex: Edges in new split faces and new cohesive edge */
725         for (e = 0, q = 0; e < supportSize; ++e) {
726           PetscInt val, edge;
727 
728           ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr);
729           if (val == 1) {
730             ierr = PetscFindInt(support[e], numSplitPoints[1], splitPoints[1], &edge);CHKERRQ(ierr);
731             if (edge < 0) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Edge %d is not a split edge", support[e]);
732             supportNew[q++] = pMaxNew[1] + edge;
733           } else if (val == -(shift + 1)) {
734             supportNew[q++] = depthOffset[1] + support[e];
735           }
736         }
737         supportNew[q] = cedge;
738         ierr          = DMPlexSetSupport(sdm, splitp, supportNew);CHKERRQ(ierr);
739         /* Cohesive edge:    Old and new split vertex, punting on support */
740         coneNew[0] = newp;
741         coneNew[1] = splitp;
742         ierr       = DMPlexSetCone(sdm, cedge, coneNew);CHKERRQ(ierr);
743       } else if (dep == dim-2) {
744         /* Split old edge:   old vertices in cone so no change */
745         /* Split new edge:   new vertices in cone */
746         for (q = 0; q < coneSize; ++q) {
747           ierr = PetscFindInt(cone[q], numSplitPoints[dim-3], splitPoints[dim-3], &v);CHKERRQ(ierr);
748 
749           coneNew[q] = pMaxNew[dim-3] + v;
750         }
751         ierr = DMPlexSetCone(sdm, splitp, coneNew);CHKERRQ(ierr);
752         /* Split old edge: Faces in positive side cells and old split faces */
753         for (e = 0, q = 0; e < supportSize; ++e) {
754           PetscInt val;
755 
756           ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr);
757           if ((val == dim-1) || (val == (shift + dim-1))) {
758             supportNew[q++] = depthOffset[dim-1] + support[e];
759           }
760         }
761         ierr = DMPlexSetSupport(sdm, newp, supportNew);CHKERRQ(ierr);
762         /* Split new edge: Faces in negative side cells and new split faces */
763         for (e = 0, q = 0; e < supportSize; ++e) {
764           PetscInt val, face;
765 
766           ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr);
767           if (val == dim-1) {
768             ierr = PetscFindInt(support[e], numSplitPoints[dim-1], splitPoints[dim-1], &face);CHKERRQ(ierr);
769             if (face < 0) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Face %d is not a split face", support[e]);
770             supportNew[q++] = pMaxNew[dim-1] + face;
771           } else if (val == -(shift + dim-1)) {
772             supportNew[q++] = depthOffset[dim-1] + support[e];
773           }
774         }
775         ierr = DMPlexSetSupport(sdm, splitp, supportNew);CHKERRQ(ierr);
776       }
777     }
778   }
779   /* Step 6b: Replace split points in negative side cones */
780   for (sp = 0; sp < numSP; ++sp) {
781     PetscInt        dep = values[sp];
782     IS              pIS;
783     PetscInt        numPoints;
784     const PetscInt *points;
785 
786     if (dep >= 0) continue;
787     ierr = DMLabelGetStratumIS(label, dep, &pIS);CHKERRQ(ierr);
788     if (!pIS) continue;
789     dep  = -dep - shift;
790     ierr = ISGetLocalSize(pIS, &numPoints);CHKERRQ(ierr);
791     ierr = ISGetIndices(pIS, &points);CHKERRQ(ierr);
792     for (p = 0; p < numPoints; ++p) {
793       const PetscInt  oldp = points[p];
794       const PetscInt  newp = depthOffset[dep] + oldp;
795       const PetscInt *cone;
796       PetscInt        coneSize, c;
797       PetscBool       replaced = PETSC_FALSE;
798 
799       /* Negative edge: replace split vertex */
800       /* Negative cell: replace split face */
801       ierr = DMPlexGetConeSize(sdm, newp, &coneSize);CHKERRQ(ierr);
802       ierr = DMPlexGetCone(sdm, newp, &cone);CHKERRQ(ierr);
803       for (c = 0; c < coneSize; ++c) {
804         const PetscInt coldp = cone[c] - depthOffset[dep-1];
805         PetscInt       csplitp, cp, val;
806 
807         ierr = DMLabelGetValue(label, coldp, &val);CHKERRQ(ierr);
808         if (val == dep-1) {
809           ierr = PetscFindInt(coldp, numSplitPoints[dep-1], splitPoints[dep-1], &cp);CHKERRQ(ierr);
810           if (cp < 0) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "Point %d is not a split point of dimension %d", oldp, dep-1);
811           csplitp  = pMaxNew[dep-1] + cp;
812           ierr     = DMPlexInsertCone(sdm, newp, c, csplitp);CHKERRQ(ierr);
813           replaced = PETSC_TRUE;
814         }
815       }
816       if (!replaced) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "The cone of point %d does not contain split points", oldp);
817     }
818     ierr = ISRestoreIndices(pIS, &points);CHKERRQ(ierr);
819     ierr = ISDestroy(&pIS);CHKERRQ(ierr);
820   }
821   /* Step 7: Stratify */
822   ierr = DMPlexStratify(sdm);CHKERRQ(ierr);
823   /* Step 8: Coordinates */
824   ierr = DMPlexShiftCoordinates_Internal(dm, depthShift, sdm);CHKERRQ(ierr);
825   ierr = DMPlexGetCoordinateSection(sdm, &coordSection);CHKERRQ(ierr);
826   ierr = DMGetCoordinatesLocal(sdm, &coordinates);CHKERRQ(ierr);
827   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
828   for (v = 0; v < (numSplitPoints ? numSplitPoints[0] : 0); ++v) {
829     const PetscInt newp   = depthOffset[0] + splitPoints[0][v];
830     const PetscInt splitp = pMaxNew[0] + v;
831     PetscInt       dof, off, soff, d;
832 
833     ierr = PetscSectionGetDof(coordSection, newp, &dof);CHKERRQ(ierr);
834     ierr = PetscSectionGetOffset(coordSection, newp, &off);CHKERRQ(ierr);
835     ierr = PetscSectionGetOffset(coordSection, splitp, &soff);CHKERRQ(ierr);
836     for (d = 0; d < dof; ++d) coords[soff+d] = coords[off+d];
837   }
838   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
839   /* Step 9: SF, if I can figure this out we can split the mesh in parallel */
840   ierr = DMPlexShiftSF_Internal(dm, depthShift, sdm);CHKERRQ(ierr);
841   /* Step 10: Labels */
842   ierr = DMPlexShiftLabels_Internal(dm, depthShift, sdm);CHKERRQ(ierr);
843   ierr = DMPlexGetNumLabels(sdm, &numLabels);CHKERRQ(ierr);
844   for (dep = 0; dep <= depth; ++dep) {
845     for (p = 0; p < numSplitPoints[dep]; ++p) {
846       const PetscInt newp   = depthOffset[dep] + splitPoints[dep][p];
847       const PetscInt splitp = pMaxNew[dep] + p;
848       PetscInt       l;
849 
850       for (l = 0; l < numLabels; ++l) {
851         DMLabel     mlabel;
852         const char *lname;
853         PetscInt    val;
854         PetscBool   isDepth;
855 
856         ierr = DMPlexGetLabelName(sdm, l, &lname);CHKERRQ(ierr);
857         ierr = PetscStrcmp(lname, "depth", &isDepth);CHKERRQ(ierr);
858         if (isDepth) continue;
859         ierr = DMPlexGetLabel(sdm, lname, &mlabel);CHKERRQ(ierr);
860         ierr = DMLabelGetValue(mlabel, newp, &val);CHKERRQ(ierr);
861         if (val >= 0) {
862           ierr = DMLabelSetValue(mlabel, splitp, val);CHKERRQ(ierr);
863 #if 0
864           /* Do not put cohesive edges into the label */
865           if (dep == 0) {
866             const PetscInt cedge = pMaxNew[1] + (depthShift[1] - depthShift[0]) + p;
867             ierr = DMLabelSetValue(mlabel, cedge, val);CHKERRQ(ierr);
868           }
869 #endif
870         }
871       }
872     }
873   }
874   for (sp = 0; sp < numSP; ++sp) {
875     const PetscInt dep = values[sp];
876 
877     if ((dep < 0) || (dep > depth)) continue;
878     if (pointIS[dep]) {ierr = ISRestoreIndices(pointIS[dep], &splitPoints[dep]);CHKERRQ(ierr);}
879     ierr = ISDestroy(&pointIS[dep]);CHKERRQ(ierr);
880   }
881   if (label) {
882     ierr = ISRestoreIndices(valueIS, &values);CHKERRQ(ierr);
883     ierr = ISDestroy(&valueIS);CHKERRQ(ierr);
884   }
885   ierr = DMPlexGetChart(sdm, NULL, &pEnd);CHKERRQ(ierr);
886   pMaxNew[0] += depthShift[0]; /* Account for shadow vertices */
887   if (depth > 1) pMaxNew[1] = pEnd - depthShift[0]; /* There is a hybrid edge for every shadow vertex */
888   if (depth > 2) pMaxNew[2] = -1; /* There are no hybrid faces */
889   ierr = DMPlexSetHybridBounds(sdm, depth >= 0 ? pMaxNew[depth] : PETSC_DETERMINE, depth>1 ? pMaxNew[depth-1] : PETSC_DETERMINE, depth>2 ? pMaxNew[1] : PETSC_DETERMINE, pMaxNew[0]);CHKERRQ(ierr);
890   ierr = PetscFree5(depthShift, depthOffset, pMaxNew, coneNew, supportNew);CHKERRQ(ierr);
891   ierr = PetscFree3(pointIS, numSplitPoints, splitPoints);CHKERRQ(ierr);
892   PetscFunctionReturn(0);
893 }
894 
895 #undef __FUNCT__
896 #define __FUNCT__ "DMPlexConstructCohesiveCells"
897 /*@C
898   DMPlexConstructCohesiveCells - Construct cohesive cells which split the face along an internal interface
899 
900   Collective on dm
901 
902   Input Parameters:
903 + dm - The original DM
904 - labelName - The label specifying the boundary faces (this could be auto-generated)
905 
906   Output Parameters:
907 - dmSplit - The new DM
908 
909   Level: developer
910 
911 .seealso: DMCreate(), DMPlexLabelCohesiveComplete()
912 @*/
913 PetscErrorCode DMPlexConstructCohesiveCells(DM dm, DMLabel label, DM *dmSplit)
914 {
915   DM             sdm;
916   PetscInt       dim;
917   PetscErrorCode ierr;
918 
919   PetscFunctionBegin;
920   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
921   PetscValidPointer(dmSplit, 4);
922   ierr = DMCreate(PetscObjectComm((PetscObject)dm), &sdm);CHKERRQ(ierr);
923   ierr = DMSetType(sdm, DMPLEX);CHKERRQ(ierr);
924   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
925   ierr = DMPlexSetDimension(sdm, dim);CHKERRQ(ierr);
926   switch (dim) {
927   case 2:
928   case 3:
929     ierr = DMPlexConstructCohesiveCells_Internal(dm, label, sdm);CHKERRQ(ierr);
930     break;
931   default:
932     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct cohesive cells for dimension %d", dim);
933   }
934   *dmSplit = sdm;
935   PetscFunctionReturn(0);
936 }
937 
938 #undef __FUNCT__
939 #define __FUNCT__ "DMPlexLabelCohesiveComplete"
940 /*@
941   DMPlexLabelCohesiveComplete - Starting with a label marking vertices on an internal surface, we add all other mesh pieces
942   to complete the surface
943 
944   Input Parameters:
945 + dm - The DM
946 - label - A DMLabel marking the surface vertices
947 
948   Output Parameter:
949 . label - A DMLabel marking all surface points
950 
951   Level: developer
952 
953 .seealso: DMPlexConstructCohesiveCells(), DMPlexLabelComplete()
954 @*/
955 PetscErrorCode DMPlexLabelCohesiveComplete(DM dm, DMLabel label)
956 {
957   IS              dimIS;
958   const PetscInt *points;
959   PetscInt        shift = 100, dim, dep, cStart, cEnd, numPoints, p, val;
960   PetscErrorCode  ierr;
961 
962   PetscFunctionBegin;
963   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
964   /* Cell orientation for face gives the side of the fault */
965   ierr = DMLabelGetStratumIS(label, dim-1, &dimIS);CHKERRQ(ierr);
966   if (!dimIS) PetscFunctionReturn(0);
967   ierr = ISGetLocalSize(dimIS, &numPoints);CHKERRQ(ierr);
968   ierr = ISGetIndices(dimIS, &points);CHKERRQ(ierr);
969   for (p = 0; p < numPoints; ++p) {
970     const PetscInt *support;
971     PetscInt        supportSize, s;
972 
973     ierr = DMPlexGetSupportSize(dm, points[p], &supportSize);CHKERRQ(ierr);
974     if (supportSize != 2) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Split face %d has %d != 2 supports", points[p], supportSize);
975     ierr = DMPlexGetSupport(dm, points[p], &support);CHKERRQ(ierr);
976     for (s = 0; s < supportSize; ++s) {
977       const PetscInt *cone, *ornt;
978       PetscInt        coneSize, c;
979       PetscBool       pos = PETSC_TRUE;
980 
981       ierr = DMPlexGetConeSize(dm, support[s], &coneSize);CHKERRQ(ierr);
982       ierr = DMPlexGetCone(dm, support[s], &cone);CHKERRQ(ierr);
983       ierr = DMPlexGetConeOrientation(dm, support[s], &ornt);CHKERRQ(ierr);
984       for (c = 0; c < coneSize; ++c) {
985         if (cone[c] == points[p]) {
986           if (ornt[c] >= 0) {
987             ierr = DMLabelSetValue(label, support[s],   shift+dim);CHKERRQ(ierr);
988           } else {
989             ierr = DMLabelSetValue(label, support[s], -(shift+dim));CHKERRQ(ierr);
990             pos  = PETSC_FALSE;
991           }
992           break;
993         }
994       }
995       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]);
996       /* Put faces touching the fault in the label */
997       for (c = 0; c < coneSize; ++c) {
998         const PetscInt point = cone[c];
999 
1000         ierr = DMLabelGetValue(label, point, &val);CHKERRQ(ierr);
1001         if (val == -1) {
1002           PetscInt *closure = NULL;
1003           PetscInt  closureSize, cl;
1004 
1005           ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1006           for (cl = 0; cl < closureSize*2; cl += 2) {
1007             const PetscInt clp = closure[cl];
1008 
1009             ierr = DMLabelGetValue(label, clp, &val);CHKERRQ(ierr);
1010             if ((val >= 0) && (val < dim-1)) {
1011               ierr = DMLabelSetValue(label, point, pos == PETSC_TRUE ? shift+dim-1 : -(shift+dim-1));CHKERRQ(ierr);
1012               break;
1013             }
1014           }
1015           ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1016         }
1017       }
1018     }
1019   }
1020   ierr = ISRestoreIndices(dimIS, &points);CHKERRQ(ierr);
1021   ierr = ISDestroy(&dimIS);CHKERRQ(ierr);
1022   /* Search for other cells/faces/edges connected to the fault by a vertex */
1023   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1024   ierr = DMLabelGetStratumIS(label, 0, &dimIS);CHKERRQ(ierr);
1025   if (!dimIS) PetscFunctionReturn(0);
1026   ierr = ISGetLocalSize(dimIS, &numPoints);CHKERRQ(ierr);
1027   ierr = ISGetIndices(dimIS, &points);CHKERRQ(ierr);
1028   for (p = 0; p < numPoints; ++p) {
1029     PetscInt *star = NULL;
1030     PetscInt  starSize, s;
1031     PetscInt  again = 1;  /* 0: Finished 1: Keep iterating after a change 2: No change */
1032 
1033     /* First mark cells connected to the fault */
1034     ierr = DMPlexGetTransitiveClosure(dm, points[p], PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
1035     while (again) {
1036       if (again > 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Could not classify all cells connected to the fault");
1037       again = 0;
1038       for (s = 0; s < starSize*2; s += 2) {
1039         const PetscInt  point = star[s];
1040         const PetscInt *cone;
1041         PetscInt        coneSize, c;
1042 
1043         if ((point < cStart) || (point >= cEnd)) continue;
1044         ierr = DMLabelGetValue(label, point, &val);CHKERRQ(ierr);
1045         if (val != -1) continue;
1046         again = 2;
1047         ierr  = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
1048         ierr  = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
1049         for (c = 0; c < coneSize; ++c) {
1050           ierr = DMLabelGetValue(label, cone[c], &val);CHKERRQ(ierr);
1051           if (val != -1) {
1052             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);
1053             if (val > 0) {
1054               ierr = DMLabelSetValue(label, point,   shift+dim);CHKERRQ(ierr);
1055             } else {
1056               ierr = DMLabelSetValue(label, point, -(shift+dim));CHKERRQ(ierr);
1057             }
1058             again = 1;
1059             break;
1060           }
1061         }
1062       }
1063     }
1064     /* Classify the rest by cell membership */
1065     for (s = 0; s < starSize*2; s += 2) {
1066       const PetscInt point = star[s];
1067 
1068       ierr = DMLabelGetValue(label, point, &val);CHKERRQ(ierr);
1069       if (val == -1) {
1070         PetscInt  *sstar = NULL;
1071         PetscInt   sstarSize, ss;
1072         PetscBool  marked = PETSC_FALSE;
1073 
1074         ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &sstarSize, &sstar);CHKERRQ(ierr);
1075         for (ss = 0; ss < sstarSize*2; ss += 2) {
1076           const PetscInt spoint = sstar[ss];
1077 
1078           if ((spoint < cStart) || (spoint >= cEnd)) continue;
1079           ierr = DMLabelGetValue(label, spoint, &val);CHKERRQ(ierr);
1080           if (val == -1) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Cell %d in star of %d does not have a valid label", spoint, point);
1081           ierr = DMPlexGetLabelValue(dm, "depth", point, &dep);CHKERRQ(ierr);
1082           if (val > 0) {
1083             ierr = DMLabelSetValue(label, point,   shift+dep);CHKERRQ(ierr);
1084           } else {
1085             ierr = DMLabelSetValue(label, point, -(shift+dep));CHKERRQ(ierr);
1086           }
1087           marked = PETSC_TRUE;
1088           break;
1089         }
1090         ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, &sstarSize, &sstar);CHKERRQ(ierr);
1091         if (!marked) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Point %d could not be classified", point);
1092       }
1093     }
1094     ierr = DMPlexRestoreTransitiveClosure(dm, points[p], PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
1095   }
1096   ierr = ISRestoreIndices(dimIS, &points);CHKERRQ(ierr);
1097   ierr = ISDestroy(&dimIS);CHKERRQ(ierr);
1098   PetscFunctionReturn(0);
1099 }
1100 
1101 #undef __FUNCT__
1102 #define __FUNCT__ "DMPlexMarkSubmesh_Uninterpolated"
1103 /* Here we need the explicit assumption that:
1104 
1105      For any marked cell, the marked vertices constitute a single face
1106 */
1107 static PetscErrorCode DMPlexMarkSubmesh_Uninterpolated(DM dm, DMLabel vertexLabel, PetscInt value, DMLabel subpointMap, PetscInt *numFaces, PetscInt *nFV, DM subdm)
1108 {
1109   IS               subvertexIS;
1110   const PetscInt  *subvertices;
1111   PetscInt        *pStart, *pEnd, *pMax, pSize;
1112   PetscInt         depth, dim, d, numSubVerticesInitial = 0, v;
1113   PetscErrorCode   ierr;
1114 
1115   PetscFunctionBegin;
1116   *numFaces = 0;
1117   *nFV      = 0;
1118   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1119   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1120   pSize = PetscMax(depth, dim) + 1;
1121   ierr = PetscMalloc3(pSize,PetscInt,&pStart,pSize,PetscInt,&pEnd,pSize,PetscInt,&pMax);CHKERRQ(ierr);
1122   ierr = DMPlexGetHybridBounds(dm, &pMax[depth], depth>1 ? &pMax[depth-1] : NULL, depth > 2 ? &pMax[1] : NULL, &pMax[0]);CHKERRQ(ierr);
1123   for (d = 0; d <= depth; ++d) {
1124     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr);
1125     if (pMax[d] >= 0) pEnd[d] = PetscMin(pEnd[d], pMax[d]);
1126   }
1127   /* Loop over initial vertices and mark all faces in the collective star() */
1128   ierr = DMLabelGetStratumIS(vertexLabel, value, &subvertexIS);CHKERRQ(ierr);
1129   if (subvertexIS) {
1130     ierr = ISGetSize(subvertexIS, &numSubVerticesInitial);CHKERRQ(ierr);
1131     ierr = ISGetIndices(subvertexIS, &subvertices);CHKERRQ(ierr);
1132   }
1133   for (v = 0; v < numSubVerticesInitial; ++v) {
1134     const PetscInt vertex = subvertices[v];
1135     PetscInt      *star   = NULL;
1136     PetscInt       starSize, s, numCells = 0, c;
1137 
1138     ierr = DMPlexGetTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
1139     for (s = 0; s < starSize*2; s += 2) {
1140       const PetscInt point = star[s];
1141       if ((point >= pStart[depth]) && (point < pEnd[depth])) star[numCells++] = point;
1142     }
1143     for (c = 0; c < numCells; ++c) {
1144       const PetscInt cell    = star[c];
1145       PetscInt      *closure = NULL;
1146       PetscInt       closureSize, cl;
1147       PetscInt       cellLoc, numCorners = 0, faceSize = 0;
1148 
1149       ierr = DMLabelGetValue(subpointMap, cell, &cellLoc);CHKERRQ(ierr);
1150       if (cellLoc == 2) continue;
1151       if (cellLoc >= 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Cell %d has dimension %d in the surface label", cell, cellLoc);
1152       ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1153       for (cl = 0; cl < closureSize*2; cl += 2) {
1154         const PetscInt point = closure[cl];
1155         PetscInt       vertexLoc;
1156 
1157         if ((point >= pStart[0]) && (point < pEnd[0])) {
1158           ++numCorners;
1159           ierr = DMLabelGetValue(vertexLabel, point, &vertexLoc);CHKERRQ(ierr);
1160           if (vertexLoc == value) closure[faceSize++] = point;
1161         }
1162       }
1163       if (!(*nFV)) {ierr = DMPlexGetNumFaceVertices(dm, dim, numCorners, nFV);CHKERRQ(ierr);}
1164       if (faceSize > *nFV) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid submesh: Too many vertices %d of an element on the surface", faceSize);
1165       if (faceSize == *nFV) {
1166         const PetscInt *cells = NULL;
1167         PetscInt        numCells, nc;
1168 
1169         ++(*numFaces);
1170         for (cl = 0; cl < faceSize; ++cl) {
1171           ierr = DMLabelSetValue(subpointMap, closure[cl], 0);CHKERRQ(ierr);
1172         }
1173         ierr = DMPlexGetJoin(dm, faceSize, closure, &numCells, &cells);CHKERRQ(ierr);
1174         for (nc = 0; nc < numCells; ++nc) {
1175           ierr = DMLabelSetValue(subpointMap, cells[nc], 2);CHKERRQ(ierr);
1176         }
1177         ierr = DMPlexRestoreJoin(dm, faceSize, closure, &numCells, &cells);CHKERRQ(ierr);
1178       }
1179       ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1180     }
1181     ierr = DMPlexRestoreTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
1182   }
1183   if (subvertexIS) {
1184     ierr = ISRestoreIndices(subvertexIS, &subvertices);CHKERRQ(ierr);
1185   }
1186   ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr);
1187   ierr = PetscFree3(pStart,pEnd,pMax);CHKERRQ(ierr);
1188   PetscFunctionReturn(0);
1189 }
1190 
1191 #undef __FUNCT__
1192 #define __FUNCT__ "DMPlexMarkSubmesh_Interpolated"
1193 static PetscErrorCode DMPlexMarkSubmesh_Interpolated(DM dm, DMLabel vertexLabel, PetscInt value, DMLabel subpointMap, DM subdm)
1194 {
1195   IS               subvertexIS;
1196   const PetscInt  *subvertices;
1197   PetscInt        *pStart, *pEnd, *pMax;
1198   PetscInt         dim, d, numSubVerticesInitial = 0, v;
1199   PetscErrorCode   ierr;
1200 
1201   PetscFunctionBegin;
1202   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1203   ierr = PetscMalloc3(dim+1,PetscInt,&pStart,dim+1,PetscInt,&pEnd,dim+1,PetscInt,&pMax);CHKERRQ(ierr);
1204   ierr = DMPlexGetHybridBounds(dm, &pMax[dim], dim>1 ? &pMax[dim-1] : NULL, dim > 2 ? &pMax[1] : NULL, &pMax[0]);CHKERRQ(ierr);
1205   for (d = 0; d <= dim; ++d) {
1206     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr);
1207     if (pMax[d] >= 0) pEnd[d] = PetscMin(pEnd[d], pMax[d]);
1208   }
1209   /* Loop over initial vertices and mark all faces in the collective star() */
1210   ierr = DMLabelGetStratumIS(vertexLabel, value, &subvertexIS);CHKERRQ(ierr);
1211   if (subvertexIS) {
1212     ierr = ISGetSize(subvertexIS, &numSubVerticesInitial);CHKERRQ(ierr);
1213     ierr = ISGetIndices(subvertexIS, &subvertices);CHKERRQ(ierr);
1214   }
1215   for (v = 0; v < numSubVerticesInitial; ++v) {
1216     const PetscInt vertex = subvertices[v];
1217     PetscInt      *star   = NULL;
1218     PetscInt       starSize, s, numFaces = 0, f;
1219 
1220     ierr = DMPlexGetTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
1221     for (s = 0; s < starSize*2; s += 2) {
1222       const PetscInt point = star[s];
1223       if ((point >= pStart[dim-1]) && (point < pEnd[dim-1])) star[numFaces++] = point;
1224     }
1225     for (f = 0; f < numFaces; ++f) {
1226       const PetscInt face    = star[f];
1227       PetscInt      *closure = NULL;
1228       PetscInt       closureSize, c;
1229       PetscInt       faceLoc;
1230 
1231       ierr = DMLabelGetValue(subpointMap, face, &faceLoc);CHKERRQ(ierr);
1232       if (faceLoc == dim-1) continue;
1233       if (faceLoc >= 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Face %d has dimension %d in the surface label", face, faceLoc);
1234       ierr = DMPlexGetTransitiveClosure(dm, face, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1235       for (c = 0; c < closureSize*2; c += 2) {
1236         const PetscInt point = closure[c];
1237         PetscInt       vertexLoc;
1238 
1239         if ((point >= pStart[0]) && (point < pEnd[0])) {
1240           ierr = DMLabelGetValue(vertexLabel, point, &vertexLoc);CHKERRQ(ierr);
1241           if (vertexLoc != value) break;
1242         }
1243       }
1244       if (c == closureSize*2) {
1245         const PetscInt *support;
1246         PetscInt        supportSize, s;
1247 
1248         for (c = 0; c < closureSize*2; c += 2) {
1249           const PetscInt point = closure[c];
1250 
1251           for (d = 0; d < dim; ++d) {
1252             if ((point >= pStart[d]) && (point < pEnd[d])) {
1253               ierr = DMLabelSetValue(subpointMap, point, d);CHKERRQ(ierr);
1254               break;
1255             }
1256           }
1257         }
1258         ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
1259         ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr);
1260         for (s = 0; s < supportSize; ++s) {
1261           ierr = DMLabelSetValue(subpointMap, support[s], dim);CHKERRQ(ierr);
1262         }
1263       }
1264       ierr = DMPlexRestoreTransitiveClosure(dm, face, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1265     }
1266     ierr = DMPlexRestoreTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
1267   }
1268   if (subvertexIS) {
1269     ierr = ISRestoreIndices(subvertexIS, &subvertices);CHKERRQ(ierr);
1270   }
1271   ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr);
1272   ierr = PetscFree3(pStart,pEnd,pMax);CHKERRQ(ierr);
1273   PetscFunctionReturn(0);
1274 }
1275 
1276 #undef __FUNCT__
1277 #define __FUNCT__ "DMPlexMarkCohesiveSubmesh_Uninterpolated"
1278 static PetscErrorCode DMPlexMarkCohesiveSubmesh_Uninterpolated(DM dm, PetscBool hasLagrange, const char labelname[], PetscInt value, DMLabel subpointMap, PetscInt *numFaces, PetscInt *nFV, PetscInt *subCells[], DM subdm)
1279 {
1280   DMLabel         label = NULL;
1281   const PetscInt *cone;
1282   PetscInt        dim, cMax, cEnd, c, p, coneSize;
1283   PetscErrorCode  ierr;
1284 
1285   PetscFunctionBegin;
1286   if (labelname) {ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr);}
1287   *numFaces = 0;
1288   *nFV = 0;
1289   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1290   ierr = DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);CHKERRQ(ierr);
1291   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
1292   if (cMax < 0) PetscFunctionReturn(0);
1293   if (label) {
1294     for (c = cMax; c < cEnd; ++c) {
1295       PetscInt val;
1296 
1297       ierr = DMLabelGetValue(label, c, &val);CHKERRQ(ierr);
1298       if (val == value) {
1299         ++(*numFaces);
1300         ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
1301       }
1302     }
1303   } else {
1304     *numFaces = cEnd - cMax;
1305     ierr = DMPlexGetConeSize(dm, cMax, &coneSize);CHKERRQ(ierr);
1306   }
1307   *nFV = hasLagrange ? coneSize/3 : coneSize/2;
1308   ierr = PetscMalloc(*numFaces *2 * sizeof(PetscInt), subCells);CHKERRQ(ierr);
1309   for (c = cMax; c < cEnd; ++c) {
1310     const PetscInt *cells;
1311     PetscInt        numCells;
1312 
1313     if (label) {
1314       PetscInt val;
1315 
1316       ierr = DMLabelGetValue(label, c, &val);CHKERRQ(ierr);
1317       if (val != value) continue;
1318     }
1319     ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr);
1320     for (p = 0; p < *nFV; ++p) {
1321       ierr = DMLabelSetValue(subpointMap, cone[p], 0);CHKERRQ(ierr);
1322     }
1323     /* Negative face */
1324     ierr = DMPlexGetJoin(dm, *nFV, cone, &numCells, &cells);CHKERRQ(ierr);
1325     if (numCells != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells");
1326     for (p = 0; p < numCells; ++p) {
1327       ierr = DMLabelSetValue(subpointMap, cells[p], 2);CHKERRQ(ierr);
1328       (*subCells)[(c-cMax)*2+p] = cells[p];
1329     }
1330     ierr = DMPlexRestoreJoin(dm, *nFV, cone, &numCells, &cells);CHKERRQ(ierr);
1331     /* Positive face is not included */
1332   }
1333   PetscFunctionReturn(0);
1334 }
1335 
1336 #undef __FUNCT__
1337 #define __FUNCT__ "DMPlexMarkCohesiveSubmesh_Interpolated"
1338 static PetscErrorCode DMPlexMarkCohesiveSubmesh_Interpolated(DM dm, PetscBool hasLagrange, const char labelname[], PetscInt value, DMLabel subpointMap, DM subdm)
1339 {
1340   DMLabel        label = NULL;
1341   PetscInt      *pStart, *pEnd;
1342   PetscInt       dim, cMax, cEnd, c, d;
1343   PetscErrorCode ierr;
1344 
1345   PetscFunctionBegin;
1346   if (labelname) {ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr);}
1347   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1348   ierr = DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);CHKERRQ(ierr);
1349   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
1350   if (cMax < 0) PetscFunctionReturn(0);
1351   ierr = PetscMalloc2(dim+1,PetscInt,&pStart,dim+1,PetscInt,&pEnd);CHKERRQ(ierr);
1352   for (d = 0; d <= dim; ++d) {
1353     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr);
1354   }
1355   for (c = cMax; c < cEnd; ++c) {
1356     const PetscInt *cone;
1357     PetscInt       *closure = NULL;
1358     PetscInt        coneSize, closureSize, cl;
1359 
1360     if (label) {
1361       PetscInt val;
1362 
1363       ierr = DMLabelGetValue(label, c, &val);CHKERRQ(ierr);
1364       if (val != value) continue;
1365     }
1366     ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
1367     if (hasLagrange) {
1368       if (coneSize != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells");
1369     } else {
1370       if (coneSize != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells");
1371     }
1372     /* Negative face */
1373     ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr);
1374     ierr = DMPlexGetTransitiveClosure(dm, cone[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1375     for (cl = 0; cl < closureSize*2; cl += 2) {
1376       const PetscInt point = closure[cl];
1377 
1378       for (d = 0; d <= dim; ++d) {
1379         if ((point >= pStart[d]) && (point < pEnd[d])) {
1380           ierr = DMLabelSetValue(subpointMap, point, d);CHKERRQ(ierr);
1381           break;
1382         }
1383       }
1384     }
1385     ierr = DMPlexRestoreTransitiveClosure(dm, cone[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1386     /* Cells -- positive face is not included */
1387     for (cl = 0; cl < 1; ++cl) {
1388       const PetscInt *support;
1389       PetscInt        supportSize, s;
1390 
1391       ierr = DMPlexGetSupportSize(dm, cone[cl], &supportSize);CHKERRQ(ierr);
1392       if (supportSize != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive faces should separate two cells");
1393       ierr = DMPlexGetSupport(dm, cone[cl], &support);CHKERRQ(ierr);
1394       for (s = 0; s < supportSize; ++s) {
1395         ierr = DMLabelSetValue(subpointMap, support[s], dim);CHKERRQ(ierr);
1396       }
1397     }
1398   }
1399   ierr = PetscFree2(pStart, pEnd);CHKERRQ(ierr);
1400   PetscFunctionReturn(0);
1401 }
1402 
1403 #undef __FUNCT__
1404 #define __FUNCT__ "DMPlexGetFaceOrientation"
1405 PetscErrorCode DMPlexGetFaceOrientation(DM dm, PetscInt cell, PetscInt numCorners, PetscInt indices[], PetscInt oppositeVertex, PetscInt origVertices[], PetscInt faceVertices[], PetscBool *posOriented)
1406 {
1407   MPI_Comm       comm;
1408   PetscBool      posOrient = PETSC_FALSE;
1409   const PetscInt debug     = 0;
1410   PetscInt       cellDim, faceSize, f;
1411   PetscErrorCode ierr;
1412 
1413   PetscFunctionBegin;
1414   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1415   ierr = DMPlexGetDimension(dm, &cellDim);CHKERRQ(ierr);
1416   if (debug) {PetscPrintf(comm, "cellDim: %d numCorners: %d\n", cellDim, numCorners);CHKERRQ(ierr);}
1417 
1418   if (cellDim == numCorners-1) {
1419     /* Simplices */
1420     faceSize  = numCorners-1;
1421     posOrient = !(oppositeVertex%2) ? PETSC_TRUE : PETSC_FALSE;
1422   } else if (cellDim == 1 && numCorners == 3) {
1423     /* Quadratic line */
1424     faceSize  = 1;
1425     posOrient = PETSC_TRUE;
1426   } else if (cellDim == 2 && numCorners == 4) {
1427     /* Quads */
1428     faceSize = 2;
1429     if ((indices[1] > indices[0]) && (indices[1] - indices[0] == 1)) {
1430       posOrient = PETSC_TRUE;
1431     } else if ((indices[0] == 3) && (indices[1] == 0)) {
1432       posOrient = PETSC_TRUE;
1433     } else {
1434       if (((indices[0] > indices[1]) && (indices[0] - indices[1] == 1)) || ((indices[0] == 0) && (indices[1] == 3))) {
1435         posOrient = PETSC_FALSE;
1436       } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid quad crossedge");
1437     }
1438   } else if (cellDim == 2 && numCorners == 6) {
1439     /* Quadratic triangle (I hate this) */
1440     /* Edges are determined by the first 2 vertices (corners of edges) */
1441     const PetscInt faceSizeTri = 3;
1442     PetscInt       sortedIndices[3], i, iFace;
1443     PetscBool      found                    = PETSC_FALSE;
1444     PetscInt       faceVerticesTriSorted[9] = {
1445       0, 3,  4, /* bottom */
1446       1, 4,  5, /* right */
1447       2, 3,  5, /* left */
1448     };
1449     PetscInt       faceVerticesTri[9] = {
1450       0, 3,  4, /* bottom */
1451       1, 4,  5, /* right */
1452       2, 5,  3, /* left */
1453     };
1454 
1455     faceSize = faceSizeTri;
1456     for (i = 0; i < faceSizeTri; ++i) sortedIndices[i] = indices[i];
1457     ierr = PetscSortInt(faceSizeTri, sortedIndices);CHKERRQ(ierr);
1458     for (iFace = 0; iFace < 3; ++iFace) {
1459       const PetscInt ii = iFace*faceSizeTri;
1460       PetscInt       fVertex, cVertex;
1461 
1462       if ((sortedIndices[0] == faceVerticesTriSorted[ii+0]) &&
1463           (sortedIndices[1] == faceVerticesTriSorted[ii+1])) {
1464         for (fVertex = 0; fVertex < faceSizeTri; ++fVertex) {
1465           for (cVertex = 0; cVertex < faceSizeTri; ++cVertex) {
1466             if (indices[cVertex] == faceVerticesTri[ii+fVertex]) {
1467               faceVertices[fVertex] = origVertices[cVertex];
1468               break;
1469             }
1470           }
1471         }
1472         found = PETSC_TRUE;
1473         break;
1474       }
1475     }
1476     if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid tri crossface");
1477     if (posOriented) *posOriented = PETSC_TRUE;
1478     PetscFunctionReturn(0);
1479   } else if (cellDim == 2 && numCorners == 9) {
1480     /* Quadratic quad (I hate this) */
1481     /* Edges are determined by the first 2 vertices (corners of edges) */
1482     const PetscInt faceSizeQuad = 3;
1483     PetscInt       sortedIndices[3], i, iFace;
1484     PetscBool      found                      = PETSC_FALSE;
1485     PetscInt       faceVerticesQuadSorted[12] = {
1486       0, 1,  4, /* bottom */
1487       1, 2,  5, /* right */
1488       2, 3,  6, /* top */
1489       0, 3,  7, /* left */
1490     };
1491     PetscInt       faceVerticesQuad[12] = {
1492       0, 1,  4, /* bottom */
1493       1, 2,  5, /* right */
1494       2, 3,  6, /* top */
1495       3, 0,  7, /* left */
1496     };
1497 
1498     faceSize = faceSizeQuad;
1499     for (i = 0; i < faceSizeQuad; ++i) sortedIndices[i] = indices[i];
1500     ierr = PetscSortInt(faceSizeQuad, sortedIndices);CHKERRQ(ierr);
1501     for (iFace = 0; iFace < 4; ++iFace) {
1502       const PetscInt ii = iFace*faceSizeQuad;
1503       PetscInt       fVertex, cVertex;
1504 
1505       if ((sortedIndices[0] == faceVerticesQuadSorted[ii+0]) &&
1506           (sortedIndices[1] == faceVerticesQuadSorted[ii+1])) {
1507         for (fVertex = 0; fVertex < faceSizeQuad; ++fVertex) {
1508           for (cVertex = 0; cVertex < faceSizeQuad; ++cVertex) {
1509             if (indices[cVertex] == faceVerticesQuad[ii+fVertex]) {
1510               faceVertices[fVertex] = origVertices[cVertex];
1511               break;
1512             }
1513           }
1514         }
1515         found = PETSC_TRUE;
1516         break;
1517       }
1518     }
1519     if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid quad crossface");
1520     if (posOriented) *posOriented = PETSC_TRUE;
1521     PetscFunctionReturn(0);
1522   } else if (cellDim == 3 && numCorners == 8) {
1523     /* Hexes
1524        A hex is two oriented quads with the normal of the first
1525        pointing up at the second.
1526 
1527           7---6
1528          /|  /|
1529         4---5 |
1530         | 3-|-2
1531         |/  |/
1532         0---1
1533 
1534         Faces are determined by the first 4 vertices (corners of faces) */
1535     const PetscInt faceSizeHex = 4;
1536     PetscInt       sortedIndices[4], i, iFace;
1537     PetscBool      found                     = PETSC_FALSE;
1538     PetscInt       faceVerticesHexSorted[24] = {
1539       0, 1, 2, 3,  /* bottom */
1540       4, 5, 6, 7,  /* top */
1541       0, 1, 4, 5,  /* front */
1542       1, 2, 5, 6,  /* right */
1543       2, 3, 6, 7,  /* back */
1544       0, 3, 4, 7,  /* left */
1545     };
1546     PetscInt       faceVerticesHex[24] = {
1547       3, 2, 1, 0,  /* bottom */
1548       4, 5, 6, 7,  /* top */
1549       0, 1, 5, 4,  /* front */
1550       1, 2, 6, 5,  /* right */
1551       2, 3, 7, 6,  /* back */
1552       3, 0, 4, 7,  /* left */
1553     };
1554 
1555     faceSize = faceSizeHex;
1556     for (i = 0; i < faceSizeHex; ++i) sortedIndices[i] = indices[i];
1557     ierr = PetscSortInt(faceSizeHex, sortedIndices);CHKERRQ(ierr);
1558     for (iFace = 0; iFace < 6; ++iFace) {
1559       const PetscInt ii = iFace*faceSizeHex;
1560       PetscInt       fVertex, cVertex;
1561 
1562       if ((sortedIndices[0] == faceVerticesHexSorted[ii+0]) &&
1563           (sortedIndices[1] == faceVerticesHexSorted[ii+1]) &&
1564           (sortedIndices[2] == faceVerticesHexSorted[ii+2]) &&
1565           (sortedIndices[3] == faceVerticesHexSorted[ii+3])) {
1566         for (fVertex = 0; fVertex < faceSizeHex; ++fVertex) {
1567           for (cVertex = 0; cVertex < faceSizeHex; ++cVertex) {
1568             if (indices[cVertex] == faceVerticesHex[ii+fVertex]) {
1569               faceVertices[fVertex] = origVertices[cVertex];
1570               break;
1571             }
1572           }
1573         }
1574         found = PETSC_TRUE;
1575         break;
1576       }
1577     }
1578     if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid hex crossface");
1579     if (posOriented) *posOriented = PETSC_TRUE;
1580     PetscFunctionReturn(0);
1581   } else if (cellDim == 3 && numCorners == 10) {
1582     /* Quadratic tet */
1583     /* Faces are determined by the first 3 vertices (corners of faces) */
1584     const PetscInt faceSizeTet = 6;
1585     PetscInt       sortedIndices[6], i, iFace;
1586     PetscBool      found                     = PETSC_FALSE;
1587     PetscInt       faceVerticesTetSorted[24] = {
1588       0, 1, 2,  6, 7, 8, /* bottom */
1589       0, 3, 4,  6, 7, 9,  /* front */
1590       1, 4, 5,  7, 8, 9,  /* right */
1591       2, 3, 5,  6, 8, 9,  /* left */
1592     };
1593     PetscInt       faceVerticesTet[24] = {
1594       0, 1, 2,  6, 7, 8, /* bottom */
1595       0, 4, 3,  6, 7, 9,  /* front */
1596       1, 5, 4,  7, 8, 9,  /* right */
1597       2, 3, 5,  8, 6, 9,  /* left */
1598     };
1599 
1600     faceSize = faceSizeTet;
1601     for (i = 0; i < faceSizeTet; ++i) sortedIndices[i] = indices[i];
1602     ierr = PetscSortInt(faceSizeTet, sortedIndices);CHKERRQ(ierr);
1603     for (iFace=0; iFace < 4; ++iFace) {
1604       const PetscInt ii = iFace*faceSizeTet;
1605       PetscInt       fVertex, cVertex;
1606 
1607       if ((sortedIndices[0] == faceVerticesTetSorted[ii+0]) &&
1608           (sortedIndices[1] == faceVerticesTetSorted[ii+1]) &&
1609           (sortedIndices[2] == faceVerticesTetSorted[ii+2]) &&
1610           (sortedIndices[3] == faceVerticesTetSorted[ii+3])) {
1611         for (fVertex = 0; fVertex < faceSizeTet; ++fVertex) {
1612           for (cVertex = 0; cVertex < faceSizeTet; ++cVertex) {
1613             if (indices[cVertex] == faceVerticesTet[ii+fVertex]) {
1614               faceVertices[fVertex] = origVertices[cVertex];
1615               break;
1616             }
1617           }
1618         }
1619         found = PETSC_TRUE;
1620         break;
1621       }
1622     }
1623     if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid tet crossface");
1624     if (posOriented) *posOriented = PETSC_TRUE;
1625     PetscFunctionReturn(0);
1626   } else if (cellDim == 3 && numCorners == 27) {
1627     /* Quadratic hexes (I hate this)
1628        A hex is two oriented quads with the normal of the first
1629        pointing up at the second.
1630 
1631          7---6
1632         /|  /|
1633        4---5 |
1634        | 3-|-2
1635        |/  |/
1636        0---1
1637 
1638        Faces are determined by the first 4 vertices (corners of faces) */
1639     const PetscInt faceSizeQuadHex = 9;
1640     PetscInt       sortedIndices[9], i, iFace;
1641     PetscBool      found                         = PETSC_FALSE;
1642     PetscInt       faceVerticesQuadHexSorted[54] = {
1643       0, 1, 2, 3,  8, 9, 10, 11,  24, /* bottom */
1644       4, 5, 6, 7,  12, 13, 14, 15,  25, /* top */
1645       0, 1, 4, 5,  8, 12, 16, 17,  22, /* front */
1646       1, 2, 5, 6,  9, 13, 17, 18,  21, /* right */
1647       2, 3, 6, 7,  10, 14, 18, 19,  23, /* back */
1648       0, 3, 4, 7,  11, 15, 16, 19,  20, /* left */
1649     };
1650     PetscInt       faceVerticesQuadHex[54] = {
1651       3, 2, 1, 0,  10, 9, 8, 11,  24, /* bottom */
1652       4, 5, 6, 7,  12, 13, 14, 15,  25, /* top */
1653       0, 1, 5, 4,  8, 17, 12, 16,  22, /* front */
1654       1, 2, 6, 5,  9, 18, 13, 17,  21, /* right */
1655       2, 3, 7, 6,  10, 19, 14, 18,  23, /* back */
1656       3, 0, 4, 7,  11, 16, 15, 19,  20 /* left */
1657     };
1658 
1659     faceSize = faceSizeQuadHex;
1660     for (i = 0; i < faceSizeQuadHex; ++i) sortedIndices[i] = indices[i];
1661     ierr = PetscSortInt(faceSizeQuadHex, sortedIndices);CHKERRQ(ierr);
1662     for (iFace = 0; iFace < 6; ++iFace) {
1663       const PetscInt ii = iFace*faceSizeQuadHex;
1664       PetscInt       fVertex, cVertex;
1665 
1666       if ((sortedIndices[0] == faceVerticesQuadHexSorted[ii+0]) &&
1667           (sortedIndices[1] == faceVerticesQuadHexSorted[ii+1]) &&
1668           (sortedIndices[2] == faceVerticesQuadHexSorted[ii+2]) &&
1669           (sortedIndices[3] == faceVerticesQuadHexSorted[ii+3])) {
1670         for (fVertex = 0; fVertex < faceSizeQuadHex; ++fVertex) {
1671           for (cVertex = 0; cVertex < faceSizeQuadHex; ++cVertex) {
1672             if (indices[cVertex] == faceVerticesQuadHex[ii+fVertex]) {
1673               faceVertices[fVertex] = origVertices[cVertex];
1674               break;
1675             }
1676           }
1677         }
1678         found = PETSC_TRUE;
1679         break;
1680       }
1681     }
1682     if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid hex crossface");
1683     if (posOriented) *posOriented = PETSC_TRUE;
1684     PetscFunctionReturn(0);
1685   } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Unknown cell type for faceOrientation().");
1686   if (!posOrient) {
1687     if (debug) {ierr = PetscPrintf(comm, "  Reversing initial face orientation\n");CHKERRQ(ierr);}
1688     for (f = 0; f < faceSize; ++f) faceVertices[f] = origVertices[faceSize-1 - f];
1689   } else {
1690     if (debug) {ierr = PetscPrintf(comm, "  Keeping initial face orientation\n");CHKERRQ(ierr);}
1691     for (f = 0; f < faceSize; ++f) faceVertices[f] = origVertices[f];
1692   }
1693   if (posOriented) *posOriented = posOrient;
1694   PetscFunctionReturn(0);
1695 }
1696 
1697 #undef __FUNCT__
1698 #define __FUNCT__ "DMPlexGetOrientedFace"
1699 /*
1700     Given a cell and a face, as a set of vertices,
1701       return the oriented face, as a set of vertices, in faceVertices
1702     The orientation is such that the face normal points out of the cell
1703 */
1704 PetscErrorCode DMPlexGetOrientedFace(DM dm, PetscInt cell, PetscInt faceSize, const PetscInt face[], PetscInt numCorners, PetscInt indices[], PetscInt origVertices[], PetscInt faceVertices[], PetscBool *posOriented)
1705 {
1706   const PetscInt *cone = NULL;
1707   PetscInt        coneSize, v, f, v2;
1708   PetscInt        oppositeVertex = -1;
1709   PetscErrorCode  ierr;
1710 
1711   PetscFunctionBegin;
1712   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
1713   ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
1714   for (v = 0, v2 = 0; v < coneSize; ++v) {
1715     PetscBool found = PETSC_FALSE;
1716 
1717     for (f = 0; f < faceSize; ++f) {
1718       if (face[f] == cone[v]) {
1719         found = PETSC_TRUE; break;
1720       }
1721     }
1722     if (found) {
1723       indices[v2]      = v;
1724       origVertices[v2] = cone[v];
1725       ++v2;
1726     } else {
1727       oppositeVertex = v;
1728     }
1729   }
1730   ierr = DMPlexGetFaceOrientation(dm, cell, numCorners, indices, oppositeVertex, origVertices, faceVertices, posOriented);CHKERRQ(ierr);
1731   PetscFunctionReturn(0);
1732 }
1733 
1734 #undef __FUNCT__
1735 #define __FUNCT__ "DMPlexInsertFace_Internal"
1736 /*
1737   DMPlexInsertFace_Internal - Puts a face into the mesh
1738 
1739   Not collective
1740 
1741   Input Parameters:
1742   + dm              - The DMPlex
1743   . numFaceVertex   - The number of vertices in the face
1744   . faceVertices    - The vertices in the face for dm
1745   . subfaceVertices - The vertices in the face for subdm
1746   . numCorners      - The number of vertices in the cell
1747   . cell            - A cell in dm containing the face
1748   . subcell         - A cell in subdm containing the face
1749   . firstFace       - First face in the mesh
1750   - newFacePoint    - Next face in the mesh
1751 
1752   Output Parameters:
1753   . newFacePoint - Contains next face point number on input, updated on output
1754 
1755   Level: developer
1756 */
1757 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)
1758 {
1759   MPI_Comm        comm;
1760   DM_Plex        *submesh = (DM_Plex*) subdm->data;
1761   const PetscInt *faces;
1762   PetscInt        numFaces, coneSize;
1763   PetscErrorCode  ierr;
1764 
1765   PetscFunctionBegin;
1766   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1767   ierr = DMPlexGetConeSize(subdm, subcell, &coneSize);CHKERRQ(ierr);
1768   if (coneSize != 1) SETERRQ2(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cone size of cell %d is %d != 1", cell, coneSize);
1769 #if 0
1770   /* Cannot use this because support() has not been constructed yet */
1771   ierr = DMPlexGetJoin(subdm, numFaceVertices, subfaceVertices, &numFaces, &faces);CHKERRQ(ierr);
1772 #else
1773   {
1774     PetscInt f;
1775 
1776     numFaces = 0;
1777     ierr     = DMGetWorkArray(subdm, 1, PETSC_INT, (void **) &faces);CHKERRQ(ierr);
1778     for (f = firstFace; f < *newFacePoint; ++f) {
1779       PetscInt dof, off, d;
1780 
1781       ierr = PetscSectionGetDof(submesh->coneSection, f, &dof);CHKERRQ(ierr);
1782       ierr = PetscSectionGetOffset(submesh->coneSection, f, &off);CHKERRQ(ierr);
1783       /* Yes, I know this is quadratic, but I expect the sizes to be <5 */
1784       for (d = 0; d < dof; ++d) {
1785         const PetscInt p = submesh->cones[off+d];
1786         PetscInt       v;
1787 
1788         for (v = 0; v < numFaceVertices; ++v) {
1789           if (subfaceVertices[v] == p) break;
1790         }
1791         if (v == numFaceVertices) break;
1792       }
1793       if (d == dof) {
1794         numFaces               = 1;
1795         ((PetscInt*) faces)[0] = f;
1796       }
1797     }
1798   }
1799 #endif
1800   if (numFaces > 1) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Vertex set had %d faces, not one", numFaces);
1801   else if (numFaces == 1) {
1802     /* Add the other cell neighbor for this face */
1803     ierr = DMPlexSetCone(subdm, subcell, faces);CHKERRQ(ierr);
1804   } else {
1805     PetscInt *indices, *origVertices, *orientedVertices, *orientedSubVertices, v, ov;
1806     PetscBool posOriented;
1807 
1808     ierr                = DMGetWorkArray(subdm, 4*numFaceVertices * sizeof(PetscInt), PETSC_INT, &orientedVertices);CHKERRQ(ierr);
1809     origVertices        = &orientedVertices[numFaceVertices];
1810     indices             = &orientedVertices[numFaceVertices*2];
1811     orientedSubVertices = &orientedVertices[numFaceVertices*3];
1812     ierr                = DMPlexGetOrientedFace(dm, cell, numFaceVertices, faceVertices, numCorners, indices, origVertices, orientedVertices, &posOriented);CHKERRQ(ierr);
1813     /* TODO: I know that routine should return a permutation, not the indices */
1814     for (v = 0; v < numFaceVertices; ++v) {
1815       const PetscInt vertex = faceVertices[v], subvertex = subfaceVertices[v];
1816       for (ov = 0; ov < numFaceVertices; ++ov) {
1817         if (orientedVertices[ov] == vertex) {
1818           orientedSubVertices[ov] = subvertex;
1819           break;
1820         }
1821       }
1822       if (ov == numFaceVertices) SETERRQ1(comm, PETSC_ERR_PLIB, "Could not find face vertex %d in orientated set", vertex);
1823     }
1824     ierr = DMPlexSetCone(subdm, *newFacePoint, orientedSubVertices);CHKERRQ(ierr);
1825     ierr = DMPlexSetCone(subdm, subcell, newFacePoint);CHKERRQ(ierr);
1826     ierr = DMRestoreWorkArray(subdm, 4*numFaceVertices * sizeof(PetscInt), PETSC_INT, &orientedVertices);CHKERRQ(ierr);
1827     ++(*newFacePoint);
1828   }
1829 #if 0
1830   ierr = DMPlexRestoreJoin(subdm, numFaceVertices, subfaceVertices, &numFaces, &faces);CHKERRQ(ierr);
1831 #else
1832   ierr = DMRestoreWorkArray(subdm, 1, PETSC_INT, (void **) &faces);CHKERRQ(ierr);
1833 #endif
1834   PetscFunctionReturn(0);
1835 }
1836 
1837 #undef __FUNCT__
1838 #define __FUNCT__ "DMPlexCreateSubmesh_Uninterpolated"
1839 static PetscErrorCode DMPlexCreateSubmesh_Uninterpolated(DM dm, const char vertexLabelName[], PetscInt value, DM subdm)
1840 {
1841   MPI_Comm        comm;
1842   DMLabel         vertexLabel, subpointMap;
1843   IS              subvertexIS,  subcellIS;
1844   const PetscInt *subVertices, *subCells;
1845   PetscInt        numSubVertices, firstSubVertex, numSubCells;
1846   PetscInt       *subface, maxConeSize, numSubFaces, firstSubFace, newFacePoint, nFV;
1847   PetscInt        vStart, vEnd, c, f;
1848   PetscErrorCode  ierr;
1849 
1850   PetscFunctionBegin;
1851   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1852   /* Create subpointMap which marks the submesh */
1853   ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr);
1854   ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr);
1855   ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr);
1856   ierr = DMPlexGetLabel(dm, vertexLabelName, &vertexLabel);CHKERRQ(ierr);
1857   ierr = DMPlexMarkSubmesh_Uninterpolated(dm, vertexLabel, value, subpointMap, &numSubFaces, &nFV, subdm);CHKERRQ(ierr);
1858   /* Setup chart */
1859   ierr = DMLabelGetStratumSize(subpointMap, 0, &numSubVertices);CHKERRQ(ierr);
1860   ierr = DMLabelGetStratumSize(subpointMap, 2, &numSubCells);CHKERRQ(ierr);
1861   ierr = DMPlexSetChart(subdm, 0, numSubCells+numSubFaces+numSubVertices);CHKERRQ(ierr);
1862   ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr);
1863   /* Set cone sizes */
1864   firstSubVertex = numSubCells;
1865   firstSubFace   = numSubCells+numSubVertices;
1866   newFacePoint   = firstSubFace;
1867   ierr = DMLabelGetStratumIS(subpointMap, 0, &subvertexIS);CHKERRQ(ierr);
1868   if (subvertexIS) {ierr = ISGetIndices(subvertexIS, &subVertices);CHKERRQ(ierr);}
1869   ierr = DMLabelGetStratumIS(subpointMap, 2, &subcellIS);CHKERRQ(ierr);
1870   if (subcellIS) {ierr = ISGetIndices(subcellIS, &subCells);CHKERRQ(ierr);}
1871   for (c = 0; c < numSubCells; ++c) {
1872     ierr = DMPlexSetConeSize(subdm, c, 1);CHKERRQ(ierr);
1873   }
1874   for (f = firstSubFace; f < firstSubFace+numSubFaces; ++f) {
1875     ierr = DMPlexSetConeSize(subdm, f, nFV);CHKERRQ(ierr);
1876   }
1877   ierr = DMSetUp(subdm);CHKERRQ(ierr);
1878   /* Create face cones */
1879   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
1880   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
1881   ierr = DMGetWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr);
1882   for (c = 0; c < numSubCells; ++c) {
1883     const PetscInt cell    = subCells[c];
1884     const PetscInt subcell = c;
1885     PetscInt      *closure = NULL;
1886     PetscInt       closureSize, cl, numCorners = 0, faceSize = 0;
1887 
1888     ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1889     for (cl = 0; cl < closureSize*2; cl += 2) {
1890       const PetscInt point = closure[cl];
1891       PetscInt       subVertex;
1892 
1893       if ((point >= vStart) && (point < vEnd)) {
1894         ++numCorners;
1895         ierr = PetscFindInt(point, numSubVertices, subVertices, &subVertex);CHKERRQ(ierr);
1896         if (subVertex >= 0) {
1897           closure[faceSize] = point;
1898           subface[faceSize] = firstSubVertex+subVertex;
1899           ++faceSize;
1900         }
1901       }
1902     }
1903     if (faceSize > nFV) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Invalid submesh: Too many vertices %d of an element on the surface", faceSize);
1904     if (faceSize == nFV) {
1905       ierr = DMPlexInsertFace_Internal(dm, subdm, faceSize, closure, subface, numCorners, cell, subcell, firstSubFace, &newFacePoint);CHKERRQ(ierr);
1906     }
1907     ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1908   }
1909   ierr = DMRestoreWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr);
1910   ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr);
1911   ierr = DMPlexStratify(subdm);CHKERRQ(ierr);
1912   /* Build coordinates */
1913   {
1914     PetscSection coordSection, subCoordSection;
1915     Vec          coordinates, subCoordinates;
1916     PetscScalar *coords, *subCoords;
1917     PetscInt     numComp, coordSize, v;
1918 
1919     ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1920     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1921     ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr);
1922     ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr);
1923     ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr);
1924     ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr);
1925     ierr = PetscSectionSetChart(subCoordSection, firstSubVertex, firstSubVertex+numSubVertices);CHKERRQ(ierr);
1926     for (v = 0; v < numSubVertices; ++v) {
1927       const PetscInt vertex    = subVertices[v];
1928       const PetscInt subvertex = firstSubVertex+v;
1929       PetscInt       dof;
1930 
1931       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
1932       ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr);
1933       ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr);
1934     }
1935     ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr);
1936     ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr);
1937     ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr);
1938     ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
1939     ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr);
1940     if (coordSize) {
1941       ierr = VecGetArray(coordinates,    &coords);CHKERRQ(ierr);
1942       ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr);
1943       for (v = 0; v < numSubVertices; ++v) {
1944         const PetscInt vertex    = subVertices[v];
1945         const PetscInt subvertex = firstSubVertex+v;
1946         PetscInt       dof, off, sdof, soff, d;
1947 
1948         ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
1949         ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr);
1950         ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr);
1951         ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr);
1952         if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof);
1953         for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d];
1954       }
1955       ierr = VecRestoreArray(coordinates,    &coords);CHKERRQ(ierr);
1956       ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr);
1957     }
1958     ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr);
1959     ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr);
1960   }
1961   /* Cleanup */
1962   if (subvertexIS) {ierr = ISRestoreIndices(subvertexIS, &subVertices);CHKERRQ(ierr);}
1963   ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr);
1964   if (subcellIS) {ierr = ISRestoreIndices(subcellIS, &subCells);CHKERRQ(ierr);}
1965   ierr = ISDestroy(&subcellIS);CHKERRQ(ierr);
1966   PetscFunctionReturn(0);
1967 }
1968 
1969 #undef __FUNCT__
1970 #define __FUNCT__ "DMPlexCreateSubmesh_Interpolated"
1971 static PetscErrorCode DMPlexCreateSubmesh_Interpolated(DM dm, const char vertexLabelName[], PetscInt value, DM subdm)
1972 {
1973   MPI_Comm         comm;
1974   DMLabel          subpointMap, vertexLabel;
1975   IS              *subpointIS;
1976   const PetscInt **subpoints;
1977   PetscInt        *numSubPoints, *firstSubPoint, *coneNew;
1978   PetscInt         totSubPoints = 0, maxConeSize, dim, p, d, v;
1979   PetscErrorCode   ierr;
1980 
1981   PetscFunctionBegin;
1982   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1983   /* Create subpointMap which marks the submesh */
1984   ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr);
1985   ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr);
1986   ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr);
1987   ierr = DMPlexGetLabel(dm, vertexLabelName, &vertexLabel);CHKERRQ(ierr);
1988   ierr = DMPlexMarkSubmesh_Interpolated(dm, vertexLabel, value, subpointMap, subdm);CHKERRQ(ierr);
1989   /* Setup chart */
1990   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1991   ierr = PetscMalloc4(dim+1,PetscInt,&numSubPoints,dim+1,PetscInt,&firstSubPoint,dim+1,IS,&subpointIS,dim+1,const PetscInt *,&subpoints);CHKERRQ(ierr);
1992   for (d = 0; d <= dim; ++d) {
1993     ierr = DMLabelGetStratumSize(subpointMap, d, &numSubPoints[d]);CHKERRQ(ierr);
1994     totSubPoints += numSubPoints[d];
1995   }
1996   ierr = DMPlexSetChart(subdm, 0, totSubPoints);CHKERRQ(ierr);
1997   ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr);
1998   /* Set cone sizes */
1999   firstSubPoint[dim] = 0;
2000   firstSubPoint[0]   = firstSubPoint[dim] + numSubPoints[dim];
2001   if (dim > 1) {firstSubPoint[dim-1] = firstSubPoint[0]     + numSubPoints[0];}
2002   if (dim > 2) {firstSubPoint[dim-2] = firstSubPoint[dim-1] + numSubPoints[dim-1];}
2003   for (d = 0; d <= dim; ++d) {
2004     ierr = DMLabelGetStratumIS(subpointMap, d, &subpointIS[d]);CHKERRQ(ierr);
2005     ierr = ISGetIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr);
2006   }
2007   for (d = 0; d <= dim; ++d) {
2008     for (p = 0; p < numSubPoints[d]; ++p) {
2009       const PetscInt  point    = subpoints[d][p];
2010       const PetscInt  subpoint = firstSubPoint[d] + p;
2011       const PetscInt *cone;
2012       PetscInt        coneSize, coneSizeNew, c, val;
2013 
2014       ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
2015       ierr = DMPlexSetConeSize(subdm, subpoint, coneSize);CHKERRQ(ierr);
2016       if (d == dim) {
2017         ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
2018         for (c = 0, coneSizeNew = 0; c < coneSize; ++c) {
2019           ierr = DMLabelGetValue(subpointMap, cone[c], &val);CHKERRQ(ierr);
2020           if (val >= 0) coneSizeNew++;
2021         }
2022         ierr = DMPlexSetConeSize(subdm, subpoint, coneSizeNew);CHKERRQ(ierr);
2023       }
2024     }
2025   }
2026   ierr = DMSetUp(subdm);CHKERRQ(ierr);
2027   /* Set cones */
2028   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
2029   ierr = PetscMalloc(maxConeSize * sizeof(PetscInt), &coneNew);CHKERRQ(ierr);
2030   for (d = 0; d <= dim; ++d) {
2031     for (p = 0; p < numSubPoints[d]; ++p) {
2032       const PetscInt  point    = subpoints[d][p];
2033       const PetscInt  subpoint = firstSubPoint[d] + p;
2034       const PetscInt *cone;
2035       PetscInt        coneSize, subconeSize, coneSizeNew, c, subc;
2036 
2037       ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
2038       ierr = DMPlexGetConeSize(subdm, subpoint, &subconeSize);CHKERRQ(ierr);
2039       ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
2040       for (c = 0, coneSizeNew = 0; c < coneSize; ++c) {
2041         ierr = PetscFindInt(cone[c], numSubPoints[d-1], subpoints[d-1], &subc);CHKERRQ(ierr);
2042         if (subc >= 0) coneNew[coneSizeNew++] = firstSubPoint[d-1] + subc;
2043       }
2044       if (coneSizeNew != subconeSize) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of cone points located %d does not match subcone size %d", coneSizeNew, subconeSize);
2045       ierr = DMPlexSetCone(subdm, subpoint, coneNew);CHKERRQ(ierr);
2046     }
2047   }
2048   ierr = PetscFree(coneNew);CHKERRQ(ierr);
2049   ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr);
2050   ierr = DMPlexStratify(subdm);CHKERRQ(ierr);
2051   /* Build coordinates */
2052   {
2053     PetscSection coordSection, subCoordSection;
2054     Vec          coordinates, subCoordinates;
2055     PetscScalar *coords, *subCoords;
2056     PetscInt     numComp, coordSize;
2057 
2058     ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
2059     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
2060     ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr);
2061     ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr);
2062     ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr);
2063     ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr);
2064     ierr = PetscSectionSetChart(subCoordSection, firstSubPoint[0], firstSubPoint[0]+numSubPoints[0]);CHKERRQ(ierr);
2065     for (v = 0; v < numSubPoints[0]; ++v) {
2066       const PetscInt vertex    = subpoints[0][v];
2067       const PetscInt subvertex = firstSubPoint[0]+v;
2068       PetscInt       dof;
2069 
2070       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
2071       ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr);
2072       ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr);
2073     }
2074     ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr);
2075     ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr);
2076     ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr);
2077     ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
2078     ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr);
2079     ierr = VecGetArray(coordinates,    &coords);CHKERRQ(ierr);
2080     ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr);
2081     for (v = 0; v < numSubPoints[0]; ++v) {
2082       const PetscInt vertex    = subpoints[0][v];
2083       const PetscInt subvertex = firstSubPoint[0]+v;
2084       PetscInt dof, off, sdof, soff, d;
2085 
2086       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
2087       ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr);
2088       ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr);
2089       ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr);
2090       if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof);
2091       for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d];
2092     }
2093     ierr = VecRestoreArray(coordinates,    &coords);CHKERRQ(ierr);
2094     ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr);
2095     ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr);
2096     ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr);
2097   }
2098   /* Cleanup */
2099   for (d = 0; d <= dim; ++d) {
2100     ierr = ISRestoreIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr);
2101     ierr = ISDestroy(&subpointIS[d]);CHKERRQ(ierr);
2102   }
2103   ierr = PetscFree4(numSubPoints,firstSubPoint,subpointIS,subpoints);CHKERRQ(ierr);
2104   PetscFunctionReturn(0);
2105 }
2106 
2107 #undef __FUNCT__
2108 #define __FUNCT__ "DMPlexCreateSubmesh"
2109 /*@C
2110   DMPlexCreateSubmesh - Extract a hypersurface from the mesh using vertices defined by a label
2111 
2112   Input Parameters:
2113 + dm           - The original mesh
2114 . vertexLabel  - The DMLabel marking vertices contained in the surface
2115 - value        - The label value to use
2116 
2117   Output Parameter:
2118 . subdm - The surface mesh
2119 
2120   Note: This function produces a DMLabel mapping original points in the submesh to their depth. This can be obtained using DMPlexGetSubpointMap().
2121 
2122   Level: developer
2123 
2124 .seealso: DMPlexGetSubpointMap(), DMPlexGetLabel(), DMLabelSetValue()
2125 @*/
2126 PetscErrorCode DMPlexCreateSubmesh(DM dm, const char vertexLabel[], PetscInt value, DM *subdm)
2127 {
2128   PetscInt       dim, depth;
2129   PetscErrorCode ierr;
2130 
2131   PetscFunctionBegin;
2132   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2133   PetscValidCharPointer(vertexLabel, 2);
2134   PetscValidPointer(subdm, 3);
2135   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
2136   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
2137   ierr = DMCreate(PetscObjectComm((PetscObject)dm), subdm);CHKERRQ(ierr);
2138   ierr = DMSetType(*subdm, DMPLEX);CHKERRQ(ierr);
2139   ierr = DMPlexSetDimension(*subdm, dim-1);CHKERRQ(ierr);
2140   if (depth == dim) {
2141     ierr = DMPlexCreateSubmesh_Interpolated(dm, vertexLabel, value, *subdm);CHKERRQ(ierr);
2142   } else {
2143     ierr = DMPlexCreateSubmesh_Uninterpolated(dm, vertexLabel, value, *subdm);CHKERRQ(ierr);
2144   }
2145   PetscFunctionReturn(0);
2146 }
2147 
2148 #undef __FUNCT__
2149 #define __FUNCT__ "DMPlexCreateCohesiveSubmesh_Uninterpolated"
2150 static PetscErrorCode DMPlexCreateCohesiveSubmesh_Uninterpolated(DM dm, PetscBool hasLagrange, const char label[], PetscInt value, DM subdm)
2151 {
2152   MPI_Comm        comm;
2153   DMLabel         subpointMap;
2154   IS              subvertexIS;
2155   const PetscInt *subVertices;
2156   PetscInt        numSubVertices, firstSubVertex, numSubCells, *subCells;
2157   PetscInt       *subface, maxConeSize, numSubFaces, firstSubFace, newFacePoint, nFV;
2158   PetscInt        cMax, c, f;
2159   PetscErrorCode  ierr;
2160 
2161   PetscFunctionBegin;
2162   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
2163   /* Create subpointMap which marks the submesh */
2164   ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr);
2165   ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr);
2166   ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr);
2167   ierr = DMPlexMarkCohesiveSubmesh_Uninterpolated(dm, hasLagrange, label, value, subpointMap, &numSubFaces, &nFV, &subCells, subdm);CHKERRQ(ierr);
2168   /* Setup chart */
2169   ierr = DMLabelGetStratumSize(subpointMap, 0, &numSubVertices);CHKERRQ(ierr);
2170   ierr = DMLabelGetStratumSize(subpointMap, 2, &numSubCells);CHKERRQ(ierr);
2171   ierr = DMPlexSetChart(subdm, 0, numSubCells+numSubFaces+numSubVertices);CHKERRQ(ierr);
2172   ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr);
2173   /* Set cone sizes */
2174   firstSubVertex = numSubCells;
2175   firstSubFace   = numSubCells+numSubVertices;
2176   newFacePoint   = firstSubFace;
2177   ierr = DMLabelGetStratumIS(subpointMap, 0, &subvertexIS);CHKERRQ(ierr);
2178   if (subvertexIS) {ierr = ISGetIndices(subvertexIS, &subVertices);CHKERRQ(ierr);}
2179   for (c = 0; c < numSubCells; ++c) {
2180     ierr = DMPlexSetConeSize(subdm, c, 1);CHKERRQ(ierr);
2181   }
2182   for (f = firstSubFace; f < firstSubFace+numSubFaces; ++f) {
2183     ierr = DMPlexSetConeSize(subdm, f, nFV);CHKERRQ(ierr);
2184   }
2185   ierr = DMSetUp(subdm);CHKERRQ(ierr);
2186   /* Create face cones */
2187   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
2188   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
2189   ierr = DMGetWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr);
2190   for (c = 0; c < numSubCells; ++c) {
2191     const PetscInt  cell    = subCells[c];
2192     const PetscInt  subcell = c;
2193     const PetscInt *cone, *cells;
2194     PetscInt        numCells, subVertex, p, v;
2195 
2196     if (cell < cMax) continue;
2197     ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
2198     for (v = 0; v < nFV; ++v) {
2199       ierr = PetscFindInt(cone[v], numSubVertices, subVertices, &subVertex);CHKERRQ(ierr);
2200       subface[v] = firstSubVertex+subVertex;
2201     }
2202     ierr = DMPlexSetCone(subdm, newFacePoint, subface);CHKERRQ(ierr);
2203     ierr = DMPlexSetCone(subdm, subcell, &newFacePoint);CHKERRQ(ierr);
2204     ierr = DMPlexGetJoin(dm, nFV, cone, &numCells, &cells);CHKERRQ(ierr);
2205     if (numCells != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells");
2206     for (p = 0; p < numCells; ++p) {
2207       PetscInt negsubcell;
2208 
2209       if (cells[p] >= cMax) continue;
2210       /* I know this is a crap search */
2211       for (negsubcell = 0; negsubcell < numSubCells; ++negsubcell) {
2212         if (subCells[negsubcell] == cells[p]) break;
2213       }
2214       if (negsubcell == numSubCells) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find negative face neighbor for cohesive cell %d", cell);
2215       ierr = DMPlexSetCone(subdm, negsubcell, &newFacePoint);CHKERRQ(ierr);
2216     }
2217     ierr = DMPlexRestoreJoin(dm, nFV, cone, &numCells, &cells);CHKERRQ(ierr);
2218     ++newFacePoint;
2219   }
2220   ierr = DMRestoreWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr);
2221   ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr);
2222   ierr = DMPlexStratify(subdm);CHKERRQ(ierr);
2223   /* Build coordinates */
2224   {
2225     PetscSection coordSection, subCoordSection;
2226     Vec          coordinates, subCoordinates;
2227     PetscScalar *coords, *subCoords;
2228     PetscInt     numComp, coordSize, v;
2229 
2230     ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
2231     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
2232     ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr);
2233     ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr);
2234     ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr);
2235     ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr);
2236     ierr = PetscSectionSetChart(subCoordSection, firstSubVertex, firstSubVertex+numSubVertices);CHKERRQ(ierr);
2237     for (v = 0; v < numSubVertices; ++v) {
2238       const PetscInt vertex    = subVertices[v];
2239       const PetscInt subvertex = firstSubVertex+v;
2240       PetscInt       dof;
2241 
2242       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
2243       ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr);
2244       ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr);
2245     }
2246     ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr);
2247     ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr);
2248     ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr);
2249     if (coordSize) {
2250       ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
2251       ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr);
2252       ierr = VecGetArray(coordinates,    &coords);CHKERRQ(ierr);
2253       ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr);
2254       for (v = 0; v < numSubVertices; ++v) {
2255         const PetscInt vertex    = subVertices[v];
2256         const PetscInt subvertex = firstSubVertex+v;
2257         PetscInt       dof, off, sdof, soff, d;
2258 
2259         ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
2260         ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr);
2261         ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr);
2262         ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr);
2263         if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof);
2264         for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d];
2265       }
2266       ierr = VecRestoreArray(coordinates,    &coords);CHKERRQ(ierr);
2267       ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr);
2268     }
2269     ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr);
2270     ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr);
2271   }
2272   /* Cleanup */
2273   if (subvertexIS) {ierr = ISRestoreIndices(subvertexIS, &subVertices);CHKERRQ(ierr);}
2274   ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr);
2275   ierr = PetscFree(subCells);CHKERRQ(ierr);
2276   PetscFunctionReturn(0);
2277 }
2278 
2279 #undef __FUNCT__
2280 #define __FUNCT__ "DMPlexCreateCohesiveSubmesh_Interpolated"
2281 static PetscErrorCode DMPlexCreateCohesiveSubmesh_Interpolated(DM dm, PetscBool hasLagrange, const char label[], PetscInt value, DM subdm)
2282 {
2283   MPI_Comm         comm;
2284   DMLabel          subpointMap;
2285   IS              *subpointIS;
2286   const PetscInt **subpoints;
2287   PetscInt        *numSubPoints, *firstSubPoint, *coneNew;
2288   PetscInt         totSubPoints = 0, maxConeSize, dim, p, d, v;
2289   PetscErrorCode   ierr;
2290 
2291   PetscFunctionBegin;
2292   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
2293   /* Create subpointMap which marks the submesh */
2294   ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr);
2295   ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr);
2296   ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr);
2297   ierr = DMPlexMarkCohesiveSubmesh_Interpolated(dm, hasLagrange, label, value, subpointMap, subdm);CHKERRQ(ierr);
2298   /* Setup chart */
2299   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
2300   ierr = PetscMalloc4(dim+1,PetscInt,&numSubPoints,dim+1,PetscInt,&firstSubPoint,dim+1,IS,&subpointIS,dim+1,const PetscInt *,&subpoints);CHKERRQ(ierr);
2301   for (d = 0; d <= dim; ++d) {
2302     ierr = DMLabelGetStratumSize(subpointMap, d, &numSubPoints[d]);CHKERRQ(ierr);
2303     totSubPoints += numSubPoints[d];
2304   }
2305   ierr = DMPlexSetChart(subdm, 0, totSubPoints);CHKERRQ(ierr);
2306   ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr);
2307   /* Set cone sizes */
2308   firstSubPoint[dim] = 0;
2309   firstSubPoint[0]   = firstSubPoint[dim] + numSubPoints[dim];
2310   if (dim > 1) {firstSubPoint[dim-1] = firstSubPoint[0]     + numSubPoints[0];}
2311   if (dim > 2) {firstSubPoint[dim-2] = firstSubPoint[dim-1] + numSubPoints[dim-1];}
2312   for (d = 0; d <= dim; ++d) {
2313     ierr = DMLabelGetStratumIS(subpointMap, d, &subpointIS[d]);CHKERRQ(ierr);
2314     ierr = ISGetIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr);
2315   }
2316   for (d = 0; d <= dim; ++d) {
2317     for (p = 0; p < numSubPoints[d]; ++p) {
2318       const PetscInt  point    = subpoints[d][p];
2319       const PetscInt  subpoint = firstSubPoint[d] + p;
2320       const PetscInt *cone;
2321       PetscInt        coneSize, coneSizeNew, c, val;
2322 
2323       ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
2324       ierr = DMPlexSetConeSize(subdm, subpoint, coneSize);CHKERRQ(ierr);
2325       if (d == dim) {
2326         ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
2327         for (c = 0, coneSizeNew = 0; c < coneSize; ++c) {
2328           ierr = DMLabelGetValue(subpointMap, cone[c], &val);CHKERRQ(ierr);
2329           if (val >= 0) coneSizeNew++;
2330         }
2331         ierr = DMPlexSetConeSize(subdm, subpoint, coneSizeNew);CHKERRQ(ierr);
2332       }
2333     }
2334   }
2335   ierr = DMSetUp(subdm);CHKERRQ(ierr);
2336   /* Set cones */
2337   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
2338   ierr = PetscMalloc(maxConeSize * sizeof(PetscInt), &coneNew);CHKERRQ(ierr);
2339   for (d = 0; d <= dim; ++d) {
2340     for (p = 0; p < numSubPoints[d]; ++p) {
2341       const PetscInt  point    = subpoints[d][p];
2342       const PetscInt  subpoint = firstSubPoint[d] + p;
2343       const PetscInt *cone;
2344       PetscInt        coneSize, subconeSize, coneSizeNew, c, subc;
2345 
2346       ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
2347       ierr = DMPlexGetConeSize(subdm, subpoint, &subconeSize);CHKERRQ(ierr);
2348       ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
2349       for (c = 0, coneSizeNew = 0; c < coneSize; ++c) {
2350         ierr = PetscFindInt(cone[c], numSubPoints[d-1], subpoints[d-1], &subc);CHKERRQ(ierr);
2351         if (subc >= 0) coneNew[coneSizeNew++] = firstSubPoint[d-1] + subc;
2352       }
2353       if (coneSizeNew != subconeSize) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of cone points located %d does not match subcone size %d", coneSizeNew, subconeSize);
2354       ierr = DMPlexSetCone(subdm, subpoint, coneNew);CHKERRQ(ierr);
2355     }
2356   }
2357   ierr = PetscFree(coneNew);CHKERRQ(ierr);
2358   ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr);
2359   ierr = DMPlexStratify(subdm);CHKERRQ(ierr);
2360   /* Build coordinates */
2361   {
2362     PetscSection coordSection, subCoordSection;
2363     Vec          coordinates, subCoordinates;
2364     PetscScalar *coords, *subCoords;
2365     PetscInt     numComp, coordSize;
2366 
2367     ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
2368     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
2369     ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr);
2370     ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr);
2371     ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr);
2372     ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr);
2373     ierr = PetscSectionSetChart(subCoordSection, firstSubPoint[0], firstSubPoint[0]+numSubPoints[0]);CHKERRQ(ierr);
2374     for (v = 0; v < numSubPoints[0]; ++v) {
2375       const PetscInt vertex    = subpoints[0][v];
2376       const PetscInt subvertex = firstSubPoint[0]+v;
2377       PetscInt       dof;
2378 
2379       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
2380       ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr);
2381       ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr);
2382     }
2383     ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr);
2384     ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr);
2385     ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr);
2386     ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
2387     ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr);
2388     ierr = VecGetArray(coordinates,    &coords);CHKERRQ(ierr);
2389     ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr);
2390     for (v = 0; v < numSubPoints[0]; ++v) {
2391       const PetscInt vertex    = subpoints[0][v];
2392       const PetscInt subvertex = firstSubPoint[0]+v;
2393       PetscInt dof, off, sdof, soff, d;
2394 
2395       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
2396       ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr);
2397       ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr);
2398       ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr);
2399       if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof);
2400       for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d];
2401     }
2402     ierr = VecRestoreArray(coordinates,    &coords);CHKERRQ(ierr);
2403     ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr);
2404     ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr);
2405     ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr);
2406   }
2407   /* Cleanup */
2408   for (d = 0; d <= dim; ++d) {
2409     ierr = ISRestoreIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr);
2410     ierr = ISDestroy(&subpointIS[d]);CHKERRQ(ierr);
2411   }
2412   ierr = PetscFree4(numSubPoints,firstSubPoint,subpointIS,subpoints);CHKERRQ(ierr);
2413   PetscFunctionReturn(0);
2414 }
2415 
2416 #undef __FUNCT__
2417 #define __FUNCT__ "DMPlexCreateCohesiveSubmesh"
2418 /*
2419   DMPlexCreateCohesiveSubmesh - Extract from a mesh with cohesive cells the hypersurface defined by one face of the cells. Optionally, a Label an be given to restrict the cells.
2420 
2421   Input Parameters:
2422 + dm          - The original mesh
2423 . hasLagrange - The mesh has Lagrange unknowns in the cohesive cells
2424 . label       - A label name, or PETSC_NULL
2425 - value  - A label value
2426 
2427   Output Parameter:
2428 . subdm - The surface mesh
2429 
2430   Note: This function produces a DMLabel mapping original points in the submesh to their depth. This can be obtained using DMPlexGetSubpointMap().
2431 
2432   Level: developer
2433 
2434 .seealso: DMPlexGetSubpointMap(), DMPlexCreateSubmesh()
2435 */
2436 PetscErrorCode DMPlexCreateCohesiveSubmesh(DM dm, PetscBool hasLagrange, const char label[], PetscInt value, DM *subdm)
2437 {
2438   PetscInt       dim, depth;
2439   PetscErrorCode ierr;
2440 
2441   PetscFunctionBegin;
2442   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2443   PetscValidPointer(subdm, 5);
2444   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
2445   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
2446   ierr = DMCreate(PetscObjectComm((PetscObject)dm), subdm);CHKERRQ(ierr);
2447   ierr = DMSetType(*subdm, DMPLEX);CHKERRQ(ierr);
2448   ierr = DMPlexSetDimension(*subdm, dim-1);CHKERRQ(ierr);
2449   if (depth == dim) {
2450     ierr = DMPlexCreateCohesiveSubmesh_Interpolated(dm, hasLagrange, label, value, *subdm);CHKERRQ(ierr);
2451   } else {
2452     ierr = DMPlexCreateCohesiveSubmesh_Uninterpolated(dm, hasLagrange, label, value, *subdm);CHKERRQ(ierr);
2453   }
2454   PetscFunctionReturn(0);
2455 }
2456 
2457 #undef __FUNCT__
2458 #define __FUNCT__ "DMPlexGetSubpointMap"
2459 /*@
2460   DMPlexGetSubpointMap - Returns a DMLabel with point dimension as values
2461 
2462   Input Parameter:
2463 . dm - The submesh DM
2464 
2465   Output Parameter:
2466 . subpointMap - The DMLabel of all the points from the original mesh in this submesh, or NULL if this is not a submesh
2467 
2468   Level: developer
2469 
2470 .seealso: DMPlexCreateSubmesh(), DMPlexCreateSubpointIS()
2471 @*/
2472 PetscErrorCode DMPlexGetSubpointMap(DM dm, DMLabel *subpointMap)
2473 {
2474   DM_Plex *mesh = (DM_Plex*) dm->data;
2475 
2476   PetscFunctionBegin;
2477   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2478   PetscValidPointer(subpointMap, 2);
2479   *subpointMap = mesh->subpointMap;
2480   PetscFunctionReturn(0);
2481 }
2482 
2483 #undef __FUNCT__
2484 #define __FUNCT__ "DMPlexSetSubpointMap"
2485 /* Note: Should normally not be called by the user, since it is set in DMPlexCreateSubmesh() */
2486 PetscErrorCode DMPlexSetSubpointMap(DM dm, DMLabel subpointMap)
2487 {
2488   DM_Plex       *mesh = (DM_Plex *) dm->data;
2489   DMLabel        tmp;
2490   PetscErrorCode ierr;
2491 
2492   PetscFunctionBegin;
2493   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2494   tmp  = mesh->subpointMap;
2495   mesh->subpointMap = subpointMap;
2496   ++mesh->subpointMap->refct;
2497   ierr = DMLabelDestroy(&tmp);CHKERRQ(ierr);
2498   PetscFunctionReturn(0);
2499 }
2500 
2501 #undef __FUNCT__
2502 #define __FUNCT__ "DMPlexCreateSubpointIS"
2503 /*@
2504   DMPlexCreateSubpointIS - Creates an IS covering the entire subdm chart with the original points as data
2505 
2506   Input Parameter:
2507 . dm - The submesh DM
2508 
2509   Output Parameter:
2510 . subpointIS - The IS of all the points from the original mesh in this submesh, or NULL if this is not a submesh
2511 
2512   Note: This is IS is guaranteed to be sorted by the construction of the submesh
2513 
2514   Level: developer
2515 
2516 .seealso: DMPlexCreateSubmesh(), DMPlexGetSubpointMap()
2517 @*/
2518 PetscErrorCode DMPlexCreateSubpointIS(DM dm, IS *subpointIS)
2519 {
2520   MPI_Comm        comm;
2521   DMLabel         subpointMap;
2522   IS              is;
2523   const PetscInt *opoints;
2524   PetscInt       *points, *depths;
2525   PetscInt        depth, depStart, depEnd, d, pStart, pEnd, p, n, off;
2526   PetscErrorCode  ierr;
2527 
2528   PetscFunctionBegin;
2529   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2530   PetscValidPointer(subpointIS, 2);
2531   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
2532   *subpointIS = NULL;
2533   ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr);
2534   if (subpointMap) {
2535     ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
2536     ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
2537     if (pStart) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Submeshes must start the point numbering at 0, not %d", pStart);
2538     ierr = DMGetWorkArray(dm, depth+1, PETSC_INT, &depths);CHKERRQ(ierr);
2539     depths[0] = depth;
2540     depths[1] = 0;
2541     for(d = 2; d <= depth; ++d) {depths[d] = depth+1 - d;}
2542     ierr = PetscMalloc(pEnd * sizeof(PetscInt), &points);CHKERRQ(ierr);
2543     for(d = 0, off = 0; d <= depth; ++d) {
2544       const PetscInt dep = depths[d];
2545 
2546       ierr = DMPlexGetDepthStratum(dm, dep, &depStart, &depEnd);CHKERRQ(ierr);
2547       ierr = DMLabelGetStratumSize(subpointMap, dep, &n);CHKERRQ(ierr);
2548       if (((d < 2) && (depth > 1)) || (d == 1)) { /* Only check vertices and cells for now since the map is broken for others */
2549         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);
2550       } else {
2551         if (!n) {
2552           if (d == 0) {
2553             /* Missing cells */
2554             for(p = 0; p < depEnd-depStart; ++p, ++off) points[off] = -1;
2555           } else {
2556             /* Missing faces */
2557             for(p = 0; p < depEnd-depStart; ++p, ++off) points[off] = PETSC_MAX_INT;
2558           }
2559         }
2560       }
2561       if (n) {
2562         ierr = DMLabelGetStratumIS(subpointMap, dep, &is);CHKERRQ(ierr);
2563         ierr = ISGetIndices(is, &opoints);CHKERRQ(ierr);
2564         for(p = 0; p < n; ++p, ++off) points[off] = opoints[p];
2565         ierr = ISRestoreIndices(is, &opoints);CHKERRQ(ierr);
2566         ierr = ISDestroy(&is);CHKERRQ(ierr);
2567       }
2568     }
2569     ierr = DMRestoreWorkArray(dm, depth+1, PETSC_INT, &depths);CHKERRQ(ierr);
2570     if (off != pEnd) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "The number of mapped submesh points %d should be %d", off, pEnd);
2571     ierr = ISCreateGeneral(comm, pEnd, points, PETSC_OWN_POINTER, subpointIS);CHKERRQ(ierr);
2572   }
2573   PetscFunctionReturn(0);
2574 }
2575