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