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