xref: /petsc/src/dm/impls/plex/plexsubmesh.c (revision 86bfadb514a0b5fe4e8d6e452466684f50283989)
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, pSize;
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   pSize = PetscMax(depth, dim) + 1;
1110   ierr = PetscMalloc3(pSize,PetscInt,&pStart,pSize,PetscInt,&pEnd,pSize,PetscInt,&pMax);CHKERRQ(ierr);
1111   ierr = DMPlexGetHybridBounds(dm, &pMax[depth], depth>1 ? &pMax[depth-1] : NULL, depth > 2 ? &pMax[1] : NULL, &pMax[0]);CHKERRQ(ierr);
1112   for (d = 0; d <= depth; ++d) {
1113     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr);
1114     if (pMax[d] >= 0) pEnd[d] = PetscMin(pEnd[d], pMax[d]);
1115   }
1116   /* Loop over initial vertices and mark all faces in the collective star() */
1117   ierr = DMLabelGetStratumIS(vertexLabel, value, &subvertexIS);CHKERRQ(ierr);
1118   if (subvertexIS) {
1119     ierr = ISGetSize(subvertexIS, &numSubVerticesInitial);CHKERRQ(ierr);
1120     ierr = ISGetIndices(subvertexIS, &subvertices);CHKERRQ(ierr);
1121   }
1122   for (v = 0; v < numSubVerticesInitial; ++v) {
1123     const PetscInt vertex = subvertices[v];
1124     PetscInt      *star   = NULL;
1125     PetscInt       starSize, s, numCells = 0, c;
1126 
1127     ierr = DMPlexGetTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
1128     for (s = 0; s < starSize*2; s += 2) {
1129       const PetscInt point = star[s];
1130       if ((point >= pStart[depth]) && (point < pEnd[depth])) star[numCells++] = point;
1131     }
1132     for (c = 0; c < numCells; ++c) {
1133       const PetscInt cell    = star[c];
1134       PetscInt      *closure = NULL;
1135       PetscInt       closureSize, cl;
1136       PetscInt       cellLoc, numCorners = 0, faceSize = 0;
1137 
1138       ierr = DMLabelGetValue(subpointMap, cell, &cellLoc);CHKERRQ(ierr);
1139       if (cellLoc == 2) continue;
1140       if (cellLoc >= 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Cell %d has dimension %d in the surface label", cell, cellLoc);
1141       ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1142       for (cl = 0; cl < closureSize*2; cl += 2) {
1143         const PetscInt point = closure[cl];
1144         PetscInt       vertexLoc;
1145 
1146         if ((point >= pStart[0]) && (point < pEnd[0])) {
1147           ++numCorners;
1148           ierr = DMLabelGetValue(vertexLabel, point, &vertexLoc);CHKERRQ(ierr);
1149           if (vertexLoc == value) closure[faceSize++] = point;
1150         }
1151       }
1152       if (!(*nFV)) {ierr = DMPlexGetNumFaceVertices(dm, dim, numCorners, nFV);CHKERRQ(ierr);}
1153       if (faceSize > *nFV) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid submesh: Too many vertices %d of an element on the surface", faceSize);
1154       if (faceSize == *nFV) {
1155         const PetscInt *cells = NULL;
1156         PetscInt        numCells, nc;
1157 
1158         ++(*numFaces);
1159         for (cl = 0; cl < faceSize; ++cl) {
1160           ierr = DMLabelSetValue(subpointMap, closure[cl], 0);CHKERRQ(ierr);
1161         }
1162         ierr = DMPlexGetJoin(dm, faceSize, closure, &numCells, &cells);CHKERRQ(ierr);
1163         for (nc = 0; nc < numCells; ++nc) {
1164           ierr = DMLabelSetValue(subpointMap, cells[nc], 2);CHKERRQ(ierr);
1165         }
1166         ierr = DMPlexRestoreJoin(dm, faceSize, closure, &numCells, &cells);CHKERRQ(ierr);
1167       }
1168       ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1169     }
1170     ierr = DMPlexRestoreTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
1171   }
1172   if (subvertexIS) {
1173     ierr = ISRestoreIndices(subvertexIS, &subvertices);CHKERRQ(ierr);
1174   }
1175   ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr);
1176   ierr = PetscFree3(pStart,pEnd,pMax);CHKERRQ(ierr);
1177   PetscFunctionReturn(0);
1178 }
1179 
1180 #undef __FUNCT__
1181 #define __FUNCT__ "DMPlexMarkSubmesh_Interpolated"
1182 static PetscErrorCode DMPlexMarkSubmesh_Interpolated(DM dm, DMLabel vertexLabel, PetscInt value, DMLabel subpointMap, DM subdm)
1183 {
1184   IS               subvertexIS;
1185   const PetscInt  *subvertices;
1186   PetscInt        *pStart, *pEnd, *pMax;
1187   PetscInt         dim, d, numSubVerticesInitial = 0, v;
1188   PetscErrorCode   ierr;
1189 
1190   PetscFunctionBegin;
1191   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1192   ierr = PetscMalloc3(dim+1,PetscInt,&pStart,dim+1,PetscInt,&pEnd,dim+1,PetscInt,&pMax);CHKERRQ(ierr);
1193   ierr = DMPlexGetHybridBounds(dm, &pMax[dim], dim>1 ? &pMax[dim-1] : NULL, dim > 2 ? &pMax[1] : NULL, &pMax[0]);CHKERRQ(ierr);
1194   for (d = 0; d <= dim; ++d) {
1195     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr);
1196     if (pMax[d] >= 0) pEnd[d] = PetscMin(pEnd[d], pMax[d]);
1197   }
1198   /* Loop over initial vertices and mark all faces in the collective star() */
1199   ierr = DMLabelGetStratumIS(vertexLabel, value, &subvertexIS);CHKERRQ(ierr);
1200   if (subvertexIS) {
1201     ierr = ISGetSize(subvertexIS, &numSubVerticesInitial);CHKERRQ(ierr);
1202     ierr = ISGetIndices(subvertexIS, &subvertices);CHKERRQ(ierr);
1203   }
1204   for (v = 0; v < numSubVerticesInitial; ++v) {
1205     const PetscInt vertex = subvertices[v];
1206     PetscInt      *star   = NULL;
1207     PetscInt       starSize, s, numFaces = 0, f;
1208 
1209     ierr = DMPlexGetTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
1210     for (s = 0; s < starSize*2; s += 2) {
1211       const PetscInt point = star[s];
1212       if ((point >= pStart[dim-1]) && (point < pEnd[dim-1])) star[numFaces++] = point;
1213     }
1214     for (f = 0; f < numFaces; ++f) {
1215       const PetscInt face    = star[f];
1216       PetscInt      *closure = NULL;
1217       PetscInt       closureSize, c;
1218       PetscInt       faceLoc;
1219 
1220       ierr = DMLabelGetValue(subpointMap, face, &faceLoc);CHKERRQ(ierr);
1221       if (faceLoc == dim-1) continue;
1222       if (faceLoc >= 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Face %d has dimension %d in the surface label", face, faceLoc);
1223       ierr = DMPlexGetTransitiveClosure(dm, face, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1224       for (c = 0; c < closureSize*2; c += 2) {
1225         const PetscInt point = closure[c];
1226         PetscInt       vertexLoc;
1227 
1228         if ((point >= pStart[0]) && (point < pEnd[0])) {
1229           ierr = DMLabelGetValue(vertexLabel, point, &vertexLoc);CHKERRQ(ierr);
1230           if (vertexLoc != value) break;
1231         }
1232       }
1233       if (c == closureSize*2) {
1234         const PetscInt *support;
1235         PetscInt        supportSize, s;
1236 
1237         for (c = 0; c < closureSize*2; c += 2) {
1238           const PetscInt point = closure[c];
1239 
1240           for (d = 0; d < dim; ++d) {
1241             if ((point >= pStart[d]) && (point < pEnd[d])) {
1242               ierr = DMLabelSetValue(subpointMap, point, d);CHKERRQ(ierr);
1243               break;
1244             }
1245           }
1246         }
1247         ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
1248         ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr);
1249         for (s = 0; s < supportSize; ++s) {
1250           ierr = DMLabelSetValue(subpointMap, support[s], dim);CHKERRQ(ierr);
1251         }
1252       }
1253       ierr = DMPlexRestoreTransitiveClosure(dm, face, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1254     }
1255     ierr = DMPlexRestoreTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
1256   }
1257   if (subvertexIS) {
1258     ierr = ISRestoreIndices(subvertexIS, &subvertices);CHKERRQ(ierr);
1259   }
1260   ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr);
1261   ierr = PetscFree3(pStart,pEnd,pMax);CHKERRQ(ierr);
1262   PetscFunctionReturn(0);
1263 }
1264 
1265 #undef __FUNCT__
1266 #define __FUNCT__ "DMPlexMarkCohesiveSubmesh_Uninterpolated"
1267 static PetscErrorCode DMPlexMarkCohesiveSubmesh_Uninterpolated(DM dm, PetscBool hasLagrange, DMLabel subpointMap, PetscInt *numFaces, PetscInt *nFV, PetscInt *subCells[], DM subdm)
1268 {
1269   const PetscInt *cone;
1270   PetscInt        dim, cMax, cEnd, c, p, coneSize;
1271   PetscErrorCode  ierr;
1272 
1273   PetscFunctionBegin;
1274   *numFaces = 0;
1275   *nFV = 0;
1276   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1277   ierr = DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);CHKERRQ(ierr);
1278   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
1279   if (cMax < 0) PetscFunctionReturn(0);
1280   ierr = DMPlexGetConeSize(dm, cMax, &coneSize);CHKERRQ(ierr);
1281   *numFaces = cEnd - cMax;
1282   *nFV      = hasLagrange ? coneSize/3 : coneSize/2;
1283   ierr = PetscMalloc(*numFaces *2 * sizeof(PetscInt), subCells);CHKERRQ(ierr);
1284   for (c = cMax; c < cEnd; ++c) {
1285     const PetscInt *cells;
1286     PetscInt        numCells;
1287 
1288     ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr);
1289     for (p = 0; p < *nFV; ++p) {
1290       ierr = DMLabelSetValue(subpointMap, cone[p], 0);CHKERRQ(ierr);
1291     }
1292     /* Negative face */
1293     ierr = DMPlexGetJoin(dm, *nFV, cone, &numCells, &cells);CHKERRQ(ierr);
1294     if (numCells != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells");
1295     for (p = 0; p < numCells; ++p) {
1296       ierr = DMLabelSetValue(subpointMap, cells[p], 2);CHKERRQ(ierr);
1297       (*subCells)[(c-cMax)*2+p] = cells[p];
1298     }
1299     ierr = DMPlexRestoreJoin(dm, *nFV, cone, &numCells, &cells);CHKERRQ(ierr);
1300     /* Positive face is not included */
1301   }
1302   PetscFunctionReturn(0);
1303 }
1304 
1305 #undef __FUNCT__
1306 #define __FUNCT__ "DMPlexMarkCohesiveSubmesh_Interpolated"
1307 static PetscErrorCode DMPlexMarkCohesiveSubmesh_Interpolated(DM dm, PetscBool hasLagrange, DMLabel subpointMap, DM subdm)
1308 {
1309   PetscInt      *pStart, *pEnd;
1310   PetscInt       dim, cMax, cEnd, c, d;
1311   PetscErrorCode ierr;
1312 
1313   PetscFunctionBegin;
1314   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1315   ierr = DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);CHKERRQ(ierr);
1316   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
1317   if (cMax < 0) PetscFunctionReturn(0);
1318   ierr = PetscMalloc2(dim+1,PetscInt,&pStart,dim+1,PetscInt,&pEnd);CHKERRQ(ierr);
1319   for (d = 0; d <= dim; ++d) {
1320     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr);
1321   }
1322   for (c = cMax; c < cEnd; ++c) {
1323     const PetscInt *cone;
1324     PetscInt       *closure = NULL;
1325     PetscInt        coneSize, closureSize, cl;
1326 
1327     ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
1328     if (hasLagrange) {
1329       if (coneSize != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells");
1330     } else {
1331       if (coneSize != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells");
1332     }
1333     /* Negative face */
1334     ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr);
1335     ierr = DMPlexGetTransitiveClosure(dm, cone[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1336     for (cl = 0; cl < closureSize*2; cl += 2) {
1337       const PetscInt point = closure[cl];
1338 
1339       for (d = 0; d <= dim; ++d) {
1340         if ((point >= pStart[d]) && (point < pEnd[d])) {
1341           ierr = DMLabelSetValue(subpointMap, point, d);CHKERRQ(ierr);
1342           break;
1343         }
1344       }
1345     }
1346     ierr = DMPlexRestoreTransitiveClosure(dm, cone[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1347     /* Cells -- positive face is not included */
1348     for (cl = 0; cl < 1; ++cl) {
1349       const PetscInt *support;
1350       PetscInt        supportSize, s;
1351 
1352       ierr = DMPlexGetSupportSize(dm, cone[cl], &supportSize);CHKERRQ(ierr);
1353       if (supportSize != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive faces should separate two cells");
1354       ierr = DMPlexGetSupport(dm, cone[cl], &support);CHKERRQ(ierr);
1355       for (s = 0; s < supportSize; ++s) {
1356         ierr = DMLabelSetValue(subpointMap, support[s], dim);CHKERRQ(ierr);
1357       }
1358     }
1359   }
1360   ierr = PetscFree2(pStart, pEnd);CHKERRQ(ierr);
1361   PetscFunctionReturn(0);
1362 }
1363 
1364 #undef __FUNCT__
1365 #define __FUNCT__ "DMPlexGetFaceOrientation"
1366 PetscErrorCode DMPlexGetFaceOrientation(DM dm, PetscInt cell, PetscInt numCorners, PetscInt indices[], PetscInt oppositeVertex, PetscInt origVertices[], PetscInt faceVertices[], PetscBool *posOriented)
1367 {
1368   MPI_Comm       comm;
1369   PetscBool      posOrient = PETSC_FALSE;
1370   const PetscInt debug     = 0;
1371   PetscInt       cellDim, faceSize, f;
1372   PetscErrorCode ierr;
1373 
1374   PetscFunctionBegin;
1375   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1376   ierr = DMPlexGetDimension(dm, &cellDim);CHKERRQ(ierr);
1377   if (debug) {PetscPrintf(comm, "cellDim: %d numCorners: %d\n", cellDim, numCorners);CHKERRQ(ierr);}
1378 
1379   if (cellDim == 1 && numCorners == 2) {
1380     /* Triangle */
1381     faceSize  = numCorners-1;
1382     posOrient = !(oppositeVertex%2) ? PETSC_TRUE : PETSC_FALSE;
1383   } else if (cellDim == 2 && numCorners == 3) {
1384     /* Triangle */
1385     faceSize  = numCorners-1;
1386     posOrient = !(oppositeVertex%2) ? PETSC_TRUE : PETSC_FALSE;
1387   } else if (cellDim == 3 && numCorners == 4) {
1388     /* Tetrahedron */
1389     faceSize  = numCorners-1;
1390     posOrient = (oppositeVertex%2) ? PETSC_TRUE : PETSC_FALSE;
1391   } else if (cellDim == 1 && numCorners == 3) {
1392     /* Quadratic line */
1393     faceSize  = 1;
1394     posOrient = PETSC_TRUE;
1395   } else if (cellDim == 2 && numCorners == 4) {
1396     /* Quads */
1397     faceSize = 2;
1398     if ((indices[1] > indices[0]) && (indices[1] - indices[0] == 1)) {
1399       posOrient = PETSC_TRUE;
1400     } else if ((indices[0] == 3) && (indices[1] == 0)) {
1401       posOrient = PETSC_TRUE;
1402     } else {
1403       if (((indices[0] > indices[1]) && (indices[0] - indices[1] == 1)) || ((indices[0] == 0) && (indices[1] == 3))) {
1404         posOrient = PETSC_FALSE;
1405       } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid quad crossedge");
1406     }
1407   } else if (cellDim == 2 && numCorners == 6) {
1408     /* Quadratic triangle (I hate this) */
1409     /* Edges are determined by the first 2 vertices (corners of edges) */
1410     const PetscInt faceSizeTri = 3;
1411     PetscInt       sortedIndices[3], i, iFace;
1412     PetscBool      found                    = PETSC_FALSE;
1413     PetscInt       faceVerticesTriSorted[9] = {
1414       0, 3,  4, /* bottom */
1415       1, 4,  5, /* right */
1416       2, 3,  5, /* left */
1417     };
1418     PetscInt       faceVerticesTri[9] = {
1419       0, 3,  4, /* bottom */
1420       1, 4,  5, /* right */
1421       2, 5,  3, /* left */
1422     };
1423 
1424     faceSize = faceSizeTri;
1425     for (i = 0; i < faceSizeTri; ++i) sortedIndices[i] = indices[i];
1426     ierr = PetscSortInt(faceSizeTri, sortedIndices);CHKERRQ(ierr);
1427     for (iFace = 0; iFace < 3; ++iFace) {
1428       const PetscInt ii = iFace*faceSizeTri;
1429       PetscInt       fVertex, cVertex;
1430 
1431       if ((sortedIndices[0] == faceVerticesTriSorted[ii+0]) &&
1432           (sortedIndices[1] == faceVerticesTriSorted[ii+1])) {
1433         for (fVertex = 0; fVertex < faceSizeTri; ++fVertex) {
1434           for (cVertex = 0; cVertex < faceSizeTri; ++cVertex) {
1435             if (indices[cVertex] == faceVerticesTri[ii+fVertex]) {
1436               faceVertices[fVertex] = origVertices[cVertex];
1437               break;
1438             }
1439           }
1440         }
1441         found = PETSC_TRUE;
1442         break;
1443       }
1444     }
1445     if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid tri crossface");
1446     if (posOriented) *posOriented = PETSC_TRUE;
1447     PetscFunctionReturn(0);
1448   } else if (cellDim == 2 && numCorners == 9) {
1449     /* Quadratic quad (I hate this) */
1450     /* Edges are determined by the first 2 vertices (corners of edges) */
1451     const PetscInt faceSizeQuad = 3;
1452     PetscInt       sortedIndices[3], i, iFace;
1453     PetscBool      found                      = PETSC_FALSE;
1454     PetscInt       faceVerticesQuadSorted[12] = {
1455       0, 1,  4, /* bottom */
1456       1, 2,  5, /* right */
1457       2, 3,  6, /* top */
1458       0, 3,  7, /* left */
1459     };
1460     PetscInt       faceVerticesQuad[12] = {
1461       0, 1,  4, /* bottom */
1462       1, 2,  5, /* right */
1463       2, 3,  6, /* top */
1464       3, 0,  7, /* left */
1465     };
1466 
1467     faceSize = faceSizeQuad;
1468     for (i = 0; i < faceSizeQuad; ++i) sortedIndices[i] = indices[i];
1469     ierr = PetscSortInt(faceSizeQuad, sortedIndices);CHKERRQ(ierr);
1470     for (iFace = 0; iFace < 4; ++iFace) {
1471       const PetscInt ii = iFace*faceSizeQuad;
1472       PetscInt       fVertex, cVertex;
1473 
1474       if ((sortedIndices[0] == faceVerticesQuadSorted[ii+0]) &&
1475           (sortedIndices[1] == faceVerticesQuadSorted[ii+1])) {
1476         for (fVertex = 0; fVertex < faceSizeQuad; ++fVertex) {
1477           for (cVertex = 0; cVertex < faceSizeQuad; ++cVertex) {
1478             if (indices[cVertex] == faceVerticesQuad[ii+fVertex]) {
1479               faceVertices[fVertex] = origVertices[cVertex];
1480               break;
1481             }
1482           }
1483         }
1484         found = PETSC_TRUE;
1485         break;
1486       }
1487     }
1488     if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid quad crossface");
1489     if (posOriented) *posOriented = PETSC_TRUE;
1490     PetscFunctionReturn(0);
1491   } else if (cellDim == 3 && numCorners == 8) {
1492     /* Hexes
1493        A hex is two oriented quads with the normal of the first
1494        pointing up at the second.
1495 
1496           7---6
1497          /|  /|
1498         4---5 |
1499         | 1-|-2
1500         |/  |/
1501         0---3
1502 
1503         Faces are determined by the first 4 vertices (corners of faces) */
1504     const PetscInt faceSizeHex = 4;
1505     PetscInt       sortedIndices[4], i, iFace;
1506     PetscBool      found                     = PETSC_FALSE;
1507     PetscInt       faceVerticesHexSorted[24] = {
1508       0, 1, 2, 3,  /* bottom */
1509       4, 5, 6, 7,  /* top */
1510       0, 3, 4, 5,  /* front */
1511       2, 3, 5, 6,  /* right */
1512       1, 2, 6, 7,  /* back */
1513       0, 1, 4, 7,  /* left */
1514     };
1515     PetscInt       faceVerticesHex[24] = {
1516       1, 2, 3, 0,  /* bottom */
1517       4, 5, 6, 7,  /* top */
1518       0, 3, 5, 4,  /* front */
1519       3, 2, 6, 5,  /* right */
1520       2, 1, 7, 6,  /* back */
1521       1, 0, 4, 7,  /* left */
1522     };
1523 
1524     faceSize = faceSizeHex;
1525     for (i = 0; i < faceSizeHex; ++i) sortedIndices[i] = indices[i];
1526     ierr = PetscSortInt(faceSizeHex, sortedIndices);CHKERRQ(ierr);
1527     for (iFace = 0; iFace < 6; ++iFace) {
1528       const PetscInt ii = iFace*faceSizeHex;
1529       PetscInt       fVertex, cVertex;
1530 
1531       if ((sortedIndices[0] == faceVerticesHexSorted[ii+0]) &&
1532           (sortedIndices[1] == faceVerticesHexSorted[ii+1]) &&
1533           (sortedIndices[2] == faceVerticesHexSorted[ii+2]) &&
1534           (sortedIndices[3] == faceVerticesHexSorted[ii+3])) {
1535         for (fVertex = 0; fVertex < faceSizeHex; ++fVertex) {
1536           for (cVertex = 0; cVertex < faceSizeHex; ++cVertex) {
1537             if (indices[cVertex] == faceVerticesHex[ii+fVertex]) {
1538               faceVertices[fVertex] = origVertices[cVertex];
1539               break;
1540             }
1541           }
1542         }
1543         found = PETSC_TRUE;
1544         break;
1545       }
1546     }
1547     if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid hex crossface");
1548     if (posOriented) *posOriented = PETSC_TRUE;
1549     PetscFunctionReturn(0);
1550   } else if (cellDim == 3 && numCorners == 10) {
1551     /* Quadratic tet */
1552     /* Faces are determined by the first 3 vertices (corners of faces) */
1553     const PetscInt faceSizeTet = 6;
1554     PetscInt       sortedIndices[6], i, iFace;
1555     PetscBool      found                     = PETSC_FALSE;
1556     PetscInt       faceVerticesTetSorted[24] = {
1557       0, 1, 2,  6, 7, 8, /* bottom */
1558       0, 3, 4,  6, 7, 9,  /* front */
1559       1, 4, 5,  7, 8, 9,  /* right */
1560       2, 3, 5,  6, 8, 9,  /* left */
1561     };
1562     PetscInt       faceVerticesTet[24] = {
1563       0, 1, 2,  6, 7, 8, /* bottom */
1564       0, 4, 3,  6, 7, 9,  /* front */
1565       1, 5, 4,  7, 8, 9,  /* right */
1566       2, 3, 5,  8, 6, 9,  /* left */
1567     };
1568 
1569     faceSize = faceSizeTet;
1570     for (i = 0; i < faceSizeTet; ++i) sortedIndices[i] = indices[i];
1571     ierr = PetscSortInt(faceSizeTet, sortedIndices);CHKERRQ(ierr);
1572     for (iFace=0; iFace < 4; ++iFace) {
1573       const PetscInt ii = iFace*faceSizeTet;
1574       PetscInt       fVertex, cVertex;
1575 
1576       if ((sortedIndices[0] == faceVerticesTetSorted[ii+0]) &&
1577           (sortedIndices[1] == faceVerticesTetSorted[ii+1]) &&
1578           (sortedIndices[2] == faceVerticesTetSorted[ii+2]) &&
1579           (sortedIndices[3] == faceVerticesTetSorted[ii+3])) {
1580         for (fVertex = 0; fVertex < faceSizeTet; ++fVertex) {
1581           for (cVertex = 0; cVertex < faceSizeTet; ++cVertex) {
1582             if (indices[cVertex] == faceVerticesTet[ii+fVertex]) {
1583               faceVertices[fVertex] = origVertices[cVertex];
1584               break;
1585             }
1586           }
1587         }
1588         found = PETSC_TRUE;
1589         break;
1590       }
1591     }
1592     if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid tet crossface");
1593     if (posOriented) *posOriented = PETSC_TRUE;
1594     PetscFunctionReturn(0);
1595   } else if (cellDim == 3 && numCorners == 27) {
1596     /* Quadratic hexes (I hate this)
1597        A hex is two oriented quads with the normal of the first
1598        pointing up at the second.
1599 
1600          7---6
1601         /|  /|
1602        4---5 |
1603        | 3-|-2
1604        |/  |/
1605        0---1
1606 
1607        Faces are determined by the first 4 vertices (corners of faces) */
1608     const PetscInt faceSizeQuadHex = 9;
1609     PetscInt       sortedIndices[9], i, iFace;
1610     PetscBool      found                         = PETSC_FALSE;
1611     PetscInt       faceVerticesQuadHexSorted[54] = {
1612       0, 1, 2, 3,  8, 9, 10, 11,  24, /* bottom */
1613       4, 5, 6, 7,  12, 13, 14, 15,  25, /* top */
1614       0, 1, 4, 5,  8, 12, 16, 17,  22, /* front */
1615       1, 2, 5, 6,  9, 13, 17, 18,  21, /* right */
1616       2, 3, 6, 7,  10, 14, 18, 19,  23, /* back */
1617       0, 3, 4, 7,  11, 15, 16, 19,  20, /* left */
1618     };
1619     PetscInt       faceVerticesQuadHex[54] = {
1620       3, 2, 1, 0,  10, 9, 8, 11,  24, /* bottom */
1621       4, 5, 6, 7,  12, 13, 14, 15,  25, /* top */
1622       0, 1, 5, 4,  8, 17, 12, 16,  22, /* front */
1623       1, 2, 6, 5,  9, 18, 13, 17,  21, /* right */
1624       2, 3, 7, 6,  10, 19, 14, 18,  23, /* back */
1625       3, 0, 4, 7,  11, 16, 15, 19,  20 /* left */
1626     };
1627 
1628     faceSize = faceSizeQuadHex;
1629     for (i = 0; i < faceSizeQuadHex; ++i) sortedIndices[i] = indices[i];
1630     ierr = PetscSortInt(faceSizeQuadHex, sortedIndices);CHKERRQ(ierr);
1631     for (iFace = 0; iFace < 6; ++iFace) {
1632       const PetscInt ii = iFace*faceSizeQuadHex;
1633       PetscInt       fVertex, cVertex;
1634 
1635       if ((sortedIndices[0] == faceVerticesQuadHexSorted[ii+0]) &&
1636           (sortedIndices[1] == faceVerticesQuadHexSorted[ii+1]) &&
1637           (sortedIndices[2] == faceVerticesQuadHexSorted[ii+2]) &&
1638           (sortedIndices[3] == faceVerticesQuadHexSorted[ii+3])) {
1639         for (fVertex = 0; fVertex < faceSizeQuadHex; ++fVertex) {
1640           for (cVertex = 0; cVertex < faceSizeQuadHex; ++cVertex) {
1641             if (indices[cVertex] == faceVerticesQuadHex[ii+fVertex]) {
1642               faceVertices[fVertex] = origVertices[cVertex];
1643               break;
1644             }
1645           }
1646         }
1647         found = PETSC_TRUE;
1648         break;
1649       }
1650     }
1651     if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid hex crossface");
1652     if (posOriented) *posOriented = PETSC_TRUE;
1653     PetscFunctionReturn(0);
1654   } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Unknown cell type for faceOrientation().");
1655   if (!posOrient) {
1656     if (debug) {ierr = PetscPrintf(comm, "  Reversing initial face orientation\n");CHKERRQ(ierr);}
1657     for (f = 0; f < faceSize; ++f) faceVertices[f] = origVertices[faceSize-1 - f];
1658   } else {
1659     if (debug) {ierr = PetscPrintf(comm, "  Keeping initial face orientation\n");CHKERRQ(ierr);}
1660     for (f = 0; f < faceSize; ++f) faceVertices[f] = origVertices[f];
1661   }
1662   if (posOriented) *posOriented = posOrient;
1663   PetscFunctionReturn(0);
1664 }
1665 
1666 #undef __FUNCT__
1667 #define __FUNCT__ "DMPlexGetOrientedFace"
1668 /*
1669     Given a cell and a face, as a set of vertices,
1670       return the oriented face, as a set of vertices, in faceVertices
1671     The orientation is such that the face normal points out of the cell
1672 */
1673 PetscErrorCode DMPlexGetOrientedFace(DM dm, PetscInt cell, PetscInt faceSize, const PetscInt face[], PetscInt numCorners, PetscInt indices[], PetscInt origVertices[], PetscInt faceVertices[], PetscBool *posOriented)
1674 {
1675   const PetscInt *cone = NULL;
1676   PetscInt        coneSize, v, f, v2;
1677   PetscInt        oppositeVertex = -1;
1678   PetscErrorCode  ierr;
1679 
1680   PetscFunctionBegin;
1681   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
1682   ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
1683   for (v = 0, v2 = 0; v < coneSize; ++v) {
1684     PetscBool found = PETSC_FALSE;
1685 
1686     for (f = 0; f < faceSize; ++f) {
1687       if (face[f] == cone[v]) {
1688         found = PETSC_TRUE; break;
1689       }
1690     }
1691     if (found) {
1692       indices[v2]      = v;
1693       origVertices[v2] = cone[v];
1694       ++v2;
1695     } else {
1696       oppositeVertex = v;
1697     }
1698   }
1699   ierr = DMPlexGetFaceOrientation(dm, cell, numCorners, indices, oppositeVertex, origVertices, faceVertices, posOriented);CHKERRQ(ierr);
1700   PetscFunctionReturn(0);
1701 }
1702 
1703 #undef __FUNCT__
1704 #define __FUNCT__ "DMPlexInsertFace_Internal"
1705 /*
1706   DMPlexInsertFace_Internal - Puts a face into the mesh
1707 
1708   Not collective
1709 
1710   Input Parameters:
1711   + dm              - The DMPlex
1712   . numFaceVertex   - The number of vertices in the face
1713   . faceVertices    - The vertices in the face for dm
1714   . subfaceVertices - The vertices in the face for subdm
1715   . numCorners      - The number of vertices in the cell
1716   . cell            - A cell in dm containing the face
1717   . subcell         - A cell in subdm containing the face
1718   . firstFace       - First face in the mesh
1719   - newFacePoint    - Next face in the mesh
1720 
1721   Output Parameters:
1722   . newFacePoint - Contains next face point number on input, updated on output
1723 
1724   Level: developer
1725 */
1726 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)
1727 {
1728   MPI_Comm        comm;
1729   DM_Plex        *submesh = (DM_Plex*) subdm->data;
1730   const PetscInt *faces;
1731   PetscInt        numFaces, coneSize;
1732   PetscErrorCode  ierr;
1733 
1734   PetscFunctionBegin;
1735   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1736   ierr = DMPlexGetConeSize(subdm, subcell, &coneSize);CHKERRQ(ierr);
1737   if (coneSize != 1) SETERRQ2(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cone size of cell %d is %d != 1", cell, coneSize);
1738 #if 0
1739   /* Cannot use this because support() has not been constructed yet */
1740   ierr = DMPlexGetJoin(subdm, numFaceVertices, subfaceVertices, &numFaces, &faces);CHKERRQ(ierr);
1741 #else
1742   {
1743     PetscInt f;
1744 
1745     numFaces = 0;
1746     ierr     = DMGetWorkArray(subdm, 1, PETSC_INT, (void **) &faces);CHKERRQ(ierr);
1747     for (f = firstFace; f < *newFacePoint; ++f) {
1748       PetscInt dof, off, d;
1749 
1750       ierr = PetscSectionGetDof(submesh->coneSection, f, &dof);CHKERRQ(ierr);
1751       ierr = PetscSectionGetOffset(submesh->coneSection, f, &off);CHKERRQ(ierr);
1752       /* Yes, I know this is quadratic, but I expect the sizes to be <5 */
1753       for (d = 0; d < dof; ++d) {
1754         const PetscInt p = submesh->cones[off+d];
1755         PetscInt       v;
1756 
1757         for (v = 0; v < numFaceVertices; ++v) {
1758           if (subfaceVertices[v] == p) break;
1759         }
1760         if (v == numFaceVertices) break;
1761       }
1762       if (d == dof) {
1763         numFaces               = 1;
1764         ((PetscInt*) faces)[0] = f;
1765       }
1766     }
1767   }
1768 #endif
1769   if (numFaces > 1) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Vertex set had %d faces, not one", numFaces);
1770   else if (numFaces == 1) {
1771     /* Add the other cell neighbor for this face */
1772     ierr = DMPlexSetCone(subdm, subcell, faces);CHKERRQ(ierr);
1773   } else {
1774     PetscInt *indices, *origVertices, *orientedVertices, *orientedSubVertices, v, ov;
1775     PetscBool posOriented;
1776 
1777     ierr                = DMGetWorkArray(subdm, 4*numFaceVertices * sizeof(PetscInt), PETSC_INT, &orientedVertices);CHKERRQ(ierr);
1778     origVertices        = &orientedVertices[numFaceVertices];
1779     indices             = &orientedVertices[numFaceVertices*2];
1780     orientedSubVertices = &orientedVertices[numFaceVertices*3];
1781     ierr                = DMPlexGetOrientedFace(dm, cell, numFaceVertices, faceVertices, numCorners, indices, origVertices, orientedVertices, &posOriented);CHKERRQ(ierr);
1782     /* TODO: I know that routine should return a permutation, not the indices */
1783     for (v = 0; v < numFaceVertices; ++v) {
1784       const PetscInt vertex = faceVertices[v], subvertex = subfaceVertices[v];
1785       for (ov = 0; ov < numFaceVertices; ++ov) {
1786         if (orientedVertices[ov] == vertex) {
1787           orientedSubVertices[ov] = subvertex;
1788           break;
1789         }
1790       }
1791       if (ov == numFaceVertices) SETERRQ1(comm, PETSC_ERR_PLIB, "Could not find face vertex %d in orientated set", vertex);
1792     }
1793     ierr = DMPlexSetCone(subdm, *newFacePoint, orientedSubVertices);CHKERRQ(ierr);
1794     ierr = DMPlexSetCone(subdm, subcell, newFacePoint);CHKERRQ(ierr);
1795     ierr = DMRestoreWorkArray(subdm, 4*numFaceVertices * sizeof(PetscInt), PETSC_INT, &orientedVertices);CHKERRQ(ierr);
1796     ++(*newFacePoint);
1797   }
1798 #if 0
1799   ierr = DMPlexRestoreJoin(subdm, numFaceVertices, subfaceVertices, &numFaces, &faces);CHKERRQ(ierr);
1800 #else
1801   ierr = DMRestoreWorkArray(subdm, 1, PETSC_INT, (void **) &faces);CHKERRQ(ierr);
1802 #endif
1803   PetscFunctionReturn(0);
1804 }
1805 
1806 #undef __FUNCT__
1807 #define __FUNCT__ "DMPlexCreateSubmesh_Uninterpolated"
1808 static PetscErrorCode DMPlexCreateSubmesh_Uninterpolated(DM dm, const char vertexLabelName[], PetscInt value, DM subdm)
1809 {
1810   MPI_Comm        comm;
1811   DMLabel         vertexLabel, subpointMap;
1812   IS              subvertexIS,  subcellIS;
1813   const PetscInt *subVertices, *subCells;
1814   PetscInt        numSubVertices, firstSubVertex, numSubCells;
1815   PetscInt       *subface, maxConeSize, numSubFaces, firstSubFace, newFacePoint, nFV;
1816   PetscInt        vStart, vEnd, c, f;
1817   PetscErrorCode  ierr;
1818 
1819   PetscFunctionBegin;
1820   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1821   /* Create subpointMap which marks the submesh */
1822   ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr);
1823   ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr);
1824   ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr);
1825   ierr = DMPlexGetLabel(dm, vertexLabelName, &vertexLabel);CHKERRQ(ierr);
1826   ierr = DMPlexMarkSubmesh_Uninterpolated(dm, vertexLabel, value, subpointMap, &numSubFaces, &nFV, subdm);CHKERRQ(ierr);
1827   /* Setup chart */
1828   ierr = DMLabelGetStratumSize(subpointMap, 0, &numSubVertices);CHKERRQ(ierr);
1829   ierr = DMLabelGetStratumSize(subpointMap, 2, &numSubCells);CHKERRQ(ierr);
1830   ierr = DMPlexSetChart(subdm, 0, numSubCells+numSubFaces+numSubVertices);CHKERRQ(ierr);
1831   ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr);
1832   /* Set cone sizes */
1833   firstSubVertex = numSubCells;
1834   firstSubFace   = numSubCells+numSubVertices;
1835   newFacePoint   = firstSubFace;
1836   ierr = DMLabelGetStratumIS(subpointMap, 0, &subvertexIS);CHKERRQ(ierr);
1837   if (subvertexIS) {ierr = ISGetIndices(subvertexIS, &subVertices);CHKERRQ(ierr);}
1838   ierr = DMLabelGetStratumIS(subpointMap, 2, &subcellIS);CHKERRQ(ierr);
1839   if (subcellIS) {ierr = ISGetIndices(subcellIS, &subCells);CHKERRQ(ierr);}
1840   for (c = 0; c < numSubCells; ++c) {
1841     ierr = DMPlexSetConeSize(subdm, c, 1);CHKERRQ(ierr);
1842   }
1843   for (f = firstSubFace; f < firstSubFace+numSubFaces; ++f) {
1844     ierr = DMPlexSetConeSize(subdm, f, nFV);CHKERRQ(ierr);
1845   }
1846   ierr = DMSetUp(subdm);CHKERRQ(ierr);
1847   /* Create face cones */
1848   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
1849   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
1850   ierr = DMGetWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr);
1851   for (c = 0; c < numSubCells; ++c) {
1852     const PetscInt cell    = subCells[c];
1853     const PetscInt subcell = c;
1854     PetscInt      *closure = NULL;
1855     PetscInt       closureSize, cl, numCorners = 0, faceSize = 0;
1856 
1857     ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1858     for (cl = 0; cl < closureSize*2; cl += 2) {
1859       const PetscInt point = closure[cl];
1860       PetscInt       subVertex;
1861 
1862       if ((point >= vStart) && (point < vEnd)) {
1863         ++numCorners;
1864         ierr = PetscFindInt(point, numSubVertices, subVertices, &subVertex);CHKERRQ(ierr);
1865         if (subVertex >= 0) {
1866           closure[faceSize] = point;
1867           subface[faceSize] = firstSubVertex+subVertex;
1868           ++faceSize;
1869         }
1870       }
1871     }
1872     if (faceSize > nFV) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Invalid submesh: Too many vertices %d of an element on the surface", faceSize);
1873     if (faceSize == nFV) {
1874       ierr = DMPlexInsertFace_Internal(dm, subdm, faceSize, closure, subface, numCorners, cell, subcell, firstSubFace, &newFacePoint);CHKERRQ(ierr);
1875     }
1876     ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1877   }
1878   ierr = DMRestoreWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr);
1879   ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr);
1880   ierr = DMPlexStratify(subdm);CHKERRQ(ierr);
1881   /* Build coordinates */
1882   {
1883     PetscSection coordSection, subCoordSection;
1884     Vec          coordinates, subCoordinates;
1885     PetscScalar *coords, *subCoords;
1886     PetscInt     numComp, coordSize, v;
1887 
1888     ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1889     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1890     ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr);
1891     ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr);
1892     ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr);
1893     ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr);
1894     ierr = PetscSectionSetChart(subCoordSection, firstSubVertex, firstSubVertex+numSubVertices);CHKERRQ(ierr);
1895     for (v = 0; v < numSubVertices; ++v) {
1896       const PetscInt vertex    = subVertices[v];
1897       const PetscInt subvertex = firstSubVertex+v;
1898       PetscInt       dof;
1899 
1900       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
1901       ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr);
1902       ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr);
1903     }
1904     ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr);
1905     ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr);
1906     ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr);
1907     ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
1908     ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr);
1909     if (coordSize) {
1910       ierr = VecGetArray(coordinates,    &coords);CHKERRQ(ierr);
1911       ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr);
1912       for (v = 0; v < numSubVertices; ++v) {
1913         const PetscInt vertex    = subVertices[v];
1914         const PetscInt subvertex = firstSubVertex+v;
1915         PetscInt       dof, off, sdof, soff, d;
1916 
1917         ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
1918         ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr);
1919         ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr);
1920         ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr);
1921         if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof);
1922         for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d];
1923       }
1924       ierr = VecRestoreArray(coordinates,    &coords);CHKERRQ(ierr);
1925       ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr);
1926     }
1927     ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr);
1928     ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr);
1929   }
1930   /* Cleanup */
1931   if (subvertexIS) {ierr = ISRestoreIndices(subvertexIS, &subVertices);CHKERRQ(ierr);}
1932   ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr);
1933   if (subcellIS) {ierr = ISRestoreIndices(subcellIS, &subCells);CHKERRQ(ierr);}
1934   ierr = ISDestroy(&subcellIS);CHKERRQ(ierr);
1935   PetscFunctionReturn(0);
1936 }
1937 
1938 #undef __FUNCT__
1939 #define __FUNCT__ "DMPlexCreateSubmesh_Interpolated"
1940 static PetscErrorCode DMPlexCreateSubmesh_Interpolated(DM dm, const char vertexLabelName[], PetscInt value, DM subdm)
1941 {
1942   MPI_Comm         comm;
1943   DMLabel          subpointMap, vertexLabel;
1944   IS              *subpointIS;
1945   const PetscInt **subpoints;
1946   PetscInt        *numSubPoints, *firstSubPoint, *coneNew;
1947   PetscInt         totSubPoints = 0, maxConeSize, dim, p, d, v;
1948   PetscErrorCode   ierr;
1949 
1950   PetscFunctionBegin;
1951   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1952   /* Create subpointMap which marks the submesh */
1953   ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr);
1954   ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr);
1955   ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr);
1956   ierr = DMPlexGetLabel(dm, vertexLabelName, &vertexLabel);CHKERRQ(ierr);
1957   ierr = DMPlexMarkSubmesh_Interpolated(dm, vertexLabel, value, subpointMap, subdm);CHKERRQ(ierr);
1958   /* Setup chart */
1959   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1960   ierr = PetscMalloc4(dim+1,PetscInt,&numSubPoints,dim+1,PetscInt,&firstSubPoint,dim+1,IS,&subpointIS,dim+1,const PetscInt *,&subpoints);CHKERRQ(ierr);
1961   for (d = 0; d <= dim; ++d) {
1962     ierr = DMLabelGetStratumSize(subpointMap, d, &numSubPoints[d]);CHKERRQ(ierr);
1963     totSubPoints += numSubPoints[d];
1964   }
1965   ierr = DMPlexSetChart(subdm, 0, totSubPoints);CHKERRQ(ierr);
1966   ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr);
1967   /* Set cone sizes */
1968   firstSubPoint[dim] = 0;
1969   firstSubPoint[0]   = firstSubPoint[dim] + numSubPoints[dim];
1970   if (dim > 1) {firstSubPoint[dim-1] = firstSubPoint[0]     + numSubPoints[0];}
1971   if (dim > 2) {firstSubPoint[dim-2] = firstSubPoint[dim-1] + numSubPoints[dim-1];}
1972   for (d = 0; d <= dim; ++d) {
1973     ierr = DMLabelGetStratumIS(subpointMap, d, &subpointIS[d]);CHKERRQ(ierr);
1974     ierr = ISGetIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr);
1975   }
1976   for (d = 0; d <= dim; ++d) {
1977     for (p = 0; p < numSubPoints[d]; ++p) {
1978       const PetscInt  point    = subpoints[d][p];
1979       const PetscInt  subpoint = firstSubPoint[d] + p;
1980       const PetscInt *cone;
1981       PetscInt        coneSize, coneSizeNew, c, val;
1982 
1983       ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
1984       ierr = DMPlexSetConeSize(subdm, subpoint, coneSize);CHKERRQ(ierr);
1985       if (d == dim) {
1986         ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
1987         for (c = 0, coneSizeNew = 0; c < coneSize; ++c) {
1988           ierr = DMLabelGetValue(subpointMap, cone[c], &val);CHKERRQ(ierr);
1989           if (val >= 0) coneSizeNew++;
1990         }
1991         ierr = DMPlexSetConeSize(subdm, subpoint, coneSizeNew);CHKERRQ(ierr);
1992       }
1993     }
1994   }
1995   ierr = DMSetUp(subdm);CHKERRQ(ierr);
1996   /* Set cones */
1997   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
1998   ierr = PetscMalloc(maxConeSize * sizeof(PetscInt), &coneNew);CHKERRQ(ierr);
1999   for (d = 0; d <= dim; ++d) {
2000     for (p = 0; p < numSubPoints[d]; ++p) {
2001       const PetscInt  point    = subpoints[d][p];
2002       const PetscInt  subpoint = firstSubPoint[d] + p;
2003       const PetscInt *cone;
2004       PetscInt        coneSize, subconeSize, coneSizeNew, c, subc;
2005 
2006       ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
2007       ierr = DMPlexGetConeSize(subdm, subpoint, &subconeSize);CHKERRQ(ierr);
2008       ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
2009       for (c = 0, coneSizeNew = 0; c < coneSize; ++c) {
2010         ierr = PetscFindInt(cone[c], numSubPoints[d-1], subpoints[d-1], &subc);CHKERRQ(ierr);
2011         if (subc >= 0) coneNew[coneSizeNew++] = firstSubPoint[d-1] + subc;
2012       }
2013       if (coneSizeNew != subconeSize) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of cone points located %d does not match subcone size %d", coneSizeNew, subconeSize);
2014       ierr = DMPlexSetCone(subdm, subpoint, coneNew);CHKERRQ(ierr);
2015     }
2016   }
2017   ierr = PetscFree(coneNew);CHKERRQ(ierr);
2018   ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr);
2019   ierr = DMPlexStratify(subdm);CHKERRQ(ierr);
2020   /* Build coordinates */
2021   {
2022     PetscSection coordSection, subCoordSection;
2023     Vec          coordinates, subCoordinates;
2024     PetscScalar *coords, *subCoords;
2025     PetscInt     numComp, coordSize;
2026 
2027     ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
2028     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
2029     ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr);
2030     ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr);
2031     ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr);
2032     ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr);
2033     ierr = PetscSectionSetChart(subCoordSection, firstSubPoint[0], firstSubPoint[0]+numSubPoints[0]);CHKERRQ(ierr);
2034     for (v = 0; v < numSubPoints[0]; ++v) {
2035       const PetscInt vertex    = subpoints[0][v];
2036       const PetscInt subvertex = firstSubPoint[0]+v;
2037       PetscInt       dof;
2038 
2039       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
2040       ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr);
2041       ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr);
2042     }
2043     ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr);
2044     ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr);
2045     ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr);
2046     ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
2047     ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr);
2048     ierr = VecGetArray(coordinates,    &coords);CHKERRQ(ierr);
2049     ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr);
2050     for (v = 0; v < numSubPoints[0]; ++v) {
2051       const PetscInt vertex    = subpoints[0][v];
2052       const PetscInt subvertex = firstSubPoint[0]+v;
2053       PetscInt dof, off, sdof, soff, d;
2054 
2055       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
2056       ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr);
2057       ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr);
2058       ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr);
2059       if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof);
2060       for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d];
2061     }
2062     ierr = VecRestoreArray(coordinates,    &coords);CHKERRQ(ierr);
2063     ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr);
2064     ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr);
2065     ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr);
2066   }
2067   /* Cleanup */
2068   for (d = 0; d <= dim; ++d) {
2069     ierr = ISRestoreIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr);
2070     ierr = ISDestroy(&subpointIS[d]);CHKERRQ(ierr);
2071   }
2072   ierr = PetscFree4(numSubPoints,firstSubPoint,subpointIS,subpoints);CHKERRQ(ierr);
2073   PetscFunctionReturn(0);
2074 }
2075 
2076 #undef __FUNCT__
2077 #define __FUNCT__ "DMPlexCreateSubmesh"
2078 /*@C
2079   DMPlexCreateSubmesh - Extract a hypersurface from the mesh using vertices defined by a label
2080 
2081   Input Parameters:
2082 + dm           - The original mesh
2083 . vertexLabel  - The DMLabel marking vertices contained in the surface
2084 - value        - The label value to use
2085 
2086   Output Parameter:
2087 . subdm - The surface mesh
2088 
2089   Note: This function produces a DMLabel mapping original points in the submesh to their depth. This can be obtained using DMPlexGetSubpointMap().
2090 
2091   Level: developer
2092 
2093 .seealso: DMPlexGetSubpointMap(), DMPlexGetLabel(), DMLabelSetValue()
2094 @*/
2095 PetscErrorCode DMPlexCreateSubmesh(DM dm, const char vertexLabel[], PetscInt value, DM *subdm)
2096 {
2097   PetscInt       dim, depth;
2098   PetscErrorCode ierr;
2099 
2100   PetscFunctionBegin;
2101   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2102   PetscValidCharPointer(vertexLabel, 2);
2103   PetscValidPointer(subdm, 3);
2104   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
2105   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
2106   ierr = DMCreate(PetscObjectComm((PetscObject)dm), subdm);CHKERRQ(ierr);
2107   ierr = DMSetType(*subdm, DMPLEX);CHKERRQ(ierr);
2108   ierr = DMPlexSetDimension(*subdm, dim-1);CHKERRQ(ierr);
2109   if (depth == dim) {
2110     ierr = DMPlexCreateSubmesh_Interpolated(dm, vertexLabel, value, *subdm);CHKERRQ(ierr);
2111   } else {
2112     ierr = DMPlexCreateSubmesh_Uninterpolated(dm, vertexLabel, value, *subdm);CHKERRQ(ierr);
2113   }
2114   PetscFunctionReturn(0);
2115 }
2116 
2117 #undef __FUNCT__
2118 #define __FUNCT__ "DMPlexCreateCohesiveSubmesh_Uninterpolated"
2119 static PetscErrorCode DMPlexCreateCohesiveSubmesh_Uninterpolated(DM dm, PetscBool hasLagrange, DM subdm)
2120 {
2121   MPI_Comm        comm;
2122   DMLabel         subpointMap;
2123   IS              subvertexIS;
2124   const PetscInt *subVertices;
2125   PetscInt        numSubVertices, firstSubVertex, numSubCells, *subCells;
2126   PetscInt       *subface, maxConeSize, numSubFaces, firstSubFace, newFacePoint, nFV;
2127   PetscInt        cMax, c, f;
2128   PetscErrorCode  ierr;
2129 
2130   PetscFunctionBegin;
2131   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
2132   /* Create subpointMap which marks the submesh */
2133   ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr);
2134   ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr);
2135   ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr);
2136   ierr = DMPlexMarkCohesiveSubmesh_Uninterpolated(dm, hasLagrange, subpointMap, &numSubFaces, &nFV, &subCells, subdm);CHKERRQ(ierr);
2137   /* Setup chart */
2138   ierr = DMLabelGetStratumSize(subpointMap, 0, &numSubVertices);CHKERRQ(ierr);
2139   ierr = DMLabelGetStratumSize(subpointMap, 2, &numSubCells);CHKERRQ(ierr);
2140   ierr = DMPlexSetChart(subdm, 0, numSubCells+numSubFaces+numSubVertices);CHKERRQ(ierr);
2141   ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr);
2142   /* Set cone sizes */
2143   firstSubVertex = numSubCells;
2144   firstSubFace   = numSubCells+numSubVertices;
2145   newFacePoint   = firstSubFace;
2146   ierr = DMLabelGetStratumIS(subpointMap, 0, &subvertexIS);CHKERRQ(ierr);
2147   if (subvertexIS) {ierr = ISGetIndices(subvertexIS, &subVertices);CHKERRQ(ierr);}
2148   for (c = 0; c < numSubCells; ++c) {
2149     ierr = DMPlexSetConeSize(subdm, c, 1);CHKERRQ(ierr);
2150   }
2151   for (f = firstSubFace; f < firstSubFace+numSubFaces; ++f) {
2152     ierr = DMPlexSetConeSize(subdm, f, nFV);CHKERRQ(ierr);
2153   }
2154   ierr = DMSetUp(subdm);CHKERRQ(ierr);
2155   /* Create face cones */
2156   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
2157   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
2158   ierr = DMGetWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr);
2159   for (c = 0; c < numSubCells; ++c) {
2160     const PetscInt  cell    = subCells[c];
2161     const PetscInt  subcell = c;
2162     const PetscInt *cone, *cells;
2163     PetscInt        numCells, subVertex, p, v;
2164 
2165     if (cell < cMax) continue;
2166     ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
2167     for (v = 0; v < nFV; ++v) {
2168       ierr = PetscFindInt(cone[v], numSubVertices, subVertices, &subVertex);CHKERRQ(ierr);
2169       subface[v] = firstSubVertex+subVertex;
2170     }
2171     ierr = DMPlexSetCone(subdm, newFacePoint, subface);CHKERRQ(ierr);
2172     ierr = DMPlexSetCone(subdm, subcell, &newFacePoint);CHKERRQ(ierr);
2173     ierr = DMPlexGetJoin(dm, nFV, cone, &numCells, &cells);CHKERRQ(ierr);
2174     if (numCells != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells");
2175     for (p = 0; p < numCells; ++p) {
2176       PetscInt negsubcell;
2177 
2178       if (cells[p] >= cMax) continue;
2179       /* I know this is a crap search */
2180       for (negsubcell = 0; negsubcell < numSubCells; ++negsubcell) {
2181         if (subCells[negsubcell] == cells[p]) break;
2182       }
2183       if (negsubcell == numSubCells) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find negative face neighbor for cohesive cell %d", cell);
2184       ierr = DMPlexSetCone(subdm, negsubcell, &newFacePoint);CHKERRQ(ierr);
2185     }
2186     ierr = DMPlexRestoreJoin(dm, nFV, cone, &numCells, &cells);CHKERRQ(ierr);
2187     ++newFacePoint;
2188   }
2189   ierr = DMRestoreWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr);
2190   ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr);
2191   ierr = DMPlexStratify(subdm);CHKERRQ(ierr);
2192   /* Build coordinates */
2193   {
2194     PetscSection coordSection, subCoordSection;
2195     Vec          coordinates, subCoordinates;
2196     PetscScalar *coords, *subCoords;
2197     PetscInt     numComp, coordSize, v;
2198 
2199     ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
2200     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
2201     ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr);
2202     ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr);
2203     ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr);
2204     ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr);
2205     ierr = PetscSectionSetChart(subCoordSection, firstSubVertex, firstSubVertex+numSubVertices);CHKERRQ(ierr);
2206     for (v = 0; v < numSubVertices; ++v) {
2207       const PetscInt vertex    = subVertices[v];
2208       const PetscInt subvertex = firstSubVertex+v;
2209       PetscInt       dof;
2210 
2211       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
2212       ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr);
2213       ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr);
2214     }
2215     ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr);
2216     ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr);
2217     ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr);
2218     if (coordSize) {
2219       ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
2220       ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr);
2221       ierr = VecGetArray(coordinates,    &coords);CHKERRQ(ierr);
2222       ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr);
2223       for (v = 0; v < numSubVertices; ++v) {
2224         const PetscInt vertex    = subVertices[v];
2225         const PetscInt subvertex = firstSubVertex+v;
2226         PetscInt       dof, off, sdof, soff, d;
2227 
2228         ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
2229         ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr);
2230         ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr);
2231         ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr);
2232         if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof);
2233         for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d];
2234       }
2235       ierr = VecRestoreArray(coordinates,    &coords);CHKERRQ(ierr);
2236       ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr);
2237     }
2238     ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr);
2239     ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr);
2240   }
2241   /* Cleanup */
2242   if (subvertexIS) {ierr = ISRestoreIndices(subvertexIS, &subVertices);CHKERRQ(ierr);}
2243   ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr);
2244   ierr = PetscFree(subCells);CHKERRQ(ierr);
2245   PetscFunctionReturn(0);
2246 }
2247 
2248 #undef __FUNCT__
2249 #define __FUNCT__ "DMPlexCreateCohesiveSubmesh_Interpolated"
2250 static PetscErrorCode DMPlexCreateCohesiveSubmesh_Interpolated(DM dm, PetscBool hasLagrange, DM subdm)
2251 {
2252   MPI_Comm         comm;
2253   DMLabel          subpointMap;
2254   IS              *subpointIS;
2255   const PetscInt **subpoints;
2256   PetscInt        *numSubPoints, *firstSubPoint, *coneNew;
2257   PetscInt         totSubPoints = 0, maxConeSize, dim, p, d, v;
2258   PetscErrorCode   ierr;
2259 
2260   PetscFunctionBegin;
2261   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
2262   /* Create subpointMap which marks the submesh */
2263   ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr);
2264   ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr);
2265   ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr);
2266   ierr = DMPlexMarkCohesiveSubmesh_Interpolated(dm, hasLagrange, subpointMap, subdm);CHKERRQ(ierr);
2267   /* Setup chart */
2268   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
2269   ierr = PetscMalloc4(dim+1,PetscInt,&numSubPoints,dim+1,PetscInt,&firstSubPoint,dim+1,IS,&subpointIS,dim+1,const PetscInt *,&subpoints);CHKERRQ(ierr);
2270   for (d = 0; d <= dim; ++d) {
2271     ierr = DMLabelGetStratumSize(subpointMap, d, &numSubPoints[d]);CHKERRQ(ierr);
2272     totSubPoints += numSubPoints[d];
2273   }
2274   ierr = DMPlexSetChart(subdm, 0, totSubPoints);CHKERRQ(ierr);
2275   ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr);
2276   /* Set cone sizes */
2277   firstSubPoint[dim] = 0;
2278   firstSubPoint[0]   = firstSubPoint[dim] + numSubPoints[dim];
2279   if (dim > 1) {firstSubPoint[dim-1] = firstSubPoint[0]     + numSubPoints[0];}
2280   if (dim > 2) {firstSubPoint[dim-2] = firstSubPoint[dim-1] + numSubPoints[dim-1];}
2281   for (d = 0; d <= dim; ++d) {
2282     ierr = DMLabelGetStratumIS(subpointMap, d, &subpointIS[d]);CHKERRQ(ierr);
2283     ierr = ISGetIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr);
2284   }
2285   for (d = 0; d <= dim; ++d) {
2286     for (p = 0; p < numSubPoints[d]; ++p) {
2287       const PetscInt  point    = subpoints[d][p];
2288       const PetscInt  subpoint = firstSubPoint[d] + p;
2289       const PetscInt *cone;
2290       PetscInt        coneSize, coneSizeNew, c, val;
2291 
2292       ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
2293       ierr = DMPlexSetConeSize(subdm, subpoint, coneSize);CHKERRQ(ierr);
2294       if (d == dim) {
2295         ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
2296         for (c = 0, coneSizeNew = 0; c < coneSize; ++c) {
2297           ierr = DMLabelGetValue(subpointMap, cone[c], &val);CHKERRQ(ierr);
2298           if (val >= 0) coneSizeNew++;
2299         }
2300         ierr = DMPlexSetConeSize(subdm, subpoint, coneSizeNew);CHKERRQ(ierr);
2301       }
2302     }
2303   }
2304   ierr = DMSetUp(subdm);CHKERRQ(ierr);
2305   /* Set cones */
2306   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
2307   ierr = PetscMalloc(maxConeSize * sizeof(PetscInt), &coneNew);CHKERRQ(ierr);
2308   for (d = 0; d <= dim; ++d) {
2309     for (p = 0; p < numSubPoints[d]; ++p) {
2310       const PetscInt  point    = subpoints[d][p];
2311       const PetscInt  subpoint = firstSubPoint[d] + p;
2312       const PetscInt *cone;
2313       PetscInt        coneSize, subconeSize, coneSizeNew, c, subc;
2314 
2315       ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
2316       ierr = DMPlexGetConeSize(subdm, subpoint, &subconeSize);CHKERRQ(ierr);
2317       ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
2318       for (c = 0, coneSizeNew = 0; c < coneSize; ++c) {
2319         ierr = PetscFindInt(cone[c], numSubPoints[d-1], subpoints[d-1], &subc);CHKERRQ(ierr);
2320         if (subc >= 0) coneNew[coneSizeNew++] = firstSubPoint[d-1] + subc;
2321       }
2322       if (coneSizeNew != subconeSize) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of cone points located %d does not match subcone size %d", coneSizeNew, subconeSize);
2323       ierr = DMPlexSetCone(subdm, subpoint, coneNew);CHKERRQ(ierr);
2324     }
2325   }
2326   ierr = PetscFree(coneNew);CHKERRQ(ierr);
2327   ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr);
2328   ierr = DMPlexStratify(subdm);CHKERRQ(ierr);
2329   /* Build coordinates */
2330   {
2331     PetscSection coordSection, subCoordSection;
2332     Vec          coordinates, subCoordinates;
2333     PetscScalar *coords, *subCoords;
2334     PetscInt     numComp, coordSize;
2335 
2336     ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
2337     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
2338     ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr);
2339     ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr);
2340     ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr);
2341     ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr);
2342     ierr = PetscSectionSetChart(subCoordSection, firstSubPoint[0], firstSubPoint[0]+numSubPoints[0]);CHKERRQ(ierr);
2343     for (v = 0; v < numSubPoints[0]; ++v) {
2344       const PetscInt vertex    = subpoints[0][v];
2345       const PetscInt subvertex = firstSubPoint[0]+v;
2346       PetscInt       dof;
2347 
2348       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
2349       ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr);
2350       ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr);
2351     }
2352     ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr);
2353     ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr);
2354     ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr);
2355     ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
2356     ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr);
2357     ierr = VecGetArray(coordinates,    &coords);CHKERRQ(ierr);
2358     ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr);
2359     for (v = 0; v < numSubPoints[0]; ++v) {
2360       const PetscInt vertex    = subpoints[0][v];
2361       const PetscInt subvertex = firstSubPoint[0]+v;
2362       PetscInt dof, off, sdof, soff, d;
2363 
2364       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
2365       ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr);
2366       ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr);
2367       ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr);
2368       if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof);
2369       for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d];
2370     }
2371     ierr = VecRestoreArray(coordinates,    &coords);CHKERRQ(ierr);
2372     ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr);
2373     ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr);
2374     ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr);
2375   }
2376   /* Cleanup */
2377   for (d = 0; d <= dim; ++d) {
2378     ierr = ISRestoreIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr);
2379     ierr = ISDestroy(&subpointIS[d]);CHKERRQ(ierr);
2380   }
2381   ierr = PetscFree4(numSubPoints,firstSubPoint,subpointIS,subpoints);CHKERRQ(ierr);
2382   PetscFunctionReturn(0);
2383 }
2384 
2385 #undef __FUNCT__
2386 #define __FUNCT__ "DMPlexCreateCohesiveSubmesh"
2387 /*
2388   DMPlexCreateCohesiveSubmesh - Extract from a mesh with cohesive cells the hypersurface defined by one face of the cells.
2389 
2390   Input Parameters:
2391 + dm          - The original mesh
2392 - hasLagrange - The mesh has Lagrange unknowns in the cohesive cells
2393 
2394   Output Parameter:
2395 . subdm - The surface mesh
2396 
2397   Note: This function produces a DMLabel mapping original points in the submesh to their depth. This can be obtained using DMPlexGetSubpointMap().
2398 
2399   Level: developer
2400 
2401 .seealso: DMPlexGetSubpointMap(), DMPlexCreateSubmesh()
2402 */
2403 PetscErrorCode DMPlexCreateCohesiveSubmesh(DM dm, PetscBool hasLagrange, DM *subdm)
2404 {
2405   PetscInt       dim, depth;
2406   PetscErrorCode ierr;
2407 
2408   PetscFunctionBegin;
2409   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2410   PetscValidPointer(subdm, 3);
2411   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
2412   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
2413   ierr = DMCreate(PetscObjectComm((PetscObject)dm), subdm);CHKERRQ(ierr);
2414   ierr = DMSetType(*subdm, DMPLEX);CHKERRQ(ierr);
2415   ierr = DMPlexSetDimension(*subdm, dim-1);CHKERRQ(ierr);
2416   if (depth == dim) {
2417     ierr = DMPlexCreateCohesiveSubmesh_Interpolated(dm, hasLagrange, *subdm);CHKERRQ(ierr);
2418   } else {
2419     ierr = DMPlexCreateCohesiveSubmesh_Uninterpolated(dm, hasLagrange, *subdm);CHKERRQ(ierr);
2420   }
2421   PetscFunctionReturn(0);
2422 }
2423 
2424 #undef __FUNCT__
2425 #define __FUNCT__ "DMPlexGetSubpointMap"
2426 /*@
2427   DMPlexGetSubpointMap - Returns a DMLabel with point dimension as values
2428 
2429   Input Parameter:
2430 . dm - The submesh DM
2431 
2432   Output Parameter:
2433 . subpointMap - The DMLabel of all the points from the original mesh in this submesh, or NULL if this is not a submesh
2434 
2435   Level: developer
2436 
2437 .seealso: DMPlexCreateSubmesh(), DMPlexCreateSubpointIS()
2438 @*/
2439 PetscErrorCode DMPlexGetSubpointMap(DM dm, DMLabel *subpointMap)
2440 {
2441   DM_Plex *mesh = (DM_Plex*) dm->data;
2442 
2443   PetscFunctionBegin;
2444   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2445   PetscValidPointer(subpointMap, 2);
2446   *subpointMap = mesh->subpointMap;
2447   PetscFunctionReturn(0);
2448 }
2449 
2450 #undef __FUNCT__
2451 #define __FUNCT__ "DMPlexSetSubpointMap"
2452 /* Note: Should normally not be called by the user, since it is set in DMPlexCreateSubmesh() */
2453 PetscErrorCode DMPlexSetSubpointMap(DM dm, DMLabel subpointMap)
2454 {
2455   DM_Plex       *mesh = (DM_Plex *) dm->data;
2456   DMLabel        tmp;
2457   PetscErrorCode ierr;
2458 
2459   PetscFunctionBegin;
2460   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2461   tmp  = mesh->subpointMap;
2462   mesh->subpointMap = subpointMap;
2463   ++mesh->subpointMap->refct;
2464   ierr = DMLabelDestroy(&tmp);CHKERRQ(ierr);
2465   PetscFunctionReturn(0);
2466 }
2467 
2468 #undef __FUNCT__
2469 #define __FUNCT__ "DMPlexCreateSubpointIS"
2470 /*@
2471   DMPlexCreateSubpointIS - Creates an IS covering the entire subdm chart with the original points as data
2472 
2473   Input Parameter:
2474 . dm - The submesh DM
2475 
2476   Output Parameter:
2477 . subpointIS - The IS of all the points from the original mesh in this submesh, or NULL if this is not a submesh
2478 
2479   Note: This is IS is guaranteed to be sorted by the construction of the submesh
2480 
2481   Level: developer
2482 
2483 .seealso: DMPlexCreateSubmesh(), DMPlexGetSubpointMap()
2484 @*/
2485 PetscErrorCode DMPlexCreateSubpointIS(DM dm, IS *subpointIS)
2486 {
2487   MPI_Comm        comm;
2488   DMLabel         subpointMap;
2489   IS              is;
2490   const PetscInt *opoints;
2491   PetscInt       *points, *depths;
2492   PetscInt        depth, depStart, depEnd, d, pStart, pEnd, p, n, off;
2493   PetscErrorCode  ierr;
2494 
2495   PetscFunctionBegin;
2496   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2497   PetscValidPointer(subpointIS, 2);
2498   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
2499   *subpointIS = NULL;
2500   ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr);
2501   if (subpointMap) {
2502     ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
2503     ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
2504     if (pStart) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Submeshes must start the point numbering at 0, not %d", pStart);
2505     ierr = DMGetWorkArray(dm, depth+1, PETSC_INT, &depths);CHKERRQ(ierr);
2506     depths[0] = depth;
2507     depths[1] = 0;
2508     for(d = 2; d <= depth; ++d) {depths[d] = depth+1 - d;}
2509     ierr = PetscMalloc(pEnd * sizeof(PetscInt), &points);CHKERRQ(ierr);
2510     for(d = 0, off = 0; d <= depth; ++d) {
2511       const PetscInt dep = depths[d];
2512 
2513       ierr = DMPlexGetDepthStratum(dm, dep, &depStart, &depEnd);CHKERRQ(ierr);
2514       ierr = DMLabelGetStratumSize(subpointMap, dep, &n);CHKERRQ(ierr);
2515       if (((d < 2) && (depth > 1)) || (d == 1)) { /* Only check vertices and cells for now since the map is broken for others */
2516         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);
2517       } else {
2518         if (!n) {
2519           if (d == 0) {
2520             /* Missing cells */
2521             for(p = 0; p < depEnd-depStart; ++p, ++off) points[off] = -1;
2522           } else {
2523             /* Missing faces */
2524             for(p = 0; p < depEnd-depStart; ++p, ++off) points[off] = PETSC_MAX_INT;
2525           }
2526         }
2527       }
2528       if (n) {
2529         ierr = DMLabelGetStratumIS(subpointMap, dep, &is);CHKERRQ(ierr);
2530         ierr = ISGetIndices(is, &opoints);CHKERRQ(ierr);
2531         for(p = 0; p < n; ++p, ++off) points[off] = opoints[p];
2532         ierr = ISRestoreIndices(is, &opoints);CHKERRQ(ierr);
2533         ierr = ISDestroy(&is);CHKERRQ(ierr);
2534       }
2535     }
2536     ierr = DMRestoreWorkArray(dm, depth+1, PETSC_INT, &depths);CHKERRQ(ierr);
2537     if (off != pEnd) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "The number of mapped submesh points %d should be %d", off, pEnd);
2538     ierr = ISCreateGeneral(comm, pEnd, points, PETSC_OWN_POINTER, subpointIS);CHKERRQ(ierr);
2539   }
2540   PetscFunctionReturn(0);
2541 }
2542