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