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