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