xref: /petsc/src/dm/impls/plex/plexinterpolate.c (revision e9fa77a147e78afc60c94377e66a31d6beab20bf)
1 #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 #include <petsc/private/hashmapi.h>
3 #include <petsc/private/hashmapij.h>
4 
5 /* HashIJKL */
6 
7 #include <petsc/private/hashmap.h>
8 
9 typedef struct _PetscHashIJKLKey { PetscInt i, j, k, l; } PetscHashIJKLKey;
10 
11 #define PetscHashIJKLKeyHash(key) \
12   PetscHashCombine(PetscHashCombine(PetscHashInt((key).i),PetscHashInt((key).j)), \
13                    PetscHashCombine(PetscHashInt((key).k),PetscHashInt((key).l)))
14 
15 #define PetscHashIJKLKeyEqual(k1,k2) \
16   (((k1).i==(k2).i) ? ((k1).j==(k2).j) ? ((k1).k==(k2).k) ? ((k1).l==(k2).l) : 0 : 0 : 0)
17 
18 PETSC_HASH_MAP(HashIJKL, PetscHashIJKLKey, PetscInt, PetscHashIJKLKeyHash, PetscHashIJKLKeyEqual, -1)
19 
20 
21 /*
22   DMPlexGetFaces_Internal - Gets groups of vertices that correspond to faces for the given cell
23   This assumes that the mesh is not interpolated from the depth of point p to the vertices
24 */
25 PetscErrorCode DMPlexGetFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[])
26 {
27   const PetscInt *cone = NULL;
28   PetscInt        coneSize;
29   PetscErrorCode  ierr;
30 
31   PetscFunctionBegin;
32   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
33   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
34   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
35   ierr = DMPlexGetRawFaces_Internal(dm, dim, coneSize, cone, numFaces, faceSize, faces);CHKERRQ(ierr);
36   PetscFunctionReturn(0);
37 }
38 
39 /*
40   DMPlexRestoreFaces_Internal - Restores the array
41 */
42 PetscErrorCode DMPlexRestoreFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[])
43 {
44   PetscErrorCode  ierr;
45 
46   PetscFunctionBegin;
47   if (faces) { ierr = DMRestoreWorkArray(dm, 0, MPIU_INT, (void *) faces);CHKERRQ(ierr); }
48   PetscFunctionReturn(0);
49 }
50 
51 /*
52   DMPlexGetRawFaces_Internal - Gets groups of vertices that correspond to faces for the given cone
53 */
54 PetscErrorCode DMPlexGetRawFaces_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[])
55 {
56   PetscInt       *facesTmp;
57   PetscInt        maxConeSize, maxSupportSize;
58   PetscErrorCode  ierr;
59 
60   PetscFunctionBegin;
61   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
62   if (faces && coneSize) PetscValidIntPointer(cone,4);
63   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
64   if (faces) {ierr = DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), MPIU_INT, &facesTmp);CHKERRQ(ierr);}
65   switch (dim) {
66   case 1:
67     switch (coneSize) {
68     case 2:
69       if (faces) {
70         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
71         *faces = facesTmp;
72       }
73       if (numFaces) *numFaces = 2;
74       if (faceSize) *faceSize = 1;
75       break;
76     default:
77       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
78     }
79     break;
80   case 2:
81     switch (coneSize) {
82     case 3:
83       if (faces) {
84         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
85         facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
86         facesTmp[4] = cone[2]; facesTmp[5] = cone[0];
87         *faces = facesTmp;
88       }
89       if (numFaces) *numFaces = 3;
90       if (faceSize) *faceSize = 2;
91       break;
92     case 4:
93       /* Vertices follow right hand rule */
94       if (faces) {
95         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
96         facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
97         facesTmp[4] = cone[2]; facesTmp[5] = cone[3];
98         facesTmp[6] = cone[3]; facesTmp[7] = cone[0];
99         *faces = facesTmp;
100       }
101       if (numFaces) *numFaces = 4;
102       if (faceSize) *faceSize = 2;
103       break;
104     default:
105       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
106     }
107     break;
108   case 3:
109     switch (coneSize) {
110     case 3:
111       if (faces) {
112         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
113         facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
114         facesTmp[4] = cone[2]; facesTmp[5] = cone[0];
115         *faces = facesTmp;
116       }
117       if (numFaces) *numFaces = 3;
118       if (faceSize) *faceSize = 2;
119       break;
120     case 4:
121       /* Vertices of first face follow right hand rule and normal points away from last vertex */
122       if (faces) {
123         facesTmp[0] = cone[0]; facesTmp[1]  = cone[1]; facesTmp[2]  = cone[2];
124         facesTmp[3] = cone[0]; facesTmp[4]  = cone[3]; facesTmp[5]  = cone[1];
125         facesTmp[6] = cone[0]; facesTmp[7]  = cone[2]; facesTmp[8]  = cone[3];
126         facesTmp[9] = cone[2]; facesTmp[10] = cone[1]; facesTmp[11] = cone[3];
127         *faces = facesTmp;
128       }
129       if (numFaces) *numFaces = 4;
130       if (faceSize) *faceSize = 3;
131       break;
132     case 8:
133       /*  7--------6
134          /|       /|
135         / |      / |
136        4--------5  |
137        |  |     |  |
138        |  |     |  |
139        |  1--------2
140        | /      | /
141        |/       |/
142        0--------3
143        */
144       if (faces) {
145         facesTmp[0]  = cone[0]; facesTmp[1]  = cone[1]; facesTmp[2]  = cone[2]; facesTmp[3]  = cone[3]; /* Bottom */
146         facesTmp[4]  = cone[4]; facesTmp[5]  = cone[5]; facesTmp[6]  = cone[6]; facesTmp[7]  = cone[7]; /* Top */
147         facesTmp[8]  = cone[0]; facesTmp[9]  = cone[3]; facesTmp[10] = cone[5]; facesTmp[11] = cone[4]; /* Front */
148         facesTmp[12] = cone[2]; facesTmp[13] = cone[1]; facesTmp[14] = cone[7]; facesTmp[15] = cone[6]; /* Back */
149         facesTmp[16] = cone[3]; facesTmp[17] = cone[2]; facesTmp[18] = cone[6]; facesTmp[19] = cone[5]; /* Right */
150         facesTmp[20] = cone[0]; facesTmp[21] = cone[4]; facesTmp[22] = cone[7]; facesTmp[23] = cone[1]; /* Left */
151         *faces = facesTmp;
152       }
153       if (numFaces) *numFaces = 6;
154       if (faceSize) *faceSize = 4;
155       break;
156     default:
157       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
158     }
159     break;
160   default:
161     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim);
162   }
163   PetscFunctionReturn(0);
164 }
165 
166 /*
167   DMPlexGetRawFacesHybrid_Internal - Gets groups of vertices that correspond to faces for the given cone using hybrid ordering (prisms)
168 */
169 static PetscErrorCode DMPlexGetRawFacesHybrid_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *numFacesNotH, PetscInt *faceSize, const PetscInt *faces[])
170 {
171   PetscInt       *facesTmp;
172   PetscInt        maxConeSize, maxSupportSize;
173   PetscErrorCode  ierr;
174 
175   PetscFunctionBegin;
176   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
177   if (faces && coneSize) PetscValidIntPointer(cone,4);
178   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
179   if (faces) {ierr = DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), MPIU_INT, &facesTmp);CHKERRQ(ierr);}
180   switch (dim) {
181   case 1:
182     switch (coneSize) {
183     case 2:
184       if (faces) {
185         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
186         *faces = facesTmp;
187       }
188       if (numFaces)     *numFaces = 2;
189       if (numFacesNotH) *numFacesNotH = 2;
190       if (faceSize)     *faceSize = 1;
191       break;
192     default:
193       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
194     }
195     break;
196   case 2:
197     switch (coneSize) {
198     case 4:
199       if (faces) {
200         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
201         facesTmp[2] = cone[2]; facesTmp[3] = cone[3];
202         facesTmp[4] = cone[0]; facesTmp[5] = cone[2];
203         facesTmp[6] = cone[1]; facesTmp[7] = cone[3];
204         *faces = facesTmp;
205       }
206       if (numFaces)     *numFaces = 4;
207       if (numFacesNotH) *numFacesNotH = 2;
208       if (faceSize)     *faceSize = 2;
209       break;
210     default:
211       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
212     }
213     break;
214   case 3:
215     switch (coneSize) {
216     case 6: /* triangular prism */
217       if (faces) {
218         facesTmp[0]  = cone[0]; facesTmp[1]  = cone[1]; facesTmp[2]  = cone[2]; facesTmp[3]  = -1;      /* Bottom */
219         facesTmp[4]  = cone[3]; facesTmp[5]  = cone[4]; facesTmp[6]  = cone[5]; facesTmp[7]  = -1;      /* Top */
220         facesTmp[8]  = cone[0]; facesTmp[9]  = cone[1]; facesTmp[10] = cone[3]; facesTmp[11] = cone[4]; /* Back left */
221         facesTmp[12] = cone[1]; facesTmp[13] = cone[2]; facesTmp[14] = cone[4]; facesTmp[15] = cone[5]; /* Back right */
222         facesTmp[16] = cone[2]; facesTmp[17] = cone[0]; facesTmp[18] = cone[5]; facesTmp[19] = cone[3]; /* Front */
223         *faces = facesTmp;
224       }
225       if (numFaces)     *numFaces = 5;
226       if (numFacesNotH) *numFacesNotH = 2;
227       if (faceSize)     *faceSize = -4;
228       break;
229     default:
230       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
231     }
232     break;
233   default:
234     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim);
235   }
236   PetscFunctionReturn(0);
237 }
238 
239 static PetscErrorCode DMPlexRestoreRawFacesHybrid_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *numFacesNotH, PetscInt *faceSize, const PetscInt *faces[])
240 {
241   PetscErrorCode  ierr;
242 
243   PetscFunctionBegin;
244   if (faces) { ierr = DMRestoreWorkArray(dm, 0, MPIU_INT, (void *) faces);CHKERRQ(ierr); }
245   PetscFunctionReturn(0);
246 }
247 
248 static PetscErrorCode DMPlexGetFacesHybrid_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *numFacesNotH, PetscInt *faceSize, const PetscInt *faces[])
249 {
250   const PetscInt *cone = NULL;
251   PetscInt        coneSize;
252   PetscErrorCode  ierr;
253 
254   PetscFunctionBegin;
255   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
256   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
257   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
258   ierr = DMPlexGetRawFacesHybrid_Internal(dm, dim, coneSize, cone, numFaces, numFacesNotH, faceSize, faces);CHKERRQ(ierr);
259   PetscFunctionReturn(0);
260 }
261 
262 /* This interpolates faces for cells at some stratum */
263 static PetscErrorCode DMPlexInterpolateFaces_Internal(DM dm, PetscInt cellDepth, DM idm)
264 {
265   DMLabel        subpointMap;
266   PetscHashIJKL  faceTable;
267   PetscInt      *pStart, *pEnd;
268   PetscInt       cellDim, depth, faceDepth = cellDepth, numPoints = 0, faceSizeAll = 0, face, c, d;
269   PetscInt       coneSizeH = 0, faceSizeAllH = 0, faceSizeAllT = 0, numCellFacesH = 0, faceT = 0, faceH, pMax = -1, dim, outerloop;
270   PetscInt       cMax, fMax, eMax, vMax;
271   PetscErrorCode ierr;
272 
273   PetscFunctionBegin;
274   ierr = DMGetDimension(dm, &cellDim);CHKERRQ(ierr);
275   /* HACK: I need a better way to determine face dimension, or an alternative to GetFaces() */
276   ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr);
277   if (subpointMap) ++cellDim;
278   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
279   ++depth;
280   ++cellDepth;
281   cellDim -= depth - cellDepth;
282   ierr = PetscMalloc2(depth+1,&pStart,depth+1,&pEnd);CHKERRQ(ierr);
283   for (d = depth-1; d >= faceDepth; --d) {
284     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d+1], &pEnd[d+1]);CHKERRQ(ierr);
285   }
286   ierr = DMPlexGetDepthStratum(dm, -1, NULL, &pStart[faceDepth]);CHKERRQ(ierr);
287   pEnd[faceDepth] = pStart[faceDepth];
288   for (d = faceDepth-1; d >= 0; --d) {
289     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr);
290   }
291   cMax = fMax = eMax = vMax = PETSC_DETERMINE;
292   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
293   if (cellDim == dim) {
294     ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
295     pMax = cMax;
296   } else if (cellDim == dim -1) {
297     ierr = DMPlexGetHybridBounds(dm, &cMax, &fMax, NULL, NULL);CHKERRQ(ierr);
298     pMax = fMax;
299   }
300   pMax = pMax < 0 ? pEnd[cellDepth] : pMax;
301   if (pMax < pEnd[cellDepth]) {
302     const PetscInt *cellFaces, *cone;
303     PetscInt        numCellFacesT, faceSize, cf;
304 
305     /* First get normal cell face size (we now allow hybrid cells to meet normal cells on either hybrid or normal faces */
306     ierr = DMPlexGetFaces_Internal(dm, cellDim, pStart[cellDepth], NULL, &faceSizeAll, NULL);CHKERRQ(ierr);
307 
308     ierr = DMPlexGetConeSize(dm, pMax, &coneSizeH);CHKERRQ(ierr);
309     ierr = DMPlexGetCone(dm, pMax, &cone);CHKERRQ(ierr);
310     ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSizeH, cone, &numCellFacesH, &numCellFacesT, &faceSize, &cellFaces);CHKERRQ(ierr);
311     if (faceSize < 0) {
312       PetscInt *sizes, minv, maxv;
313 
314       /* count vertices of hybrid and non-hybrid faces */
315       ierr = PetscCalloc1(numCellFacesH, &sizes);CHKERRQ(ierr);
316       for (cf = 0; cf < numCellFacesT; ++cf) { /* These are the non-hybrid faces */
317         const PetscInt *cellFace = &cellFaces[-cf*faceSize];
318         PetscInt       f;
319 
320         for (f = 0; f < -faceSize; ++f) sizes[cf] += (cellFace[f] >= 0 ? 1 : 0);
321       }
322       ierr = PetscSortInt(numCellFacesT, sizes);CHKERRQ(ierr);
323       minv = sizes[0];
324       maxv = sizes[PetscMax(numCellFacesT-1, 0)];
325       if (minv != maxv) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Different number of vertices for non-hybrid face %D != %D", minv, maxv);
326       faceSizeAllT = minv;
327       ierr = PetscMemzero(sizes, numCellFacesH*sizeof(PetscInt));CHKERRQ(ierr);
328       for (cf = numCellFacesT; cf < numCellFacesH; ++cf) { /* These are the hybrid faces */
329         const PetscInt *cellFace = &cellFaces[-cf*faceSize];
330         PetscInt       f;
331 
332         for (f = 0; f < -faceSize; ++f) sizes[cf-numCellFacesT] += (cellFace[f] >= 0 ? 1 : 0);
333       }
334       ierr = PetscSortInt(numCellFacesH - numCellFacesT, sizes);CHKERRQ(ierr);
335       minv = sizes[0];
336       maxv = sizes[PetscMax(numCellFacesH - numCellFacesT-1, 0)];
337       ierr = PetscFree(sizes);CHKERRQ(ierr);
338       if (minv != maxv) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Different number of vertices for hybrid face %D != %D", minv, maxv);
339       faceSizeAllH = minv;
340     } else { /* the size of the faces in hybrid cells is the same */
341       faceSizeAll = faceSizeAllH = faceSizeAllT = faceSize;
342     }
343     ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSizeH, cone, &numCellFacesH, &numCellFacesT, &faceSize, &cellFaces);CHKERRQ(ierr);
344   } else if (pEnd[cellDepth] > pStart[cellDepth]) {
345     ierr = DMPlexGetFaces_Internal(dm, cellDim, pStart[cellDepth], NULL, &faceSizeAll, NULL);CHKERRQ(ierr);
346   }
347   if (faceSizeAll > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSizeAll);
348 
349   /* With hybrid grids, we first iterate on hybrid cells and start numbering the non-hybrid faces
350      Then, faces for non-hybrid cells are numbered.
351      This is to guarantee consistent orientations (all 0) of all the points in the cone of the hybrid cells */
352   ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr);
353   for (outerloop = 0, face = pStart[faceDepth]; outerloop < 2; outerloop++) {
354     PetscInt start, end;
355 
356     start = outerloop == 0 ? pMax : pStart[cellDepth];
357     end = outerloop == 0 ? pEnd[cellDepth] : pMax;
358     for (c = start; c < end; ++c) {
359       const PetscInt *cellFaces;
360       PetscInt        numCellFaces, faceSize, faceSizeInc, faceSizeCheck, cf;
361 
362       if (c < pMax) {
363         ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
364         if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll);
365         faceSizeCheck = faceSizeAll;
366       } else { /* Hybrid cell */
367         const PetscInt *cone;
368         PetscInt        numCellFacesN, coneSize;
369 
370         ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
371         if (coneSize != coneSizeH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid coneSize %D != %D", coneSize, coneSizeH);
372         ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr);
373         ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr);
374         if (numCellFaces != numCellFacesH) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected numCellFaces %D != %D for hybrid cell %D", numCellFaces, numCellFacesH, c);
375         faceSize = PetscMax(faceSize, -faceSize);
376         if (faceSize > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSize);
377         numCellFaces = numCellFacesN; /* process only non-hybrid faces */
378         faceSizeCheck = faceSizeAllT;
379       }
380       faceSizeInc = faceSize;
381       for (cf = 0; cf < numCellFaces; ++cf) {
382         const PetscInt   *cellFace = &cellFaces[cf*faceSizeInc];
383         PetscInt          faceSizeH = faceSize;
384         PetscHashIJKLKey  key;
385         PetscHashIter     iter;
386         PetscBool         missing;
387 
388         if (faceSizeInc == 2) {
389           key.i = PetscMin(cellFace[0], cellFace[1]);
390           key.j = PetscMax(cellFace[0], cellFace[1]);
391           key.k = PETSC_MAX_INT;
392           key.l = PETSC_MAX_INT;
393         } else {
394           key.i = cellFace[0];
395           key.j = cellFace[1];
396           key.k = cellFace[2];
397           key.l = faceSize > 3 ? (cellFace[3] < 0 ? faceSizeH = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT;
398           ierr  = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr);
399         }
400         /* this check is redundant for non-hybrid meshes */
401         if (faceSizeH != faceSizeCheck) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected number of vertices for face %D of point %D -> %D != %D", cf, c, faceSizeH, faceSizeCheck);
402         ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr);
403         if (missing) {
404           ierr = PetscHashIJKLIterSet(faceTable, iter, face++);CHKERRQ(ierr);
405           if (c >= pMax) ++faceT;
406         }
407       }
408       if (c < pMax) {
409         ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
410       } else {
411         ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSizeH, NULL, NULL, NULL, NULL, &cellFaces);CHKERRQ(ierr);
412       }
413     }
414   }
415   pEnd[faceDepth] = face;
416 
417   /* Second pass for hybrid meshes: number hybrid faces */
418   for (c = pMax; c < pEnd[cellDepth]; ++c) {
419     const PetscInt *cellFaces, *cone;
420     PetscInt        numCellFaces, numCellFacesN, faceSize, cf, coneSize;
421 
422     ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
423     ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr);
424     ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr);
425     if (numCellFaces != numCellFacesH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid numCellFaces %D != %D", numCellFaces, numCellFacesH);
426     faceSize = PetscMax(faceSize, -faceSize);
427     for (cf = numCellFacesN; cf < numCellFaces; ++cf) { /* These are the hybrid faces */
428       const PetscInt   *cellFace = &cellFaces[cf*faceSize];
429       PetscHashIJKLKey  key;
430       PetscHashIter     iter;
431       PetscBool         missing;
432       PetscInt          faceSizeH = faceSize;
433 
434       if (faceSize == 2) {
435         key.i = PetscMin(cellFace[0], cellFace[1]);
436         key.j = PetscMax(cellFace[0], cellFace[1]);
437         key.k = PETSC_MAX_INT;
438         key.l = PETSC_MAX_INT;
439       } else {
440         key.i = cellFace[0];
441         key.j = cellFace[1];
442         key.k = cellFace[2];
443         key.l = faceSize > 3 ? (cellFace[3] < 0 ? faceSizeH = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT;
444         ierr  = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr);
445       }
446       if (faceSizeH != faceSizeAllH) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected number of vertices for hybrid face %D of point %D -> %D != %D", cf, c, faceSizeH, faceSizeAllH);
447       ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr);
448       if (missing) {ierr = PetscHashIJKLIterSet(faceTable, iter, face++);CHKERRQ(ierr);}
449     }
450     ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr);
451   }
452   faceH = face - pEnd[faceDepth];
453   if (faceH) {
454     if (fMax == PETSC_DETERMINE) fMax = pEnd[faceDepth];
455     else if (eMax == PETSC_DETERMINE) eMax = pEnd[faceDepth];
456     else SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of unassigned hybrid facets %D for cellDim %D and dimension %D", faceH, cellDim, dim);
457   }
458   pEnd[faceDepth] = face;
459   ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr);
460   /* Count new points */
461   for (d = 0; d <= depth; ++d) {
462     numPoints += pEnd[d]-pStart[d];
463   }
464   ierr = DMPlexSetChart(idm, 0, numPoints);CHKERRQ(ierr);
465   /* Set cone sizes */
466   for (d = 0; d <= depth; ++d) {
467     PetscInt coneSize, p;
468 
469     if (d == faceDepth) {
470       /* Now we have two cases: */
471       if (faceSizeAll == faceSizeAllT) {
472         /* I see no way to do this if we admit faces of different shapes */
473         for (p = pStart[d]; p < pEnd[d]-faceH; ++p) {
474           ierr = DMPlexSetConeSize(idm, p, faceSizeAll);CHKERRQ(ierr);
475         }
476         for (p = pEnd[d]-faceH; p < pEnd[d]; ++p) {
477           ierr = DMPlexSetConeSize(idm, p, faceSizeAllH);CHKERRQ(ierr);
478         }
479       } else if (faceSizeAll == faceSizeAllH) {
480         for (p = pStart[d]; p < pStart[d]+faceT; ++p) {
481           ierr = DMPlexSetConeSize(idm, p, faceSizeAllT);CHKERRQ(ierr);
482         }
483         for (p = pStart[d]+faceT; p < pEnd[d]-faceH; ++p) {
484           ierr = DMPlexSetConeSize(idm, p, faceSizeAll);CHKERRQ(ierr);
485         }
486         for (p = pEnd[d]-faceH; p < pEnd[d]; ++p) {
487           ierr = DMPlexSetConeSize(idm, p, faceSizeAllH);CHKERRQ(ierr);
488         }
489       } else SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent faces sizes N: %D T: %D H: %D", faceSizeAll, faceSizeAllT, faceSizeAllH);
490     } else if (d == cellDepth) {
491       for (p = pStart[d]; p < pEnd[d]; ++p) {
492         /* Number of cell faces may be different from number of cell vertices*/
493         if (p < pMax) {
494           ierr = DMPlexGetFaces_Internal(dm, cellDim, p, &coneSize, NULL, NULL);CHKERRQ(ierr);
495         } else {
496           ierr = DMPlexGetFacesHybrid_Internal(dm, cellDim, p, &coneSize, NULL, NULL, NULL);CHKERRQ(ierr);
497         }
498         ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr);
499       }
500     } else {
501       for (p = pStart[d]; p < pEnd[d]; ++p) {
502         ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
503         ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr);
504       }
505     }
506   }
507   ierr = DMSetUp(idm);CHKERRQ(ierr);
508   /* Get face cones from subsets of cell vertices */
509   if (faceSizeAll > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSizeAll);
510   ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr);
511   for (d = depth; d > cellDepth; --d) {
512     const PetscInt *cone;
513     PetscInt        p;
514 
515     for (p = pStart[d]; p < pEnd[d]; ++p) {
516       ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
517       ierr = DMPlexSetCone(idm, p, cone);CHKERRQ(ierr);
518       ierr = DMPlexGetConeOrientation(dm, p, &cone);CHKERRQ(ierr);
519       ierr = DMPlexSetConeOrientation(idm, p, cone);CHKERRQ(ierr);
520     }
521   }
522   for (outerloop = 0, face = pStart[faceDepth]; outerloop < 2; outerloop++) {
523     PetscInt start, end;
524 
525     start = outerloop == 0 ? pMax : pStart[cellDepth];
526     end = outerloop == 0 ? pEnd[cellDepth] : pMax;
527     for (c = start; c < end; ++c) {
528       const PetscInt *cellFaces;
529       PetscInt        numCellFaces, faceSize, faceSizeInc, cf;
530 
531       if (c < pMax) {
532         ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
533         if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll);
534       } else {
535         const PetscInt *cone;
536         PetscInt        numCellFacesN, coneSize;
537 
538         ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
539         if (coneSize != coneSizeH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid coneSize %D != %D", coneSize, coneSizeH);
540         ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr);
541         ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr);
542         if (numCellFaces != numCellFacesH) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected numCellFaces %D != %D for hybrid cell %D", numCellFaces, numCellFacesH, c);
543         faceSize = PetscMax(faceSize, -faceSize);
544         if (faceSize > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSize);
545         numCellFaces = numCellFacesN; /* process only non-hybrid faces */
546       }
547       faceSizeInc = faceSize;
548       for (cf = 0; cf < numCellFaces; ++cf) {
549         const PetscInt  *cellFace = &cellFaces[cf*faceSizeInc];
550         PetscHashIJKLKey key;
551         PetscHashIter    iter;
552         PetscBool        missing;
553 
554         if (faceSizeInc == 2) {
555           key.i = PetscMin(cellFace[0], cellFace[1]);
556           key.j = PetscMax(cellFace[0], cellFace[1]);
557           key.k = PETSC_MAX_INT;
558           key.l = PETSC_MAX_INT;
559         } else {
560           key.i = cellFace[0];
561           key.j = cellFace[1];
562           key.k = cellFace[2];
563           key.l = faceSizeInc > 3 ? (cellFace[3] < 0 ? faceSize = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT;
564           ierr  = PetscSortInt(faceSizeInc, (PetscInt *) &key);CHKERRQ(ierr);
565         }
566         ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr);
567         if (missing) {
568           ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr);
569           ierr = PetscHashIJKLIterSet(faceTable, iter, face);CHKERRQ(ierr);
570           ierr = DMPlexInsertCone(idm, c, cf, face++);CHKERRQ(ierr);
571         } else {
572           const PetscInt *cone;
573           PetscInt        coneSize, ornt, i, j, f;
574 
575           ierr = PetscHashIJKLIterGet(faceTable, iter, &f);CHKERRQ(ierr);
576           ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr);
577           /* Orient face: Do not allow reverse orientation at the first vertex */
578           ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr);
579           ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr);
580           if (coneSize != faceSize) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of face vertices %D for face %D should be %D", coneSize, f, faceSize);
581           /* - First find the initial vertex */
582           for (i = 0; i < faceSize; ++i) if (cellFace[0] == cone[i]) break;
583           /* - Try forward comparison */
584           for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+j)%faceSize]) break;
585           if (j == faceSize) {
586             if ((faceSize == 2) && (i == 1)) ornt = -2;
587             else                             ornt = i;
588           } else {
589             /* - Try backward comparison */
590             for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+faceSize-j)%faceSize]) break;
591             if (j == faceSize) {
592               if (i == 0) ornt = -faceSize;
593               else        ornt = -i;
594             } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine orientation of face %D in cell %D", f, c);
595           }
596           ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr);
597         }
598       }
599       if (c < pMax) {
600         ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
601       } else {
602         ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSizeH, NULL, NULL, NULL, NULL, &cellFaces);CHKERRQ(ierr);
603       }
604     }
605   }
606   /* Second pass for hybrid meshes: orient hybrid faces */
607   for (c = pMax; c < pEnd[cellDepth]; ++c) {
608     const PetscInt *cellFaces, *cone;
609     PetscInt        numCellFaces, numCellFacesN, faceSize, cf, coneSize;
610 
611     ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
612     ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr);
613     ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr);
614     if (numCellFaces != numCellFacesH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid numCellFaces %D != %D", numCellFaces, numCellFacesH);
615     faceSize = PetscMax(faceSize, -faceSize);
616     for (cf = numCellFacesN; cf < numCellFaces; ++cf) { /* These are the hybrid faces */
617       const PetscInt   *cellFace = &cellFaces[cf*faceSize];
618       PetscHashIJKLKey key;
619       PetscHashIter    iter;
620       PetscBool        missing;
621       PetscInt         faceSizeH = faceSize;
622 
623       if (faceSize == 2) {
624         key.i = PetscMin(cellFace[0], cellFace[1]);
625         key.j = PetscMax(cellFace[0], cellFace[1]);
626         key.k = PETSC_MAX_INT;
627         key.l = PETSC_MAX_INT;
628       } else {
629         key.i = cellFace[0];
630         key.j = cellFace[1];
631         key.k = cellFace[2];
632         key.l = faceSize > 3 ? (cellFace[3] < 0 ? faceSizeH = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT;
633         ierr  = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr);
634       }
635       if (faceSizeH != faceSizeAllH) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected number of vertices for hybrid face %D of point %D -> %D != %D", cf, c, faceSizeH, faceSizeAllH);
636       ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr);
637       if (missing) {
638         ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr);
639         ierr = PetscHashIJKLIterSet(faceTable, iter, face);CHKERRQ(ierr);
640         ierr = DMPlexInsertCone(idm, c, cf, face++);CHKERRQ(ierr);
641       } else {
642         PetscInt        fv[4] = {0, 1, 2, 3};
643         const PetscInt *cone;
644         PetscInt        coneSize, ornt, i, j, f;
645 
646         ierr = PetscHashIJKLIterGet(faceTable, iter, &f);CHKERRQ(ierr);
647         ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr);
648         /* Orient face: Do not allow reverse orientation at the first vertex */
649         ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr);
650         ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr);
651         if (coneSize != faceSizeH) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of face vertices %D for face %D should be %D", coneSize, f, faceSizeH);
652         /* Hybrid faces are stored as tensor products of edges, so to compare them to normal faces, we have to flip */
653         if (faceSize == 4 && c >= pMax && faceSizeAll != faceSizeAllT && f < pEnd[faceDepth] - faceH) {fv[2] = 3; fv[3] = 2;}
654         /* - First find the initial vertex */
655         for (i = 0; i < faceSizeH; ++i) if (cellFace[fv[0]] == cone[i]) break;
656         /* - Try forward comparison */
657         for (j = 0; j < faceSizeH; ++j) if (cellFace[fv[j]] != cone[(i+j)%faceSizeH]) break;
658         if (j == faceSizeH) {
659           if ((faceSizeH == 2) && (i == 1)) ornt = -2;
660           else                              ornt = i;
661         } else {
662           /* - Try backward comparison */
663           for (j = 0; j < faceSizeH; ++j) if (cellFace[fv[j]] != cone[(i+faceSizeH-j)%faceSizeH]) break;
664           if (j == faceSizeH) {
665             if (i == 0) ornt = -faceSizeH;
666             else        ornt = -i;
667           } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine orientation of face %D in cell %D", f, c);
668         }
669         ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr);
670       }
671     }
672     ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr);
673   }
674   if (face != pEnd[faceDepth]) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of faces %D should be %D", face-pStart[faceDepth], pEnd[faceDepth]-pStart[faceDepth]);
675   ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr);
676   ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr);
677   ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr);
678   ierr = DMPlexSetHybridBounds(idm, cMax, fMax, eMax, vMax);CHKERRQ(ierr);
679   ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr);
680   ierr = DMPlexStratify(idm);CHKERRQ(ierr);
681   PetscFunctionReturn(0);
682 }
683 
684 PetscErrorCode DMPlexOrientCell(DM dm, PetscInt p, PetscInt masterConeSize, const PetscInt masterCone[])
685 {
686   PetscInt coneSize;
687   PetscInt start1=0;
688   PetscBool reverse1=PETSC_FALSE;
689   const PetscInt *cone=NULL;
690   PetscErrorCode ierr;
691 
692   PetscFunctionBegin;
693   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
694   if (!coneSize) PetscFunctionReturn(0); /* do nothing for points with no cone */
695   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
696   ierr = DMPlexFixFaceOrientations_Orient_Private(coneSize, masterConeSize, masterCone, cone, &start1, &reverse1);CHKERRQ(ierr);
697 #if defined(PETSC_USE_DEBUG)
698   if (PetscUnlikely(cone[start1] != masterCone[0])) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "The algorithm above is wrong as cone[%d] = %d != %d = masterCone[0]", start1, cone[start1], masterCone[0]);
699 #endif
700   ierr = DMPlexOrientCell_Internal(dm, p, start1, reverse1);CHKERRQ(ierr);
701 #if defined(PETSC_USE_DEBUG)
702   {
703     PetscInt c;
704     ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
705     for (c = 0; c < 2; c++) {
706       if (PetscUnlikely(cone[c] != masterCone[c])) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_PLIB, "The algorithm above is wrong as cone[%d] = %d != %d = masterCone[%d]", c, cone[c], masterCone[c], c);
707     }
708   }
709 #endif
710   PetscFunctionReturn(0);
711 }
712 
713 PetscErrorCode DMPlexOrientCell_Internal(DM dm, PetscInt p, PetscInt start1, PetscBool reverse1)
714 {
715   PetscInt i, j, k, maxConeSize, coneSize, coneConeSize, supportSize, supportConeSize;
716   PetscInt start0, start;
717   PetscBool reverse0, reverse;
718   PetscInt newornt;
719   const PetscInt *cone=NULL, *support=NULL, *supportCone=NULL, *ornts=NULL;
720   PetscInt *newcone=NULL, *newornts=NULL;
721   PetscErrorCode ierr;
722 
723   PetscFunctionBegin;
724   if (!start1 && !reverse1) PetscFunctionReturn(0);
725   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
726   if (!coneSize) PetscFunctionReturn(0); /* do nothing for points with no cone */
727   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
728   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
729   /* permute p's cone and orientations */
730   ierr = DMPlexGetConeOrientation(dm, p, &ornts);CHKERRQ(ierr);
731   ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &newcone);CHKERRQ(ierr);
732   ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &newornts);CHKERRQ(ierr);
733   ierr = DMPlexFixFaceOrientations_Permute_Private(coneSize, cone, start1, reverse1, newcone);CHKERRQ(ierr);
734   ierr = DMPlexFixFaceOrientations_Permute_Private(coneSize, ornts, start1, reverse1, newornts);CHKERRQ(ierr);
735   /* if direction of p (face) is flipped, flip also p's cone points (edges) */
736   if (reverse1) {
737     for (i=0; i<coneSize; i++) {
738       ierr = DMPlexGetConeSize(dm, cone[i], &coneConeSize);CHKERRQ(ierr);
739       ierr = DMPlexFixFaceOrientations_Translate_Private(newornts[i], &start0, &reverse0);CHKERRQ(ierr);
740       ierr = DMPlexFixFaceOrientations_Combine_Private(coneConeSize, start0, reverse0, 1, PETSC_FALSE, &start, &reverse);CHKERRQ(ierr);
741       ierr = DMPlexFixFaceOrientations_TranslateBack_Private(coneConeSize, start, reverse, &newornts[i]);CHKERRQ(ierr);
742     }
743   }
744   ierr = DMPlexSetConeOrientation(dm, p, newornts);CHKERRQ(ierr);
745   /* fix oriention of p within cones of p's support points */
746   ierr = DMPlexGetSupport(dm, p, &support);CHKERRQ(ierr);
747   ierr = DMPlexGetSupportSize(dm, p, &supportSize);CHKERRQ(ierr);
748   for (j=0; j<supportSize; j++) {
749     ierr = DMPlexGetCone(dm, support[j], &supportCone);CHKERRQ(ierr);
750     ierr = DMPlexGetConeSize(dm, support[j], &supportConeSize);CHKERRQ(ierr);
751     ierr = DMPlexGetConeOrientation(dm, support[j], &ornts);CHKERRQ(ierr);
752     for (k=0; k<supportConeSize; k++) {
753       if (supportCone[k] != p) continue;
754       ierr = DMPlexFixFaceOrientations_Translate_Private(ornts[k], &start0, &reverse0);CHKERRQ(ierr);
755       ierr = DMPlexFixFaceOrientations_Combine_Private(coneSize, start0, reverse0, start1, reverse1, &start, &reverse);CHKERRQ(ierr);
756       ierr = DMPlexFixFaceOrientations_TranslateBack_Private(coneSize, start, reverse, &newornt);CHKERRQ(ierr);
757       ierr = DMPlexInsertConeOrientation(dm, support[j], k, newornt);CHKERRQ(ierr);
758     }
759   }
760   /* rewrite cone */
761   ierr = DMPlexSetCone(dm, p, newcone);CHKERRQ(ierr);
762   ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &newcone);CHKERRQ(ierr);
763   ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &newornts);CHKERRQ(ierr);
764   PetscFunctionReturn(0);
765 }
766 
767 static PetscErrorCode SortRmineRremoteByRemote_Private(PetscSF sf, PetscInt *rmine1[], PetscInt *rremote1[])
768 {
769   PetscInt            nleaves;
770   PetscInt            nranks;
771   const PetscMPIInt  *ranks=NULL;
772   const PetscInt     *roffset=NULL, *rmine=NULL, *rremote=NULL;
773   PetscInt            n, o, r;
774   PetscErrorCode      ierr;
775 
776   PetscFunctionBegin;
777   ierr = PetscSFGetRanks(sf, &nranks, &ranks, &roffset, &rmine, &rremote);CHKERRQ(ierr);
778   nleaves = roffset[nranks];
779   ierr = PetscMalloc2(nleaves, rmine1, nleaves, rremote1);CHKERRQ(ierr);
780   for (r=0; r<nranks; r++) {
781     /* simultaneously sort rank-wise portions of rmine & rremote by values in rremote
782        - to unify order with the other side */
783     o = roffset[r];
784     n = roffset[r+1] - o;
785     ierr = PetscMemcpy(&(*rmine1)[o], &rmine[o], n*sizeof(PetscInt));CHKERRQ(ierr);
786     ierr = PetscMemcpy(&(*rremote1)[o], &rremote[o], n*sizeof(PetscInt));CHKERRQ(ierr);
787     ierr = PetscSortIntWithArray(n, &(*rremote1)[o], &(*rmine1)[o]);CHKERRQ(ierr);
788   }
789   PetscFunctionReturn(0);
790 }
791 
792 PetscErrorCode DMPlexOrientInterface(DM dm)
793 {
794   PetscSF           sf=NULL;
795   PetscInt          (*roots)[2], (*leaves)[2];
796   PetscMPIInt       (*rootsRanks)[2], (*leavesRanks)[2];
797   const PetscInt    *locals=NULL;
798   const PetscSFNode *remotes=NULL;
799   PetscInt           nroots, nleaves, p, c;
800   PetscInt           nranks, n, o, r;
801   const PetscMPIInt *ranks=NULL;
802   const PetscInt    *roffset=NULL;
803   PetscInt          *rmine1=NULL, *rremote1=NULL; /* rmine and rremote copies simultaneously sorted by rank and rremote */
804   const PetscInt    *cone=NULL;
805   PetscInt           coneSize, ind0;
806   MPI_Comm           comm;
807   PetscMPIInt        rank;
808   PetscInt           debug = 0;
809   PetscErrorCode     ierr;
810 
811   PetscFunctionBegin;
812   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
813   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &locals, &remotes);CHKERRQ(ierr);
814   if (nroots < 0) PetscFunctionReturn(0);
815   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
816   ierr = PetscSFGetRanks(sf, &nranks, &ranks, &roffset, NULL, NULL);CHKERRQ(ierr);
817 #if defined(PETSC_USE_DEBUG)
818   ierr = DMViewFromOptions(dm, NULL, "-before_fix_dm_view");CHKERRQ(ierr);
819   ierr = DMPlexCheckPointSF(dm);CHKERRQ(ierr);
820 #endif
821   ierr = SortRmineRremoteByRemote_Private(sf, &rmine1, &rremote1);CHKERRQ(ierr);
822   ierr = PetscMalloc4(nroots, &roots, nroots, &leaves, nroots, &rootsRanks, nroots, &leavesRanks);CHKERRQ(ierr);
823   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
824   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
825   if (debug && rank == 0) {ierr = PetscSynchronizedPrintf(comm, "Roots\n");CHKERRQ(ierr);}
826   for (p = 0; p < nroots; ++p) {
827     ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
828     ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
829     /* Translate all points to root numbering */
830     for (c = 0; c < 2; c++) {
831       if (coneSize > 1) {
832         ierr = PetscFindInt(cone[c], nleaves, locals, &ind0);CHKERRQ(ierr);
833         if (ind0 < 0) {
834           roots[p][c] = cone[c];
835           rootsRanks[p][c] = rank;
836         } else {
837           roots[p][c] = remotes[ind0].index;
838           rootsRanks[p][c] = remotes[ind0].rank;
839         }
840       } else {
841         roots[p][c] = -1;
842         rootsRanks[p][c] = -1;
843       }
844     }
845   }
846   if (debug) {
847     for (p = 0; p < nroots; ++p) {
848       ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
849       ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
850       if (coneSize > 1) {
851         ierr = PetscSynchronizedPrintf(comm, "[%d]  %D: cone=[%D %D] roots=[%D %D] rootsRanks=[%D %D]\n", rank, p, cone[0], cone[1], roots[p][0], roots[p][1], rootsRanks[p][0], rootsRanks[p][1]);CHKERRQ(ierr);
852       }
853     }
854     ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);
855   }
856   for (p = 0; p < nroots; ++p) {
857     for (c = 0; c < 2; c++) {
858       leaves[p][c] = -2;
859       leavesRanks[p][c] = -2;
860     }
861   }
862   ierr = PetscSFBcastBegin(sf, MPIU_2INT, roots, leaves);CHKERRQ(ierr);
863   ierr = PetscSFBcastBegin(sf, MPI_2INT, rootsRanks, leavesRanks);CHKERRQ(ierr);
864   ierr = PetscSFBcastEnd(sf, MPIU_2INT, roots, leaves);CHKERRQ(ierr);
865   ierr = PetscSFBcastEnd(sf, MPI_2INT, rootsRanks, leavesRanks);CHKERRQ(ierr);
866   if (debug) {ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);}
867   if (debug && rank == 0) {ierr = PetscSynchronizedPrintf(comm, "Referred leaves\n");CHKERRQ(ierr);}
868   for (p = 0; p < nroots; ++p) {
869     if (leaves[p][0] < 0) continue;
870     ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
871     ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
872     if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d]  %D: cone=[%D %D] leaves=[%D %D] roots=[%D %D] leavesRanks=[%D %D] rootsRanks=[%D %D]\n", rank, p, cone[0], cone[1], leaves[p][0], leaves[p][1], roots[p][0], roots[p][1], leavesRanks[p][0], leavesRanks[p][1], rootsRanks[p][0], rootsRanks[p][1]);CHKERRQ(ierr);}
873     if ((leaves[p][0] != roots[p][0]) || (leaves[p][1] != roots[p][1]) || (leavesRanks[p][0] != rootsRanks[p][0]) || (leavesRanks[p][0] != rootsRanks[p][0])) {
874       PetscInt masterCone[2];
875       /* Translate these two cone points back to leave numbering */
876       for (c = 0; c < 2; c++) {
877         if (leavesRanks[p][c] == rank) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "this should never happen - remote rank of point %D is the same rank",leavesRanks[p][c]);
878         /* Find index of rank leavesRanks[p][c] among remote ranks */
879         /* No need for PetscMPIIntCast because these integers were originally cast from PetscMPIInt. */
880         ierr = PetscFindMPIInt((PetscMPIInt)leavesRanks[p][c], nranks, ranks, &r);CHKERRQ(ierr);
881         if (PetscUnlikely(r < 0)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "this should never happen - rank %D not found among remote ranks",leavesRanks[p][c]);
882         /* Find point leaves[p][c] among remote points aimed at rank leavesRanks[p][c] */
883         o = roffset[r];
884         n = roffset[r+1] - o;
885         ierr = PetscFindInt(leaves[p][c], n, &rremote1[o], &ind0);CHKERRQ(ierr);
886         if (PetscUnlikely(ind0 < 0)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "No cone point of %D is connected to (%D, %D) - it seems there is missing connection in point SF!",p,ranks[r],leaves[p][c]);
887         /* Get the corresponding local point */
888         masterCone[c] = rmine1[o+ind0];CHKERRQ(ierr);
889       }
890       if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d]  %D: masterCone=[%D %D]\n", rank, p, masterCone[0], masterCone[1]);CHKERRQ(ierr);}
891       /* Vaclav's note: Here we only compare first 2 points of the cone. Full cone size would lead to stronger self-checking. */
892       ierr = DMPlexOrientCell(dm, p, 2, masterCone);CHKERRQ(ierr);
893     }
894   }
895 #if defined(PETSC_USE_DEBUG)
896   ierr = DMViewFromOptions(dm, NULL, "-after_fix_dm_view");CHKERRQ(ierr);
897   for (r = 0; r < nleaves; ++r) {
898     p = locals[r];
899     ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
900     if (!coneSize) continue;
901     ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
902     for (c = 0; c < 2; c++) {
903       ierr = PetscFindInt(cone[c], nleaves, locals, &ind0);CHKERRQ(ierr);
904       if (ind0 < 0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point SF contains %D but is missing its cone point cone[%D] = %D!", p, c, cone[c]);
905       if (leaves[p][c] != remotes[ind0].index || leavesRanks[p][c] != remotes[ind0].rank) {
906         if (leavesRanks[p][c] == rank) {
907           PetscInt ind1;
908           ierr = PetscFindInt(leaves[p][c], nleaves, locals, &ind1);CHKERRQ(ierr);
909           if (ind1 < 0) {
910             SETERRQ8(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D = locals[%d]: cone[%D]=%D --> (%D, %D) differs from the enforced (%D, %D). The latter was not even found among the local SF points - it is probably broken!", p, r, c, cone[c], remotes[ind0].rank, remotes[ind0].index, leavesRanks[p][c], leaves[p][c]);
911           } else {
912             SETERRQ9(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D = locals[%d]: cone[%D]=%D --> (%D, %D) differs from the enforced %D --> (%D, %D). Is the algorithm above or the point SF broken?", p, r, c, cone[c], remotes[ind0].rank, remotes[ind0].index, leaves[p][c], remotes[ind1].rank, remotes[ind1].index);
913           }
914         } else {
915           SETERRQ8(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D = locals[%d]: cone[%D]=%D --> (%D, %D) differs from the enforced (%D, %D). Is the algorithm above or the point SF broken?", p, r, c, cone[c], remotes[ind0].rank, remotes[ind0].index, leavesRanks[p][c], leaves[p][c]);
916         }
917       }
918     }
919   }
920 #endif
921   if (debug) {ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);}
922   ierr = PetscFree4(roots, leaves, rootsRanks, leavesRanks);CHKERRQ(ierr);
923   ierr = PetscFree2(rmine1, rremote1);CHKERRQ(ierr);
924   PetscFunctionReturn(0);
925 }
926 
927 static PetscErrorCode IntArrayViewFromOptions(MPI_Comm comm, const char opt[], const char name[], const char idxname[], const char valname[], PetscInt n, const PetscInt a[])
928 {
929   PetscInt       idx;
930   PetscMPIInt    rank;
931   PetscBool      flg;
932   PetscErrorCode ierr;
933 
934   PetscFunctionBegin;
935   ierr = PetscOptionsHasName(NULL, NULL, opt, &flg);CHKERRQ(ierr);
936   if (!flg) PetscFunctionReturn(0);
937   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
938   ierr = PetscSynchronizedPrintf(comm, "[%d]%s:\n", rank, name);CHKERRQ(ierr);
939   for (idx = 0; idx < n; ++idx) {ierr = PetscSynchronizedPrintf(comm, "[%d]%s %D %s %D\n", rank, idxname, idx, valname, a[idx]);CHKERRQ(ierr);}
940   ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);
941   PetscFunctionReturn(0);
942 }
943 
944 static PetscErrorCode SFNodeArrayViewFromOptions(MPI_Comm comm, const char opt[], const char name[], const char idxname[], PetscInt n, const PetscSFNode a[])
945 {
946   PetscInt       idx;
947   PetscMPIInt    rank;
948   PetscBool      flg;
949   PetscErrorCode ierr;
950 
951   PetscFunctionBegin;
952   ierr = PetscOptionsHasName(NULL, NULL, opt, &flg);CHKERRQ(ierr);
953   if (!flg) PetscFunctionReturn(0);
954   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
955   ierr = PetscSynchronizedPrintf(comm, "[%d]%s:\n", rank, name);CHKERRQ(ierr);
956   if (idxname) {
957     for (idx = 0; idx < n; ++idx) {ierr = PetscSynchronizedPrintf(comm, "[%d]%s %D rank %D index %D\n", rank, idxname, idx, a[idx].rank, a[idx].index);CHKERRQ(ierr);}
958   } else {
959     for (idx = 0; idx < n; ++idx) {ierr = PetscSynchronizedPrintf(comm, "[%d]rank %D index %D\n", rank, a[idx].rank, a[idx].index);CHKERRQ(ierr);}
960   }
961   ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);
962   PetscFunctionReturn(0);
963 }
964 
965 static PetscErrorCode DMPlexMapToLocalPoint(PetscHMapIJ roothash, const PetscInt localPoints[], PetscMPIInt rank, PetscSFNode remotePoint, PetscInt *localPoint)
966 {
967   PetscErrorCode ierr;
968 
969   PetscFunctionBegin;
970   if (remotePoint.rank == rank) {
971     *localPoint = remotePoint.index;
972   } else {
973     PetscHashIJKey key;
974     PetscInt       root;
975 
976     key.i = remotePoint.index;
977     key.j = remotePoint.rank;
978     ierr = PetscHMapIJGet(roothash, key, &root);CHKERRQ(ierr);
979     if (root >= 0) {
980       *localPoint = localPoints[root];
981     } else PetscFunctionReturn(1);
982   }
983   PetscFunctionReturn(0);
984 }
985 
986 /*@
987   DMPlexInterpolatePointSF - Insert interpolated points in the overlap into the PointSF in parallel, following local interpolation
988 
989   Collective on DM
990 
991   Input Parameters:
992 + dm      - The interpolated DM
993 - pointSF - The initial SF without interpolated points
994 
995   Output Parameter:
996 . pointSF - The SF including interpolated points
997 
998   Level: intermediate
999 
1000    Note: All debugging for this process can be turned on with the options: -dm_interp_pre_view -petscsf_interp_pre_view -petscsection_interp_candidate_view -petscsection_interp_candidate_remote_view -petscsection_interp_claim_view -petscsf_interp_pre_view -dmplex_interp_debug
1001 
1002 .keywords: mesh
1003 .seealso: DMPlexInterpolate(), DMPlexUninterpolate()
1004 @*/
1005 PetscErrorCode DMPlexInterpolatePointSF(DM dm, PetscSF pointSF)
1006 {
1007   /*
1008        Okay, the algorithm is:
1009          - Take each point in the overlap (root)
1010          - Look at the neighboring points in the overlap (candidates)
1011          - Send these candidate points to neighbors
1012          - Neighbor checks for edge between root and candidate
1013          - If edge is found, it replaces candidate point with edge point
1014          - Send back the overwritten candidates (claims)
1015          - Original guy checks for edges, different from original candidate, and gets its own edge
1016          - This pair is put into SF
1017 
1018        We need a new algorithm that tolerates groups larger than 2.
1019          - Take each point in the overlap (root)
1020          - Find all collections of points in the overlap which make faces (do early join)
1021          - Send collections as candidates (add size as first number)
1022            - Make sure to send collection to all owners of all overlap points in collection
1023          - Neighbor check for face in collections
1024          - If face is found, it replaces candidate point with face point
1025          - Send back the overwritten candidates (claims)
1026          - Original guy checks for faces, different from original candidate, and gets its own face
1027          - This pair is put into SF
1028   */
1029   PetscHMapI         leafhash;
1030   PetscHMapIJ        roothash;
1031   const PetscInt    *localPoints, *rootdegree;
1032   const PetscSFNode *remotePoints;
1033   PetscSFNode       *candidates, *candidatesRemote, *claims;
1034   PetscSection       candidateSection, candidateSectionRemote, claimSection;
1035   PetscInt           numLeaves, l, numRoots, r, candidatesSize, candidatesRemoteSize;
1036   PetscMPIInt        size, rank;
1037   PetscHashIJKey     key;
1038   PetscBool          debug = PETSC_FALSE;
1039   PetscErrorCode     ierr;
1040 
1041   PetscFunctionBegin;
1042   ierr = PetscOptionsHasName(NULL, ((PetscObject) dm)->prefix, "-dmplex_interp_debug", &debug);CHKERRQ(ierr);
1043   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRQ(ierr);
1044   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
1045   ierr = PetscSFGetGraph(pointSF, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr);
1046   if (size < 2 || numRoots < 0) PetscFunctionReturn(0);
1047   ierr = PetscObjectViewFromOptions((PetscObject) dm, NULL, "-dm_interp_pre_view");CHKERRQ(ierr);
1048   ierr = PetscObjectViewFromOptions((PetscObject) pointSF, NULL, "-petscsf_interp_pre_view");CHKERRQ(ierr);
1049   ierr = PetscLogEventBegin(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr);
1050   /* Build hashes of points in the SF for efficient lookup */
1051   ierr = PetscHMapICreate(&leafhash);CHKERRQ(ierr);
1052   ierr = PetscHMapIJCreate(&roothash);CHKERRQ(ierr);
1053   for (l = 0; l < numLeaves; ++l) {
1054     key.i = remotePoints[l].index;
1055     key.j = remotePoints[l].rank;
1056     ierr = PetscHMapISet(leafhash, localPoints[l], l);CHKERRQ(ierr);
1057     ierr = PetscHMapIJSet(roothash, key, l);CHKERRQ(ierr);
1058   }
1059   /* Compute root degree to identify shared points */
1060   ierr = PetscSFComputeDegreeBegin(pointSF, &rootdegree);CHKERRQ(ierr);
1061   ierr = PetscSFComputeDegreeEnd(pointSF, &rootdegree);CHKERRQ(ierr);
1062   ierr = IntArrayViewFromOptions(PetscObjectComm((PetscObject) dm), "-interp_root_degree_view", "Root degree", "point", "degree", numRoots, rootdegree);CHKERRQ(ierr);
1063   /* Build a section / SFNode array of candidate points (face bd points) in the cone(support(leaf)),
1064      where each candidate is defined by a set of remote points (roots) for the other points that define the face. */
1065   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSection);CHKERRQ(ierr);
1066   ierr = PetscSectionSetChart(candidateSection, 0, numRoots);CHKERRQ(ierr);
1067   {
1068     PetscHMapIJ facehash;
1069 
1070     ierr = PetscHMapIJCreate(&facehash);CHKERRQ(ierr);
1071     for (l = 0; l < numLeaves; ++l) {
1072       const PetscInt    localPoint = localPoints[l];
1073       const PetscInt   *support;
1074       PetscInt          supportSize, s;
1075 
1076       if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]  Checking local point %D\n", rank, localPoint);CHKERRQ(ierr);}
1077       ierr = DMPlexGetSupportSize(dm, localPoint, &supportSize);CHKERRQ(ierr);
1078       ierr = DMPlexGetSupport(dm, localPoint, &support);CHKERRQ(ierr);
1079       for (s = 0; s < supportSize; ++s) {
1080         const PetscInt  face = support[s];
1081         const PetscInt *cone;
1082         PetscInt        coneSize, c, f, root;
1083         PetscBool       isFace = PETSC_TRUE;
1084 
1085         /* Only add face once */
1086         if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]    Support point %D\n", rank, face);CHKERRQ(ierr);}
1087         key.i = localPoint;
1088         key.j = face;
1089         ierr = PetscHMapIJGet(facehash, key, &f);CHKERRQ(ierr);
1090         if (f >= 0) continue;
1091         ierr = DMPlexGetConeSize(dm, face, &coneSize);CHKERRQ(ierr);
1092         ierr = DMPlexGetCone(dm, face, &cone);CHKERRQ(ierr);
1093         /* If a cone point does not map to leaves on any proc, then do not put face in SF */
1094         for (c = 0; c < coneSize; ++c) {
1095           if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]      Cone point %D\n", rank, cone[c]);CHKERRQ(ierr);}
1096           ierr = PetscHMapIGet(leafhash, cone[c], &root);CHKERRQ(ierr);
1097           if (!rootdegree[cone[c]] && (root < 0)) {isFace = PETSC_FALSE; break;}
1098         }
1099         if (isFace) {
1100           if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]    Found shared face %D\n", rank, face);CHKERRQ(ierr);}
1101           ierr = PetscHMapIJSet(facehash, key, l);CHKERRQ(ierr);
1102           ierr = PetscSectionAddDof(candidateSection, localPoint, coneSize);CHKERRQ(ierr);
1103         }
1104       }
1105     }
1106     if (debug) {ierr = PetscSynchronizedFlush(PetscObjectComm((PetscObject) dm), NULL);CHKERRQ(ierr);}
1107     ierr = PetscHMapIJClear(facehash);CHKERRQ(ierr);
1108     ierr = PetscSectionSetUp(candidateSection);CHKERRQ(ierr);
1109     ierr = PetscSectionGetStorageSize(candidateSection, &candidatesSize);CHKERRQ(ierr);
1110     ierr = PetscMalloc1(candidatesSize, &candidates);CHKERRQ(ierr);
1111     for (l = 0; l < numLeaves; ++l) {
1112       const PetscInt    localPoint = localPoints[l];
1113       const PetscInt   *support;
1114       PetscInt          supportSize, s, offset, idx = 0;
1115 
1116       if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]  Checking local point %D\n", rank, localPoint);CHKERRQ(ierr);}
1117       ierr = PetscSectionGetOffset(candidateSection, localPoint, &offset);CHKERRQ(ierr);
1118       ierr = DMPlexGetSupportSize(dm, localPoint, &supportSize);CHKERRQ(ierr);
1119       ierr = DMPlexGetSupport(dm, localPoint, &support);CHKERRQ(ierr);
1120       for (s = 0; s < supportSize; ++s) {
1121         const PetscInt  face = support[s];
1122         const PetscInt *cone;
1123         PetscInt        coneSize, c, f, root;
1124         PetscBool       isFace = PETSC_TRUE;
1125 
1126         /* Only add face once */
1127         if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]    Support point %D\n", rank, face);CHKERRQ(ierr);}
1128         key.i = localPoint;
1129         key.j = face;
1130         ierr = PetscHMapIJGet(facehash, key, &f);CHKERRQ(ierr);
1131         if (f >= 0) continue;
1132         ierr = DMPlexGetConeSize(dm, face, &coneSize);CHKERRQ(ierr);
1133         ierr = DMPlexGetCone(dm, face, &cone);CHKERRQ(ierr);
1134         /* If a cone point does not map to leaves on any proc, then do not put face in SF */
1135         for (c = 0; c < coneSize; ++c) {
1136           if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]      Cone point %D\n", rank, cone[c]);CHKERRQ(ierr);}
1137           ierr = PetscHMapIGet(leafhash, cone[c], &root);CHKERRQ(ierr);
1138           if (!rootdegree[cone[c]] && (root < 0)) {isFace = PETSC_FALSE; break;}
1139         }
1140         if (isFace) {
1141           if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]    Adding shared face %D at idx %D\n", rank, face, idx);CHKERRQ(ierr);}
1142           ierr = PetscHMapIJSet(facehash, key, l);CHKERRQ(ierr);
1143           candidates[offset+idx].rank    = -1;
1144           candidates[offset+idx++].index = coneSize-1;
1145           for (c = 0; c < coneSize; ++c) {
1146             if (cone[c] == localPoint) continue;
1147             if (rootdegree[cone[c]]) {
1148               candidates[offset+idx].rank    = rank;
1149               candidates[offset+idx++].index = cone[c];
1150             } else {
1151               ierr = PetscHMapIGet(leafhash, cone[c], &root);CHKERRQ(ierr);
1152               if (root < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot locate local point %D in SF", cone[c]);
1153               candidates[offset+idx++] = remotePoints[root];
1154             }
1155           }
1156         }
1157       }
1158     }
1159     if (debug) {ierr = PetscSynchronizedFlush(PetscObjectComm((PetscObject) dm), NULL);CHKERRQ(ierr);}
1160     ierr = PetscHMapIJDestroy(&facehash);CHKERRQ(ierr);
1161     ierr = PetscObjectViewFromOptions((PetscObject) candidateSection, NULL, "-petscsection_interp_candidate_view");CHKERRQ(ierr);
1162     ierr = SFNodeArrayViewFromOptions(PetscObjectComm((PetscObject) dm), "-petscsection_interp_candidate_view", "Candidates", NULL, candidatesSize, candidates);CHKERRQ(ierr);
1163   }
1164   /* Gather candidate section / array pair into the root partition via inverse(multi(pointSF)). */
1165   /*   Note that this section is indexed by offsets into leaves, not by point number */
1166   {
1167     PetscSF   sfMulti, sfInverse, sfCandidates;
1168     PetscInt *remoteOffsets;
1169 
1170     ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr);
1171     ierr = PetscSFCreateInverseSF(sfMulti, &sfInverse);CHKERRQ(ierr);
1172     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSectionRemote);CHKERRQ(ierr);
1173     ierr = PetscSFDistributeSection(sfInverse, candidateSection, &remoteOffsets, candidateSectionRemote);CHKERRQ(ierr);
1174     ierr = PetscSFCreateSectionSF(sfInverse, candidateSection, remoteOffsets, candidateSectionRemote, &sfCandidates);CHKERRQ(ierr);
1175     ierr = PetscSectionGetStorageSize(candidateSectionRemote, &candidatesRemoteSize);CHKERRQ(ierr);
1176     ierr = PetscMalloc1(candidatesRemoteSize, &candidatesRemote);CHKERRQ(ierr);
1177     ierr = PetscSFBcastBegin(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr);
1178     ierr = PetscSFBcastEnd(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr);
1179     ierr = PetscSFDestroy(&sfInverse);CHKERRQ(ierr);
1180     ierr = PetscSFDestroy(&sfCandidates);CHKERRQ(ierr);
1181     ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
1182 
1183     ierr = PetscObjectViewFromOptions((PetscObject) candidateSectionRemote, NULL, "-petscsection_interp_candidate_remote_view");CHKERRQ(ierr);
1184     ierr = SFNodeArrayViewFromOptions(PetscObjectComm((PetscObject) dm), "-petscsection_interp_candidate_remote_view", "Remote Candidates", NULL, candidatesRemoteSize, candidatesRemote);CHKERRQ(ierr);
1185   }
1186   /* */
1187   {
1188     PetscInt idx;
1189     /* There is a section point for every leaf attached to a given root point */
1190     for (r = 0, idx = 0; r < numRoots; ++r) {
1191       PetscInt deg;
1192       for (deg = 0; deg < rootdegree[r]; ++deg, ++idx) {
1193         PetscInt offset, dof, d;
1194 
1195         ierr = PetscSectionGetDof(candidateSectionRemote, idx, &dof);CHKERRQ(ierr);
1196         ierr = PetscSectionGetOffset(candidateSectionRemote, idx, &offset);CHKERRQ(ierr);
1197         for (d = 0; d < dof; ++d) {
1198           const PetscInt  sizeInd   = offset+d;
1199           const PetscInt  numPoints = candidatesRemote[sizeInd].index;
1200           const PetscInt *join      = NULL;
1201           PetscInt        points[1024], p, joinSize;
1202 
1203           points[0] = r;
1204           for (p = 0; p < numPoints; ++p) {
1205             ierr = DMPlexMapToLocalPoint(roothash, localPoints, rank, candidatesRemote[offset+(++d)], &points[p+1]);
1206             if (ierr) {d += numPoints-1 - p; break;} /* We got a point not in our overlap */
1207             if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]  Checking local candidate %D\n", rank, points[p+1]);CHKERRQ(ierr);}
1208           }
1209           if (ierr) continue;
1210           ierr = DMPlexGetJoin(dm, numPoints+1, points, &joinSize, &join);CHKERRQ(ierr);
1211           if (joinSize == 1) {
1212             if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]    Adding face %D at idx %D\n", rank, join[0], sizeInd);CHKERRQ(ierr);}
1213             candidatesRemote[sizeInd].rank  = rank;
1214             candidatesRemote[sizeInd].index = join[0];
1215           }
1216           ierr = DMPlexRestoreJoin(dm, numPoints+1, points, &joinSize, &join);CHKERRQ(ierr);
1217         }
1218       }
1219     }
1220     if (debug) {ierr = PetscSynchronizedFlush(PetscObjectComm((PetscObject) dm), NULL);CHKERRQ(ierr);}
1221   }
1222   /* Push claims back to receiver via the MultiSF and derive new pointSF mapping on receiver */
1223   {
1224     PetscSF         sfMulti, sfClaims, sfPointNew;
1225     PetscSFNode    *remotePointsNew;
1226     PetscHMapI      claimshash;
1227     PetscInt       *remoteOffsets, *localPointsNew;
1228     PetscInt        claimsSize, pStart, pEnd, root, numLocalNew, p, d;
1229 
1230     ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr);
1231     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &claimSection);CHKERRQ(ierr);
1232     ierr = PetscSFDistributeSection(sfMulti, candidateSectionRemote, &remoteOffsets, claimSection);CHKERRQ(ierr);
1233     ierr = PetscSFCreateSectionSF(sfMulti, candidateSectionRemote, remoteOffsets, claimSection, &sfClaims);CHKERRQ(ierr);
1234     ierr = PetscSectionGetStorageSize(claimSection, &claimsSize);CHKERRQ(ierr);
1235     ierr = PetscMalloc1(claimsSize, &claims);CHKERRQ(ierr);
1236     ierr = PetscSFBcastBegin(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr);
1237     ierr = PetscSFBcastEnd(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr);
1238     ierr = PetscSFDestroy(&sfClaims);CHKERRQ(ierr);
1239     ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
1240     ierr = PetscObjectViewFromOptions((PetscObject) claimSection, NULL, "-petscsection_interp_claim_view");CHKERRQ(ierr);
1241     ierr = SFNodeArrayViewFromOptions(PetscObjectComm((PetscObject) dm), "-petscsection_interp_claim_view", "Claims", NULL, claimsSize, claims);CHKERRQ(ierr);
1242     /* Walk the original section of local supports and add an SF entry for each updated item */
1243     ierr = PetscHMapICreate(&claimshash);CHKERRQ(ierr);
1244     for (p = 0; p < numRoots; ++p) {
1245       PetscInt dof, offset;
1246 
1247       if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]  Checking root for claims %D\n", rank, p);CHKERRQ(ierr);}
1248       ierr = PetscSectionGetDof(candidateSection, p, &dof);CHKERRQ(ierr);
1249       ierr = PetscSectionGetOffset(candidateSection, p, &offset);CHKERRQ(ierr);
1250       for (d = 0; d < dof;) {
1251         if (claims[offset+d].rank >= 0) {
1252           const PetscInt  faceInd   = offset+d;
1253           const PetscInt  numPoints = candidates[faceInd].index;
1254           const PetscInt *join      = NULL;
1255           PetscInt        joinSize, points[1024], c;
1256 
1257           if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]    Found claim for remote point (%D, %D)\n", rank, claims[faceInd].rank, claims[faceInd].index);CHKERRQ(ierr);}
1258           points[0] = p;
1259           if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]      point %D\n", rank, points[0]);CHKERRQ(ierr);}
1260           for (c = 0, ++d; c < numPoints; ++c, ++d) {
1261             key.i = candidates[offset+d].index;
1262             key.j = candidates[offset+d].rank;
1263             ierr = PetscHMapIJGet(roothash, key, &root);CHKERRQ(ierr);
1264             points[c+1] = localPoints[root];
1265             if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]      point %D\n", rank, points[c+1]);CHKERRQ(ierr);}
1266           }
1267           ierr = DMPlexGetJoin(dm, numPoints+1, points, &joinSize, &join);CHKERRQ(ierr);
1268           if (joinSize == 1) {
1269             if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]    Found local face %D\n", rank, join[0]);CHKERRQ(ierr);}
1270             ierr = PetscHMapISet(claimshash, join[0], faceInd);CHKERRQ(ierr);
1271           }
1272           ierr = DMPlexRestoreJoin(dm, numPoints+1, points, &joinSize, &join);CHKERRQ(ierr);
1273         } else d += claims[offset+d].index+1;
1274       }
1275     }
1276     if (debug) {ierr = PetscSynchronizedFlush(PetscObjectComm((PetscObject) dm), NULL);CHKERRQ(ierr);}
1277     /* Create new pointSF from hashed claims */
1278     ierr = PetscHMapIGetSize(claimshash, &numLocalNew);CHKERRQ(ierr);
1279     ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
1280     ierr = PetscMalloc1(numLeaves + numLocalNew, &localPointsNew);CHKERRQ(ierr);
1281     ierr = PetscMalloc1(numLeaves + numLocalNew, &remotePointsNew);CHKERRQ(ierr);
1282     for (p = 0; p < numLeaves; ++p) {
1283       localPointsNew[p] = localPoints[p];
1284       remotePointsNew[p].index = remotePoints[p].index;
1285       remotePointsNew[p].rank  = remotePoints[p].rank;
1286     }
1287     p = numLeaves;
1288     ierr = PetscHMapIGetKeys(claimshash, &p, localPointsNew);CHKERRQ(ierr);
1289     ierr = PetscSortInt(numLocalNew, &localPointsNew[numLeaves]);CHKERRQ(ierr);
1290     for (p = numLeaves; p < numLeaves + numLocalNew; ++p) {
1291       PetscInt offset;
1292       ierr = PetscHMapIGet(claimshash, localPointsNew[p], &offset);CHKERRQ(ierr);
1293       remotePointsNew[p] = claims[offset];
1294     }
1295     ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), &sfPointNew);CHKERRQ(ierr);
1296     ierr = PetscSFSetGraph(sfPointNew, pEnd-pStart, numLeaves+numLocalNew, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr);
1297     ierr = DMSetPointSF(dm, sfPointNew);CHKERRQ(ierr);
1298     ierr = PetscSFDestroy(&sfPointNew);CHKERRQ(ierr);
1299     ierr = PetscHMapIDestroy(&claimshash);CHKERRQ(ierr);
1300   }
1301   ierr = PetscHMapIDestroy(&leafhash);CHKERRQ(ierr);
1302   ierr = PetscHMapIJDestroy(&roothash);CHKERRQ(ierr);
1303   ierr = PetscSectionDestroy(&candidateSection);CHKERRQ(ierr);
1304   ierr = PetscSectionDestroy(&candidateSectionRemote);CHKERRQ(ierr);
1305   ierr = PetscSectionDestroy(&claimSection);CHKERRQ(ierr);
1306   ierr = PetscFree(candidates);CHKERRQ(ierr);
1307   ierr = PetscFree(candidatesRemote);CHKERRQ(ierr);
1308   ierr = PetscFree(claims);CHKERRQ(ierr);
1309   ierr = PetscLogEventEnd(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr);
1310   PetscFunctionReturn(0);
1311 }
1312 
1313 /*@C
1314   DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc.
1315 
1316   Collective on DM
1317 
1318   Input Parameters:
1319 + dm - The DMPlex object with only cells and vertices
1320 - dmInt - The interpolated DM
1321 
1322   Output Parameter:
1323 . dmInt - The complete DMPlex object
1324 
1325   Level: intermediate
1326 
1327   Notes:
1328     It does not copy over the coordinates.
1329 
1330 .keywords: mesh
1331 .seealso: DMPlexUninterpolate(), DMPlexCreateFromCellList(), DMPlexCopyCoordinates()
1332 @*/
1333 PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt)
1334 {
1335   DM             idm, odm = dm;
1336   PetscSF        sfPoint;
1337   PetscInt       depth, dim, d;
1338   const char    *name;
1339   PetscBool      flg=PETSC_TRUE;
1340   PetscErrorCode ierr;
1341 
1342   PetscFunctionBegin;
1343   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1344   PetscValidPointer(dmInt, 2);
1345   ierr = PetscLogEventBegin(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr);
1346   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1347   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1348   if ((depth == dim) || (dim <= 1)) {
1349     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
1350     idm  = dm;
1351   } else {
1352     for (d = 1; d < dim; ++d) {
1353       /* Create interpolated mesh */
1354       ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr);
1355       ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr);
1356       ierr = DMSetDimension(idm, dim);CHKERRQ(ierr);
1357       if (depth > 0) {
1358         ierr = DMPlexInterpolateFaces_Internal(odm, 1, idm);CHKERRQ(ierr);
1359         ierr = DMGetPointSF(odm, &sfPoint);CHKERRQ(ierr);
1360         ierr = DMPlexInterpolatePointSF(idm, sfPoint);CHKERRQ(ierr);
1361       }
1362       if (odm != dm) {ierr = DMDestroy(&odm);CHKERRQ(ierr);}
1363       odm = idm;
1364     }
1365     ierr = PetscObjectGetName((PetscObject) dm,  &name);CHKERRQ(ierr);
1366     ierr = PetscObjectSetName((PetscObject) idm,  name);CHKERRQ(ierr);
1367     ierr = DMPlexCopyCoordinates(dm, idm);CHKERRQ(ierr);
1368     ierr = DMCopyLabels(dm, idm);CHKERRQ(ierr);
1369     ierr = PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_interpolate_orient_interfaces", &flg, NULL);CHKERRQ(ierr);
1370     if (flg) {ierr = DMPlexOrientInterface(idm);CHKERRQ(ierr);}
1371   }
1372   {
1373     PetscBool            isper;
1374     const PetscReal      *maxCell, *L;
1375     const DMBoundaryType *bd;
1376 
1377     ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr);
1378     ierr = DMSetPeriodicity(idm,isper,maxCell,L,bd);CHKERRQ(ierr);
1379   }
1380   *dmInt = idm;
1381   ierr = PetscLogEventEnd(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr);
1382   PetscFunctionReturn(0);
1383 }
1384 
1385 /*@
1386   DMPlexCopyCoordinates - Copy coordinates from one mesh to another with the same vertices
1387 
1388   Collective on DM
1389 
1390   Input Parameter:
1391 . dmA - The DMPlex object with initial coordinates
1392 
1393   Output Parameter:
1394 . dmB - The DMPlex object with copied coordinates
1395 
1396   Level: intermediate
1397 
1398   Note: This is typically used when adding pieces other than vertices to a mesh
1399 
1400 .keywords: mesh
1401 .seealso: DMCopyLabels(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMGetCoordinateSection()
1402 @*/
1403 PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB)
1404 {
1405   Vec            coordinatesA, coordinatesB;
1406   VecType        vtype;
1407   PetscSection   coordSectionA, coordSectionB;
1408   PetscScalar   *coordsA, *coordsB;
1409   PetscInt       spaceDim, Nf, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d;
1410   PetscInt       cStartA, cEndA, cStartB, cEndB, cS, cE;
1411   PetscBool      lc = PETSC_FALSE;
1412   PetscErrorCode ierr;
1413 
1414   PetscFunctionBegin;
1415   PetscValidHeaderSpecific(dmA, DM_CLASSID, 1);
1416   PetscValidHeaderSpecific(dmB, DM_CLASSID, 2);
1417   if (dmA == dmB) PetscFunctionReturn(0);
1418   ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr);
1419   ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr);
1420   if ((vEndA-vStartA) != (vEndB-vStartB)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of vertices in first DM %d != %d in the second DM", vEndA-vStartA, vEndB-vStartB);
1421   ierr = DMPlexGetHeightStratum(dmA, 0, &cStartA, &cEndA);CHKERRQ(ierr);
1422   ierr = DMPlexGetHeightStratum(dmB, 0, &cStartB, &cEndB);CHKERRQ(ierr);
1423   ierr = DMGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr);
1424   ierr = DMGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr);
1425   if (coordSectionA == coordSectionB) PetscFunctionReturn(0);
1426   ierr = PetscSectionGetNumFields(coordSectionA, &Nf);CHKERRQ(ierr);
1427   if (!Nf) PetscFunctionReturn(0);
1428   if (Nf > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of coordinate fields must be 1, not %D", Nf);
1429   if (!coordSectionB) {
1430     PetscInt dim;
1431 
1432     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) coordSectionA), &coordSectionB);CHKERRQ(ierr);
1433     ierr = DMGetCoordinateDim(dmA, &dim);CHKERRQ(ierr);
1434     ierr = DMSetCoordinateSection(dmB, dim, coordSectionB);CHKERRQ(ierr);
1435     ierr = PetscObjectDereference((PetscObject) coordSectionB);CHKERRQ(ierr);
1436   }
1437   ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr);
1438   ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr);
1439   ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr);
1440   ierr = PetscSectionGetChart(coordSectionA, &cS, &cE);CHKERRQ(ierr);
1441   if (cStartA <= cS && cS < cEndA) { /* localized coordinates */
1442     if ((cEndA-cStartA) != (cEndB-cStartB)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of cells in first DM %D != %D in the second DM", cEndA-cStartA, cEndB-cStartB);
1443     cS = cS - cStartA + cStartB;
1444     cE = vEndB;
1445     lc = PETSC_TRUE;
1446   } else {
1447     cS = vStartB;
1448     cE = vEndB;
1449   }
1450   ierr = PetscSectionSetChart(coordSectionB, cS, cE);CHKERRQ(ierr);
1451   for (v = vStartB; v < vEndB; ++v) {
1452     ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr);
1453     ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr);
1454   }
1455   if (lc) { /* localized coordinates */
1456     PetscInt c;
1457 
1458     for (c = cS-cStartB; c < cEndB-cStartB; c++) {
1459       PetscInt dof;
1460 
1461       ierr = PetscSectionGetDof(coordSectionA, c + cStartA, &dof);CHKERRQ(ierr);
1462       ierr = PetscSectionSetDof(coordSectionB, c + cStartB, dof);CHKERRQ(ierr);
1463       ierr = PetscSectionSetFieldDof(coordSectionB, c + cStartB, 0, dof);CHKERRQ(ierr);
1464     }
1465   }
1466   ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr);
1467   ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr);
1468   ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr);
1469   ierr = VecCreate(PETSC_COMM_SELF, &coordinatesB);CHKERRQ(ierr);
1470   ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr);
1471   ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr);
1472   ierr = VecGetBlockSize(coordinatesA, &d);CHKERRQ(ierr);
1473   ierr = VecSetBlockSize(coordinatesB, d);CHKERRQ(ierr);
1474   ierr = VecGetType(coordinatesA, &vtype);CHKERRQ(ierr);
1475   ierr = VecSetType(coordinatesB, vtype);CHKERRQ(ierr);
1476   ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr);
1477   ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr);
1478   for (v = 0; v < vEndB-vStartB; ++v) {
1479     PetscInt offA, offB;
1480 
1481     ierr = PetscSectionGetOffset(coordSectionA, v + vStartA, &offA);CHKERRQ(ierr);
1482     ierr = PetscSectionGetOffset(coordSectionB, v + vStartB, &offB);CHKERRQ(ierr);
1483     for (d = 0; d < spaceDim; ++d) {
1484       coordsB[offB+d] = coordsA[offA+d];
1485     }
1486   }
1487   if (lc) { /* localized coordinates */
1488     PetscInt c;
1489 
1490     for (c = cS-cStartB; c < cEndB-cStartB; c++) {
1491       PetscInt dof, offA, offB;
1492 
1493       ierr = PetscSectionGetOffset(coordSectionA, c + cStartA, &offA);CHKERRQ(ierr);
1494       ierr = PetscSectionGetOffset(coordSectionB, c + cStartB, &offB);CHKERRQ(ierr);
1495       ierr = PetscSectionGetDof(coordSectionA, c + cStartA, &dof);CHKERRQ(ierr);
1496       ierr = PetscMemcpy(coordsB + offB,coordsA + offA,dof*sizeof(*coordsB));CHKERRQ(ierr);
1497     }
1498   }
1499   ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr);
1500   ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr);
1501   ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr);
1502   ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr);
1503   PetscFunctionReturn(0);
1504 }
1505 
1506 /*@
1507   DMPlexUninterpolate - Take in a mesh with all intermediate faces, edges, etc. and return a cell-vertex mesh
1508 
1509   Collective on DM
1510 
1511   Input Parameter:
1512 . dm - The complete DMPlex object
1513 
1514   Output Parameter:
1515 . dmUnint - The DMPlex object with only cells and vertices
1516 
1517   Level: intermediate
1518 
1519   Notes:
1520     It does not copy over the coordinates.
1521 
1522 .keywords: mesh
1523 .seealso: DMPlexInterpolate(), DMPlexCreateFromCellList(), DMPlexCopyCoordinates()
1524 @*/
1525 PetscErrorCode DMPlexUninterpolate(DM dm, DM *dmUnint)
1526 {
1527   DM             udm;
1528   PetscInt       dim, vStart, vEnd, cStart, cEnd, cMax, c, maxConeSize = 0, *cone;
1529   PetscErrorCode ierr;
1530 
1531   PetscFunctionBegin;
1532   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1533   PetscValidPointer(dmUnint, 2);
1534   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1535   if (dim <= 1) {
1536     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
1537     *dmUnint = dm;
1538     PetscFunctionReturn(0);
1539   }
1540   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
1541   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1542   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
1543   ierr = DMCreate(PetscObjectComm((PetscObject) dm), &udm);CHKERRQ(ierr);
1544   ierr = DMSetType(udm, DMPLEX);CHKERRQ(ierr);
1545   ierr = DMSetDimension(udm, dim);CHKERRQ(ierr);
1546   ierr = DMPlexSetChart(udm, cStart, vEnd);CHKERRQ(ierr);
1547   for (c = cStart; c < cEnd; ++c) {
1548     PetscInt *closure = NULL, closureSize, cl, coneSize = 0;
1549 
1550     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1551     for (cl = 0; cl < closureSize*2; cl += 2) {
1552       const PetscInt p = closure[cl];
1553 
1554       if ((p >= vStart) && (p < vEnd)) ++coneSize;
1555     }
1556     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1557     ierr = DMPlexSetConeSize(udm, c, coneSize);CHKERRQ(ierr);
1558     maxConeSize = PetscMax(maxConeSize, coneSize);
1559   }
1560   ierr = DMSetUp(udm);CHKERRQ(ierr);
1561   ierr = PetscMalloc1(maxConeSize, &cone);CHKERRQ(ierr);
1562   for (c = cStart; c < cEnd; ++c) {
1563     PetscInt *closure = NULL, closureSize, cl, coneSize = 0;
1564 
1565     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1566     for (cl = 0; cl < closureSize*2; cl += 2) {
1567       const PetscInt p = closure[cl];
1568 
1569       if ((p >= vStart) && (p < vEnd)) cone[coneSize++] = p;
1570     }
1571     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1572     ierr = DMPlexSetCone(udm, c, cone);CHKERRQ(ierr);
1573   }
1574   ierr = PetscFree(cone);CHKERRQ(ierr);
1575   ierr = DMPlexSetHybridBounds(udm, cMax, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1576   ierr = DMPlexSymmetrize(udm);CHKERRQ(ierr);
1577   ierr = DMPlexStratify(udm);CHKERRQ(ierr);
1578   /* Reduce SF */
1579   {
1580     PetscSF            sfPoint, sfPointUn;
1581     const PetscSFNode *remotePoints;
1582     const PetscInt    *localPoints;
1583     PetscSFNode       *remotePointsUn;
1584     PetscInt          *localPointsUn;
1585     PetscInt           vEnd, numRoots, numLeaves, l;
1586     PetscInt           numLeavesUn = 0, n = 0;
1587     PetscErrorCode     ierr;
1588 
1589     /* Get original SF information */
1590     ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
1591     ierr = DMGetPointSF(udm, &sfPointUn);CHKERRQ(ierr);
1592     ierr = DMPlexGetDepthStratum(dm, 0, NULL, &vEnd);CHKERRQ(ierr);
1593     ierr = PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr);
1594     /* Allocate space for cells and vertices */
1595     for (l = 0; l < numLeaves; ++l) if (localPoints[l] < vEnd) numLeavesUn++;
1596     /* Fill in leaves */
1597     if (vEnd >= 0) {
1598       ierr = PetscMalloc1(numLeavesUn, &remotePointsUn);CHKERRQ(ierr);
1599       ierr = PetscMalloc1(numLeavesUn, &localPointsUn);CHKERRQ(ierr);
1600       for (l = 0; l < numLeaves; l++) {
1601         if (localPoints[l] < vEnd) {
1602           localPointsUn[n]        = localPoints[l];
1603           remotePointsUn[n].rank  = remotePoints[l].rank;
1604           remotePointsUn[n].index = remotePoints[l].index;
1605           ++n;
1606         }
1607       }
1608       if (n != numLeavesUn) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent number of leaves %d != %d", n, numLeavesUn);
1609       ierr = PetscSFSetGraph(sfPointUn, vEnd, numLeavesUn, localPointsUn, PETSC_OWN_POINTER, remotePointsUn, PETSC_OWN_POINTER);CHKERRQ(ierr);
1610     }
1611   }
1612   {
1613     PetscBool            isper;
1614     const PetscReal      *maxCell, *L;
1615     const DMBoundaryType *bd;
1616 
1617     ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr);
1618     ierr = DMSetPeriodicity(udm,isper,maxCell,L,bd);CHKERRQ(ierr);
1619   }
1620 
1621   *dmUnint = udm;
1622   PetscFunctionReturn(0);
1623 }
1624