Actual source code: plexinterpolate.c

petsc-master 2019-12-03
Report Typos and Errors
  1:  #include <petsc/private/dmpleximpl.h>
  2:  #include <petsc/private/hashmapi.h>
  3:  #include <petsc/private/hashmapij.h>

  5: const char * const DMPlexInterpolatedFlags[] = {"none", "partial", "mixed", "full", "DMPlexInterpolatedFlag", "DMPLEX_INTERPOLATED_", 0};

  7: /* HashIJKL */

  9:  #include <petsc/private/hashmap.h>

 11: typedef struct _PetscHashIJKLKey { PetscInt i, j, k, l; } PetscHashIJKLKey;

 13: #define PetscHashIJKLKeyHash(key) \
 14:   PetscHashCombine(PetscHashCombine(PetscHashInt((key).i),PetscHashInt((key).j)), \
 15:                    PetscHashCombine(PetscHashInt((key).k),PetscHashInt((key).l)))

 17: #define PetscHashIJKLKeyEqual(k1,k2) \
 18:   (((k1).i==(k2).i) ? ((k1).j==(k2).j) ? ((k1).k==(k2).k) ? ((k1).l==(k2).l) : 0 : 0 : 0)

 20: PETSC_HASH_MAP(HashIJKL, PetscHashIJKLKey, PetscInt, PetscHashIJKLKeyHash, PetscHashIJKLKeyEqual, -1)

 22: static PetscSFNode _PetscInvalidSFNode = {-1, -1};

 24: typedef struct _PetscHashIJKLRemoteKey { PetscSFNode i, j, k, l; } PetscHashIJKLRemoteKey;

 26: #define PetscHashIJKLRemoteKeyHash(key) \
 27:   PetscHashCombine(PetscHashCombine(PetscHashInt((key).i.rank + (key).i.index),PetscHashInt((key).j.rank + (key).j.index)), \
 28:                    PetscHashCombine(PetscHashInt((key).k.rank + (key).k.index),PetscHashInt((key).l.rank + (key).l.index)))

 30: #define PetscHashIJKLRemoteKeyEqual(k1,k2) \
 31:   (((k1).i.rank==(k2).i.rank) ? ((k1).i.index==(k2).i.index) ? ((k1).j.rank==(k2).j.rank) ? ((k1).j.index==(k2).j.index) ? ((k1).k.rank==(k2).k.rank) ? ((k1).k.index==(k2).k.index) ? ((k1).l.rank==(k2).l.rank) ? ((k1).l.index==(k2).l.index) : 0 : 0 : 0 : 0 : 0 : 0 : 0)

 33: PETSC_HASH_MAP(HashIJKLRemote, PetscHashIJKLRemoteKey, PetscSFNode, PetscHashIJKLRemoteKeyHash, PetscHashIJKLRemoteKeyEqual, _PetscInvalidSFNode)

 35: static PetscErrorCode PetscSortSFNode(PetscInt n, PetscSFNode A[])
 36: {
 37:   PetscInt i;

 40:   for (i = 1; i < n; ++i) {
 41:     PetscSFNode x = A[i];
 42:     PetscInt    j;

 44:     for (j = i-1; j >= 0; --j) {
 45:       if ((A[j].rank > x.rank) || (A[j].rank == x.rank && A[j].index > x.index)) break;
 46:       A[j+1] = A[j];
 47:     }
 48:     A[j+1] = x;
 49:   }
 50:   return(0);
 51: }

 53: /*
 54:   DMPlexGetFaces_Internal - Gets groups of vertices that correspond to faces for the given cell
 55:   This assumes that the mesh is not interpolated from the depth of point p to the vertices
 56: */
 57: PetscErrorCode DMPlexGetFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[])
 58: {
 59:   const PetscInt *cone = NULL;
 60:   PetscInt        coneSize;
 61:   PetscErrorCode  ierr;

 65:   DMPlexGetConeSize(dm, p, &coneSize);
 66:   DMPlexGetCone(dm, p, &cone);
 67:   DMPlexGetRawFaces_Internal(dm, dim, coneSize, cone, numFaces, faceSize, faces);
 68:   return(0);
 69: }

 71: /*
 72:   DMPlexRestoreFaces_Internal - Restores the array
 73: */
 74: PetscErrorCode DMPlexRestoreFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[])
 75: {
 76:   PetscErrorCode  ierr;

 79:   if (faces) { DMRestoreWorkArray(dm, 0, MPIU_INT, (void *) faces); }
 80:   return(0);
 81: }

 83: /*
 84:   DMPlexGetRawFaces_Internal - Gets groups of vertices that correspond to faces for the given cone
 85: */
 86: PetscErrorCode DMPlexGetRawFaces_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[])
 87: {
 88:   PetscInt       *facesTmp;
 89:   PetscInt        maxConeSize, maxSupportSize;
 90:   PetscErrorCode  ierr;

 95:   DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);
 96:   if (faces) {DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), MPIU_INT, &facesTmp);}
 97:   switch (dim) {
 98:   case 1:
 99:     switch (coneSize) {
100:     case 2:
101:       if (faces) {
102:         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
103:         *faces = facesTmp;
104:       }
105:       if (numFaces) *numFaces = 2;
106:       if (faceSize) *faceSize = 1;
107:       break;
108:     default:
109:       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
110:     }
111:     break;
112:   case 2:
113:     switch (coneSize) {
114:     case 3:
115:       if (faces) {
116:         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
117:         facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
118:         facesTmp[4] = cone[2]; facesTmp[5] = cone[0];
119:         *faces = facesTmp;
120:       }
121:       if (numFaces) *numFaces = 3;
122:       if (faceSize) *faceSize = 2;
123:       break;
124:     case 4:
125:       /* Vertices follow right hand rule */
126:       if (faces) {
127:         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
128:         facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
129:         facesTmp[4] = cone[2]; facesTmp[5] = cone[3];
130:         facesTmp[6] = cone[3]; facesTmp[7] = cone[0];
131:         *faces = facesTmp;
132:       }
133:       if (numFaces) *numFaces = 4;
134:       if (faceSize) *faceSize = 2;
135:       break;
136:     default:
137:       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
138:     }
139:     break;
140:   case 3:
141:     switch (coneSize) {
142:     case 3:
143:       if (faces) {
144:         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
145:         facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
146:         facesTmp[4] = cone[2]; facesTmp[5] = cone[0];
147:         *faces = facesTmp;
148:       }
149:       if (numFaces) *numFaces = 3;
150:       if (faceSize) *faceSize = 2;
151:       break;
152:     case 4:
153:       /* Vertices of first face follow right hand rule and normal points away from last vertex */
154:       if (faces) {
155:         facesTmp[0] = cone[0]; facesTmp[1]  = cone[1]; facesTmp[2]  = cone[2];
156:         facesTmp[3] = cone[0]; facesTmp[4]  = cone[3]; facesTmp[5]  = cone[1];
157:         facesTmp[6] = cone[0]; facesTmp[7]  = cone[2]; facesTmp[8]  = cone[3];
158:         facesTmp[9] = cone[2]; facesTmp[10] = cone[1]; facesTmp[11] = cone[3];
159:         *faces = facesTmp;
160:       }
161:       if (numFaces) *numFaces = 4;
162:       if (faceSize) *faceSize = 3;
163:       break;
164:     case 8:
165:       /*  7--------6
166:          /|       /|
167:         / |      / |
168:        4--------5  |
169:        |  |     |  |
170:        |  |     |  |
171:        |  1--------2
172:        | /      | /
173:        |/       |/
174:        0--------3
175:        */
176:       if (faces) {
177:         facesTmp[0]  = cone[0]; facesTmp[1]  = cone[1]; facesTmp[2]  = cone[2]; facesTmp[3]  = cone[3]; /* Bottom */
178:         facesTmp[4]  = cone[4]; facesTmp[5]  = cone[5]; facesTmp[6]  = cone[6]; facesTmp[7]  = cone[7]; /* Top */
179:         facesTmp[8]  = cone[0]; facesTmp[9]  = cone[3]; facesTmp[10] = cone[5]; facesTmp[11] = cone[4]; /* Front */
180:         facesTmp[12] = cone[2]; facesTmp[13] = cone[1]; facesTmp[14] = cone[7]; facesTmp[15] = cone[6]; /* Back */
181:         facesTmp[16] = cone[3]; facesTmp[17] = cone[2]; facesTmp[18] = cone[6]; facesTmp[19] = cone[5]; /* Right */
182:         facesTmp[20] = cone[0]; facesTmp[21] = cone[4]; facesTmp[22] = cone[7]; facesTmp[23] = cone[1]; /* Left */
183:         *faces = facesTmp;
184:       }
185:       if (numFaces) *numFaces = 6;
186:       if (faceSize) *faceSize = 4;
187:       break;
188:     default:
189:       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
190:     }
191:     break;
192:   default:
193:     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim);
194:   }
195:   return(0);
196: }

198: /*
199:   DMPlexGetRawFacesHybrid_Internal - Gets groups of vertices that correspond to faces for the given cone using hybrid ordering (prisms)
200: */
201: static PetscErrorCode DMPlexGetRawFacesHybrid_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *numFacesNotH, PetscInt *faceSize, const PetscInt *faces[])
202: {
203:   PetscInt       *facesTmp;
204:   PetscInt        maxConeSize, maxSupportSize;
205:   PetscErrorCode  ierr;

210:   DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);
211:   if (faces) {DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), MPIU_INT, &facesTmp);}
212:   switch (dim) {
213:   case 1:
214:     switch (coneSize) {
215:     case 2:
216:       if (faces) {
217:         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
218:         *faces = facesTmp;
219:       }
220:       if (numFaces)     *numFaces = 2;
221:       if (numFacesNotH) *numFacesNotH = 2;
222:       if (faceSize)     *faceSize = 1;
223:       break;
224:     default:
225:       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
226:     }
227:     break;
228:   case 2:
229:     switch (coneSize) {
230:     case 4:
231:       if (faces) {
232:         facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
233:         facesTmp[2] = cone[2]; facesTmp[3] = cone[3];
234:         facesTmp[4] = cone[0]; facesTmp[5] = cone[2];
235:         facesTmp[6] = cone[1]; facesTmp[7] = cone[3];
236:         *faces = facesTmp;
237:       }
238:       if (numFaces)     *numFaces = 4;
239:       if (numFacesNotH) *numFacesNotH = 2;
240:       if (faceSize)     *faceSize = 2;
241:       break;
242:     default:
243:       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
244:     }
245:     break;
246:   case 3:
247:     switch (coneSize) {
248:     case 6: /* triangular prism */
249:       if (faces) {
250:         facesTmp[0]  = cone[0]; facesTmp[1]  = cone[1]; facesTmp[2]  = cone[2]; facesTmp[3]  = -1;      /* Bottom */
251:         facesTmp[4]  = cone[3]; facesTmp[5]  = cone[4]; facesTmp[6]  = cone[5]; facesTmp[7]  = -1;      /* Top */
252:         facesTmp[8]  = cone[0]; facesTmp[9]  = cone[1]; facesTmp[10] = cone[3]; facesTmp[11] = cone[4]; /* Back left */
253:         facesTmp[12] = cone[1]; facesTmp[13] = cone[2]; facesTmp[14] = cone[4]; facesTmp[15] = cone[5]; /* Back right */
254:         facesTmp[16] = cone[2]; facesTmp[17] = cone[0]; facesTmp[18] = cone[5]; facesTmp[19] = cone[3]; /* Front */
255:         *faces = facesTmp;
256:       }
257:       if (numFaces)     *numFaces = 5;
258:       if (numFacesNotH) *numFacesNotH = 2;
259:       if (faceSize)     *faceSize = -4;
260:       break;
261:     default:
262:       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
263:     }
264:     break;
265:   default:
266:     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim);
267:   }
268:   return(0);
269: }

271: static PetscErrorCode DMPlexRestoreRawFacesHybrid_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *numFacesNotH, PetscInt *faceSize, const PetscInt *faces[])
272: {
273:   PetscErrorCode  ierr;

276:   if (faces) { DMRestoreWorkArray(dm, 0, MPIU_INT, (void *) faces); }
277:   return(0);
278: }

280: static PetscErrorCode DMPlexGetFacesHybrid_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *numFacesNotH, PetscInt *faceSize, const PetscInt *faces[])
281: {
282:   const PetscInt *cone = NULL;
283:   PetscInt        coneSize;
284:   PetscErrorCode  ierr;

288:   DMPlexGetConeSize(dm, p, &coneSize);
289:   DMPlexGetCone(dm, p, &cone);
290:   DMPlexGetRawFacesHybrid_Internal(dm, dim, coneSize, cone, numFaces, numFacesNotH, faceSize, faces);
291:   return(0);
292: }

294: /* This interpolates faces for cells at some stratum */
295: static PetscErrorCode DMPlexInterpolateFaces_Internal(DM dm, PetscInt cellDepth, DM idm)
296: {
297:   DMLabel        subpointMap;
298:   PetscHashIJKL  faceTable;
299:   PetscInt      *pStart, *pEnd;
300:   PetscInt       cellDim, depth, faceDepth = cellDepth, numPoints = 0, faceSizeAll = 0, face, c, d;
301:   PetscInt       coneSizeH = 0, faceSizeAllH = 0, faceSizeAllT = 0, numCellFacesH = 0, faceT = 0, faceH, pMax = -1, dim, outerloop;
302:   PetscInt       cMax, fMax, eMax, vMax;

306:   DMGetDimension(dm, &cellDim);
307:   /* HACK: I need a better way to determine face dimension, or an alternative to GetFaces() */
308:   DMPlexGetSubpointMap(dm, &subpointMap);
309:   if (subpointMap) ++cellDim;
310:   DMPlexGetDepth(dm, &depth);
311:   ++depth;
312:   ++cellDepth;
313:   cellDim -= depth - cellDepth;
314:   PetscMalloc2(depth+1,&pStart,depth+1,&pEnd);
315:   for (d = depth-1; d >= faceDepth; --d) {
316:     DMPlexGetDepthStratum(dm, d, &pStart[d+1], &pEnd[d+1]);
317:   }
318:   DMPlexGetDepthStratum(dm, -1, NULL, &pStart[faceDepth]);
319:   pEnd[faceDepth] = pStart[faceDepth];
320:   for (d = faceDepth-1; d >= 0; --d) {
321:     DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);
322:   }
323:   cMax = fMax = eMax = vMax = PETSC_DETERMINE;
324:   DMGetDimension(dm, &dim);
325:   if (cellDim == dim) {
326:     DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);
327:     pMax = cMax;
328:   } else if (cellDim == dim -1) {
329:     DMPlexGetHybridBounds(dm, &cMax, &fMax, NULL, NULL);
330:     pMax = fMax;
331:   }
332:   pMax = pMax < 0 ? pEnd[cellDepth] : pMax;
333:   if (pMax < pEnd[cellDepth]) {
334:     const PetscInt *cellFaces, *cone;
335:     PetscInt        numCellFacesT, faceSize, cf;

337:     /* First get normal cell face size (we now allow hybrid cells to meet normal cells on either hybrid or normal faces */
338:     if (pStart[cellDepth] < pMax) {DMPlexGetFaces_Internal(dm, cellDim, pStart[cellDepth], NULL, &faceSizeAll, NULL);}

340:     DMPlexGetConeSize(dm, pMax, &coneSizeH);
341:     DMPlexGetCone(dm, pMax, &cone);
342:     DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSizeH, cone, &numCellFacesH, &numCellFacesT, &faceSize, &cellFaces);
343:     if (faceSize < 0) {
344:       PetscInt *sizes, minv, maxv;

346:       /* count vertices of hybrid and non-hybrid faces */
347:       PetscCalloc1(numCellFacesH, &sizes);
348:       for (cf = 0; cf < numCellFacesT; ++cf) { /* These are the non-hybrid faces */
349:         const PetscInt *cellFace = &cellFaces[-cf*faceSize];
350:         PetscInt       f;

352:         for (f = 0; f < -faceSize; ++f) sizes[cf] += (cellFace[f] >= 0 ? 1 : 0);
353:       }
354:       PetscSortInt(numCellFacesT, sizes);
355:       minv = sizes[0];
356:       maxv = sizes[PetscMax(numCellFacesT-1, 0)];
357:       if (minv != maxv) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Different number of vertices for non-hybrid face %D != %D", minv, maxv);
358:       faceSizeAllT = minv;
359:       PetscArrayzero(sizes, numCellFacesH);
360:       for (cf = numCellFacesT; cf < numCellFacesH; ++cf) { /* These are the hybrid faces */
361:         const PetscInt *cellFace = &cellFaces[-cf*faceSize];
362:         PetscInt       f;

364:         for (f = 0; f < -faceSize; ++f) sizes[cf-numCellFacesT] += (cellFace[f] >= 0 ? 1 : 0);
365:       }
366:       PetscSortInt(numCellFacesH - numCellFacesT, sizes);
367:       minv = sizes[0];
368:       maxv = sizes[PetscMax(numCellFacesH - numCellFacesT-1, 0)];
369:       PetscFree(sizes);
370:       if (minv != maxv) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Different number of vertices for hybrid face %D != %D", minv, maxv);
371:       faceSizeAllH = minv;
372:       if (!faceSizeAll) faceSizeAll = faceSizeAllT;
373:     } else { /* the size of the faces in hybrid cells is the same */
374:       faceSizeAll = faceSizeAllH = faceSizeAllT = faceSize;
375:     }
376:     DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSizeH, cone, &numCellFacesH, &numCellFacesT, &faceSize, &cellFaces);
377:   } else if (pEnd[cellDepth] > pStart[cellDepth]) {
378:     DMPlexGetFaces_Internal(dm, cellDim, pStart[cellDepth], NULL, &faceSizeAll, NULL);
379:     faceSizeAllH = faceSizeAllT = faceSizeAll;
380:   }
381:   if (faceSizeAll > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSizeAll);

383:   /* With hybrid grids, we first iterate on hybrid cells and start numbering the non-hybrid faces
384:      Then, faces for non-hybrid cells are numbered.
385:      This is to guarantee consistent orientations (all 0) of all the points in the cone of the hybrid cells */
386:   PetscHashIJKLCreate(&faceTable);
387:   for (outerloop = 0, face = pStart[faceDepth]; outerloop < 2; outerloop++) {
388:     PetscInt start, end;

390:     start = outerloop == 0 ? pMax : pStart[cellDepth];
391:     end = outerloop == 0 ? pEnd[cellDepth] : pMax;
392:     for (c = start; c < end; ++c) {
393:       const PetscInt *cellFaces;
394:       PetscInt        numCellFaces, faceSize, faceSizeInc, faceSizeCheck, cf;

396:       if (c < pMax) {
397:         DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);
398:         if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll);
399:         faceSizeCheck = faceSizeAll;
400:       } else { /* Hybrid cell */
401:         const PetscInt *cone;
402:         PetscInt        numCellFacesN, coneSize;

404:         DMPlexGetConeSize(dm, c, &coneSize);
405:         if (coneSize != coneSizeH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid coneSize %D != %D", coneSize, coneSizeH);
406:         DMPlexGetCone(dm, c, &cone);
407:         DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);
408:         if (numCellFaces != numCellFacesH) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected numCellFaces %D != %D for hybrid cell %D", numCellFaces, numCellFacesH, c);
409:         faceSize = PetscMax(faceSize, -faceSize);
410:         if (faceSize > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSize);
411:         numCellFaces = numCellFacesN; /* process only non-hybrid faces */
412:         faceSizeCheck = faceSizeAllT;
413:       }
414:       faceSizeInc = faceSize;
415:       for (cf = 0; cf < numCellFaces; ++cf) {
416:         const PetscInt   *cellFace = &cellFaces[cf*faceSizeInc];
417:         PetscInt          faceSizeH = faceSize;
418:         PetscHashIJKLKey  key;
419:         PetscHashIter     iter;
420:         PetscBool         missing;

422:         if (faceSizeInc == 2) {
423:           key.i = PetscMin(cellFace[0], cellFace[1]);
424:           key.j = PetscMax(cellFace[0], cellFace[1]);
425:           key.k = PETSC_MAX_INT;
426:           key.l = PETSC_MAX_INT;
427:         } else {
428:           key.i = cellFace[0];
429:           key.j = cellFace[1];
430:           key.k = cellFace[2];
431:           key.l = faceSize > 3 ? (cellFace[3] < 0 ? faceSizeH = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT;
432:           PetscSortInt(faceSize, (PetscInt *) &key);
433:         }
434:         /* this check is redundant for non-hybrid meshes */
435:         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);
436:         PetscHashIJKLPut(faceTable, key, &iter, &missing);
437:         if (missing) {
438:           PetscHashIJKLIterSet(faceTable, iter, face++);
439:           if (c >= pMax) ++faceT;
440:         }
441:       }
442:       if (c < pMax) {
443:         DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);
444:       } else {
445:         DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSizeH, NULL, NULL, NULL, NULL, &cellFaces);
446:       }
447:     }
448:   }
449:   pEnd[faceDepth] = face;

451:   /* Second pass for hybrid meshes: number hybrid faces */
452:   for (c = pMax; c < pEnd[cellDepth]; ++c) {
453:     const PetscInt *cellFaces, *cone;
454:     PetscInt        numCellFaces, numCellFacesN, faceSize, cf, coneSize;

456:     DMPlexGetConeSize(dm, c, &coneSize);
457:     DMPlexGetCone(dm, c, &cone);
458:     DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);
459:     if (numCellFaces != numCellFacesH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid numCellFaces %D != %D", numCellFaces, numCellFacesH);
460:     faceSize = PetscMax(faceSize, -faceSize);
461:     for (cf = numCellFacesN; cf < numCellFaces; ++cf) { /* These are the hybrid faces */
462:       const PetscInt   *cellFace = &cellFaces[cf*faceSize];
463:       PetscHashIJKLKey  key;
464:       PetscHashIter     iter;
465:       PetscBool         missing;
466:       PetscInt          faceSizeH = faceSize;

468:       if (faceSize == 2) {
469:         key.i = PetscMin(cellFace[0], cellFace[1]);
470:         key.j = PetscMax(cellFace[0], cellFace[1]);
471:         key.k = PETSC_MAX_INT;
472:         key.l = PETSC_MAX_INT;
473:       } else {
474:         key.i = cellFace[0];
475:         key.j = cellFace[1];
476:         key.k = cellFace[2];
477:         key.l = faceSize > 3 ? (cellFace[3] < 0 ? faceSizeH = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT;
478:         PetscSortInt(faceSize, (PetscInt *) &key);
479:       }
480:       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);
481:       PetscHashIJKLPut(faceTable, key, &iter, &missing);
482:       if (missing) {PetscHashIJKLIterSet(faceTable, iter, face++);}
483:     }
484:     DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);
485:   }
486:   faceH = face - pEnd[faceDepth];
487:   if (faceH) {
488:     if (fMax == PETSC_DETERMINE) fMax = pEnd[faceDepth];
489:     else if (eMax == PETSC_DETERMINE) eMax = pEnd[faceDepth];
490:     else SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of unassigned hybrid facets %D for cellDim %D and dimension %D", faceH, cellDim, dim);
491:   }
492:   pEnd[faceDepth] = face;
493:   PetscHashIJKLDestroy(&faceTable);
494:   /* Count new points */
495:   for (d = 0; d <= depth; ++d) {
496:     numPoints += pEnd[d]-pStart[d];
497:   }
498:   DMPlexSetChart(idm, 0, numPoints);
499:   /* Set cone sizes */
500:   for (d = 0; d <= depth; ++d) {
501:     PetscInt coneSize, p;

503:     if (d == faceDepth) {
504:       /* Now we have two cases: */
505:       if (faceSizeAll == faceSizeAllT) {
506:         /* I see no way to do this if we admit faces of different shapes */
507:         for (p = pStart[d]; p < pEnd[d]-faceH; ++p) {
508:           DMPlexSetConeSize(idm, p, faceSizeAll);
509:         }
510:         for (p = pEnd[d]-faceH; p < pEnd[d]; ++p) {
511:           DMPlexSetConeSize(idm, p, faceSizeAllH);
512:         }
513:       } else if (faceSizeAll == faceSizeAllH) {
514:         for (p = pStart[d]; p < pStart[d]+faceT; ++p) {
515:           DMPlexSetConeSize(idm, p, faceSizeAllT);
516:         }
517:         for (p = pStart[d]+faceT; p < pEnd[d]-faceH; ++p) {
518:           DMPlexSetConeSize(idm, p, faceSizeAll);
519:         }
520:         for (p = pEnd[d]-faceH; p < pEnd[d]; ++p) {
521:           DMPlexSetConeSize(idm, p, faceSizeAllH);
522:         }
523:       } else SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent faces sizes N: %D T: %D H: %D", faceSizeAll, faceSizeAllT, faceSizeAllH);
524:     } else if (d == cellDepth) {
525:       for (p = pStart[d]; p < pEnd[d]; ++p) {
526:         /* Number of cell faces may be different from number of cell vertices*/
527:         if (p < pMax) {
528:           DMPlexGetFaces_Internal(dm, cellDim, p, &coneSize, NULL, NULL);
529:         } else {
530:           DMPlexGetFacesHybrid_Internal(dm, cellDim, p, &coneSize, NULL, NULL, NULL);
531:         }
532:         DMPlexSetConeSize(idm, p, coneSize);
533:       }
534:     } else {
535:       for (p = pStart[d]; p < pEnd[d]; ++p) {
536:         DMPlexGetConeSize(dm, p, &coneSize);
537:         DMPlexSetConeSize(idm, p, coneSize);
538:       }
539:     }
540:   }
541:   DMSetUp(idm);
542:   /* Get face cones from subsets of cell vertices */
543:   if (faceSizeAll > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSizeAll);
544:   PetscHashIJKLCreate(&faceTable);
545:   for (d = depth; d > cellDepth; --d) {
546:     const PetscInt *cone;
547:     PetscInt        p;

549:     for (p = pStart[d]; p < pEnd[d]; ++p) {
550:       DMPlexGetCone(dm, p, &cone);
551:       DMPlexSetCone(idm, p, cone);
552:       DMPlexGetConeOrientation(dm, p, &cone);
553:       DMPlexSetConeOrientation(idm, p, cone);
554:     }
555:   }
556:   for (outerloop = 0, face = pStart[faceDepth]; outerloop < 2; outerloop++) {
557:     PetscInt start, end;

559:     start = outerloop == 0 ? pMax : pStart[cellDepth];
560:     end = outerloop == 0 ? pEnd[cellDepth] : pMax;
561:     for (c = start; c < end; ++c) {
562:       const PetscInt *cellFaces;
563:       PetscInt        numCellFaces, faceSize, faceSizeInc, cf;

565:       if (c < pMax) {
566:         DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);
567:         if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll);
568:       } else {
569:         const PetscInt *cone;
570:         PetscInt        numCellFacesN, coneSize;

572:         DMPlexGetConeSize(dm, c, &coneSize);
573:         if (coneSize != coneSizeH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid coneSize %D != %D", coneSize, coneSizeH);
574:         DMPlexGetCone(dm, c, &cone);
575:         DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);
576:         if (numCellFaces != numCellFacesH) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected numCellFaces %D != %D for hybrid cell %D", numCellFaces, numCellFacesH, c);
577:         faceSize = PetscMax(faceSize, -faceSize);
578:         if (faceSize > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSize);
579:         numCellFaces = numCellFacesN; /* process only non-hybrid faces */
580:       }
581:       faceSizeInc = faceSize;
582:       for (cf = 0; cf < numCellFaces; ++cf) {
583:         const PetscInt  *cellFace = &cellFaces[cf*faceSizeInc];
584:         PetscHashIJKLKey key;
585:         PetscHashIter    iter;
586:         PetscBool        missing;

588:         if (faceSizeInc == 2) {
589:           key.i = PetscMin(cellFace[0], cellFace[1]);
590:           key.j = PetscMax(cellFace[0], cellFace[1]);
591:           key.k = PETSC_MAX_INT;
592:           key.l = PETSC_MAX_INT;
593:         } else {
594:           key.i = cellFace[0];
595:           key.j = cellFace[1];
596:           key.k = cellFace[2];
597:           key.l = faceSizeInc > 3 ? (cellFace[3] < 0 ? faceSize = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT;
598:           PetscSortInt(faceSizeInc, (PetscInt *) &key);
599:         }
600:         PetscHashIJKLPut(faceTable, key, &iter, &missing);
601:         if (missing) {
602:           DMPlexSetCone(idm, face, cellFace);
603:           PetscHashIJKLIterSet(faceTable, iter, face);
604:           DMPlexInsertCone(idm, c, cf, face++);
605:         } else {
606:           const PetscInt *cone;
607:           PetscInt        coneSize, ornt, i, j, f;

609:           PetscHashIJKLIterGet(faceTable, iter, &f);
610:           DMPlexInsertCone(idm, c, cf, f);
611:           /* Orient face: Do not allow reverse orientation at the first vertex */
612:           DMPlexGetConeSize(idm, f, &coneSize);
613:           DMPlexGetCone(idm, f, &cone);
614:           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);
615:           /* - First find the initial vertex */
616:           for (i = 0; i < faceSize; ++i) if (cellFace[0] == cone[i]) break;
617:           /* - Try forward comparison */
618:           for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+j)%faceSize]) break;
619:           if (j == faceSize) {
620:             if ((faceSize == 2) && (i == 1)) ornt = -2;
621:             else                             ornt = i;
622:           } else {
623:             /* - Try backward comparison */
624:             for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+faceSize-j)%faceSize]) break;
625:             if (j == faceSize) {
626:               if (i == 0) ornt = -faceSize;
627:               else        ornt = -i;
628:             } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine orientation of face %D in cell %D", f, c);
629:           }
630:           DMPlexInsertConeOrientation(idm, c, cf, ornt);
631:         }
632:       }
633:       if (c < pMax) {
634:         DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);
635:       } else {
636:         DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSizeH, NULL, NULL, NULL, NULL, &cellFaces);
637:       }
638:     }
639:   }
640:   /* Second pass for hybrid meshes: orient hybrid faces */
641:   for (c = pMax; c < pEnd[cellDepth]; ++c) {
642:     const PetscInt *cellFaces, *cone;
643:     PetscInt        numCellFaces, numCellFacesN, faceSize, cf, coneSize;

645:     DMPlexGetConeSize(dm, c, &coneSize);
646:     DMPlexGetCone(dm, c, &cone);
647:     DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);
648:     if (numCellFaces != numCellFacesH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid numCellFaces %D != %D", numCellFaces, numCellFacesH);
649:     faceSize = PetscMax(faceSize, -faceSize);
650:     for (cf = numCellFacesN; cf < numCellFaces; ++cf) { /* These are the hybrid faces */
651:       const PetscInt   *cellFace = &cellFaces[cf*faceSize];
652:       PetscHashIJKLKey key;
653:       PetscHashIter    iter;
654:       PetscBool        missing;
655:       PetscInt         faceSizeH = faceSize;

657:       if (faceSize == 2) {
658:         key.i = PetscMin(cellFace[0], cellFace[1]);
659:         key.j = PetscMax(cellFace[0], cellFace[1]);
660:         key.k = PETSC_MAX_INT;
661:         key.l = PETSC_MAX_INT;
662:       } else {
663:         key.i = cellFace[0];
664:         key.j = cellFace[1];
665:         key.k = cellFace[2];
666:         key.l = faceSize > 3 ? (cellFace[3] < 0 ? faceSizeH = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT;
667:         PetscSortInt(faceSize, (PetscInt *) &key);
668:       }
669:       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);
670:       PetscHashIJKLPut(faceTable, key, &iter, &missing);
671:       if (missing) {
672:         DMPlexSetCone(idm, face, cellFace);
673:         PetscHashIJKLIterSet(faceTable, iter, face);
674:         DMPlexInsertCone(idm, c, cf, face++);
675:       } else {
676:         PetscInt        fv[4] = {0, 1, 2, 3};
677:         const PetscInt *cone;
678:         PetscInt        coneSize, ornt, i, j, f;
679:         PetscBool       q2h = PETSC_FALSE;

681:         PetscHashIJKLIterGet(faceTable, iter, &f);
682:         DMPlexInsertCone(idm, c, cf, f);
683:         /* Orient face: Do not allow reverse orientation at the first vertex */
684:         DMPlexGetConeSize(idm, f, &coneSize);
685:         DMPlexGetCone(idm, f, &cone);
686:         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);
687:         /* Hybrid faces are stored as tensor products of edges, so to compare them to normal faces, we have to flip */
688:         if (faceSize == 4 && c >= pMax && faceSizeAll != faceSizeAllT && f < pEnd[faceDepth] - faceH) {q2h = PETSC_TRUE; fv[2] = 3; fv[3] = 2;}
689:         /* - First find the initial vertex */
690:         for (i = 0; i < faceSizeH; ++i) if (cellFace[fv[0]] == cone[i]) break;
691:         if (q2h) { /* Matt's case: hybrid faces meeting with non-hybrid faces. This is a case that is not (and will not be) supported in general by the refinements */
692:           /* - Try forward comparison */
693:           for (j = 0; j < faceSizeH; ++j) if (cellFace[fv[j]] != cone[(i+j)%faceSizeH]) break;
694:           if (j == faceSizeH) {
695:             if ((faceSizeH == 2) && (i == 1)) ornt = -2;
696:             else                              ornt = i;
697:           } else {
698:             /* - Try backward comparison */
699:             for (j = 0; j < faceSizeH; ++j) if (cellFace[fv[j]] != cone[(i+faceSizeH-j)%faceSizeH]) break;
700:             if (j == faceSizeH) {
701:               if (i == 0) ornt = -faceSizeH;
702:               else        ornt = -i;
703:             } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine orientation of face %D in cell %D", f, c);
704:           }
705:         } else {
706:           /* when matching hybrid faces in 3D, only few cases are possible.
707:              Face traversal however can no longer follow the usual convention, this seems a serious issue to me */
708:           PetscInt tquad_map[4][4] = { {PETSC_MIN_INT,            0,PETSC_MIN_INT,PETSC_MIN_INT},
709:                                        {           -1,PETSC_MIN_INT,PETSC_MIN_INT,PETSC_MIN_INT},
710:                                        {PETSC_MIN_INT,PETSC_MIN_INT,PETSC_MIN_INT,            1},
711:                                        {PETSC_MIN_INT,PETSC_MIN_INT,           -2,PETSC_MIN_INT} };
712:           PetscInt i2;

714:           /* find the second vertex */
715:           for (i2 = 0; i2 < faceSizeH; ++i2) if (cellFace[fv[1]] == cone[i2]) break;
716:           switch (faceSizeH) {
717:           case 2:
718:             ornt = i ? -2 : 0;
719:             break;
720:           case 4:
721:             ornt = tquad_map[i][i2];
722:             break;
723:           default:
724:             SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unhandled face size %D for face %D in cell %D", faceSizeH, f, c);

726:           }
727:         }
728:         DMPlexInsertConeOrientation(idm, c, cf, ornt);
729:       }
730:     }
731:     DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);
732:   }
733:   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]);
734:   PetscFree2(pStart,pEnd);
735:   PetscHashIJKLDestroy(&faceTable);
736:   PetscFree2(pStart,pEnd);
737:   DMPlexSetHybridBounds(idm, cMax, fMax, eMax, vMax);
738:   DMPlexSymmetrize(idm);
739:   DMPlexStratify(idm);
740:   return(0);
741: }

743: static PetscErrorCode SortRmineRremoteByRemote_Private(PetscSF sf, PetscInt *rmine1[], PetscInt *rremote1[])
744: {
745:   PetscInt            nleaves;
746:   PetscInt            nranks;
747:   const PetscMPIInt  *ranks=NULL;
748:   const PetscInt     *roffset=NULL, *rmine=NULL, *rremote=NULL;
749:   PetscInt            n, o, r;
750:   PetscErrorCode      ierr;

753:   PetscSFGetRootRanks(sf, &nranks, &ranks, &roffset, &rmine, &rremote);
754:   nleaves = roffset[nranks];
755:   PetscMalloc2(nleaves, rmine1, nleaves, rremote1);
756:   for (r=0; r<nranks; r++) {
757:     /* simultaneously sort rank-wise portions of rmine & rremote by values in rremote
758:        - to unify order with the other side */
759:     o = roffset[r];
760:     n = roffset[r+1] - o;
761:     PetscArraycpy(&(*rmine1)[o], &rmine[o], n);
762:     PetscArraycpy(&(*rremote1)[o], &rremote[o], n);
763:     PetscSortIntWithArray(n, &(*rremote1)[o], &(*rmine1)[o]);
764:   }
765:   return(0);
766: }

768: PetscErrorCode DMPlexOrientInterface_Internal(DM dm)
769: {
770:   /* Here we only compare first 2 points of the cone. Full cone size would lead to stronger self-checking. */
771:   PetscInt          masterCone[2];
772:   PetscInt          (*roots)[2], (*leaves)[2];
773:   PetscMPIInt       (*rootsRanks)[2], (*leavesRanks)[2];

775:   PetscSF           sf=NULL;
776:   const PetscInt    *locals=NULL;
777:   const PetscSFNode *remotes=NULL;
778:   PetscInt           nroots, nleaves, p, c;
779:   PetscInt           nranks, n, o, r;
780:   const PetscMPIInt *ranks=NULL;
781:   const PetscInt    *roffset=NULL;
782:   PetscInt          *rmine1=NULL, *rremote1=NULL; /* rmine and rremote copies simultaneously sorted by rank and rremote */
783:   const PetscInt    *cone=NULL;
784:   PetscInt           coneSize, ind0;
785:   MPI_Comm           comm;
786:   PetscMPIInt        rank,size;
787:   PetscInt           debug = 0;
788:   PetscErrorCode     ierr;

791:   DMGetPointSF(dm, &sf);
792:   PetscSFGetGraph(sf, &nroots, &nleaves, &locals, &remotes);
793:   if (nroots < 0) return(0);
794:   PetscSFSetUp(sf);
795:   PetscSFGetRootRanks(sf, &nranks, &ranks, &roffset, NULL, NULL);
796:   DMViewFromOptions(dm, NULL, "-before_fix_dm_view");
797: #if defined(PETSC_USE_DEBUG)
798:   DMPlexCheckPointSF(dm);
799: #endif
800:   SortRmineRremoteByRemote_Private(sf, &rmine1, &rremote1);
801:   PetscMalloc4(nroots, &roots, nroots, &leaves, nroots, &rootsRanks, nroots, &leavesRanks);
802:   PetscObjectGetComm((PetscObject) dm, &comm);
803:   MPI_Comm_rank(comm, &rank);
804:   MPI_Comm_size(comm, &size);
805:   for (p = 0; p < nroots; ++p) {
806:     DMPlexGetConeSize(dm, p, &coneSize);
807:     DMPlexGetCone(dm, p, &cone);
808:     if (coneSize < 2) {
809:       for (c = 0; c < 2; c++) {
810:         roots[p][c] = -1;
811:         rootsRanks[p][c] = -1;
812:       }
813:       continue;
814:     }
815:     /* Translate all points to root numbering */
816:     for (c = 0; c < 2; c++) {
817:       PetscFindInt(cone[c], nleaves, locals, &ind0);
818:       if (ind0 < 0) {
819:         roots[p][c] = cone[c];
820:         rootsRanks[p][c] = rank;
821:       } else {
822:         roots[p][c] = remotes[ind0].index;
823:         rootsRanks[p][c] = remotes[ind0].rank;
824:       }
825:     }
826:   }
827:   for (p = 0; p < nroots; ++p) {
828:     for (c = 0; c < 2; c++) {
829:       leaves[p][c] = -2;
830:       leavesRanks[p][c] = -2;
831:     }
832:   }
833:   PetscSFBcastBegin(sf, MPIU_2INT, roots, leaves);
834:   PetscSFBcastBegin(sf, MPI_2INT, rootsRanks, leavesRanks);
835:   PetscSFBcastEnd(sf, MPIU_2INT, roots, leaves);
836:   PetscSFBcastEnd(sf, MPI_2INT, rootsRanks, leavesRanks);
837:   if (debug) {PetscSynchronizedFlush(comm, NULL);}
838:   if (debug && rank == 0) {PetscSynchronizedPrintf(comm, "Referenced roots\n");}
839:   for (p = 0; p < nroots; ++p) {
840:     if (leaves[p][0] < 0) continue;
841:     DMPlexGetConeSize(dm, p, &coneSize);
842:     DMPlexGetCone(dm, p, &cone);
843:     if (debug) {PetscSynchronizedPrintf(comm, "[%d]  %4D: cone=[%4D %4D] roots=[(%d,%4D) (%d,%4D)] leaves=[(%d,%4D) (%d,%4D)]", rank, p, cone[0], cone[1], rootsRanks[p][0], roots[p][0], rootsRanks[p][1], roots[p][1], leavesRanks[p][0], leaves[p][0], leavesRanks[p][1], leaves[p][1]);}
844:     if ((leaves[p][0] != roots[p][0]) || (leaves[p][1] != roots[p][1]) || (leavesRanks[p][0] != rootsRanks[p][0]) || (leavesRanks[p][1] != rootsRanks[p][1])) {
845:       /* Translate these two leaves to my cone points; masterCone means desired order p's cone points */
846:       for (c = 0; c < 2; c++) {
847:         if (leavesRanks[p][c] == rank) {
848:           /* A local leave is just taken as it is */
849:           masterCone[c] = leaves[p][c];
850:           continue;
851:         }
852:         /* Find index of rank leavesRanks[p][c] among remote ranks */
853:         /* No need for PetscMPIIntCast because these integers were originally cast from PetscMPIInt. */
854:         PetscFindMPIInt((PetscMPIInt)leavesRanks[p][c], nranks, ranks, &r);
855:         if (PetscUnlikely(r < 0)) SETERRQ7(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D cone[%D]=%D root (%d,%D) leave (%d,%D): leave rank not found among remote ranks",p,c,cone[c],rootsRanks[p][c],roots[p][c],leavesRanks[p][c],leaves[p][c]);
856:         if (PetscUnlikely(ranks[r] < 0 || ranks[r] >= size)) SETERRQ5(PETSC_COMM_SELF, PETSC_ERR_PLIB, "p=%D c=%D commsize=%d: ranks[%D] = %d makes no sense",p,c,size,r,ranks[r]);
857:         /* Find point leaves[p][c] among remote points aimed at rank leavesRanks[p][c] */
858:         o = roffset[r];
859:         n = roffset[r+1] - o;
860:         PetscFindInt(leaves[p][c], n, &rremote1[o], &ind0);
861:         if (PetscUnlikely(ind0 < 0)) SETERRQ7(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D cone[%D]=%D root (%d,%D) leave (%d,%D): corresponding remote point not found - it seems there is missing connection in point SF!",p,c,cone[c],rootsRanks[p][c],roots[p][c],leavesRanks[p][c],leaves[p][c]);
862:         /* Get the corresponding local point */
863:         masterCone[c] = rmine1[o+ind0];
864:       }
865:       if (debug) {PetscSynchronizedPrintf(comm, " masterCone=[%4D %4D]\n", masterCone[0], masterCone[1]);}
866:       /* Set the desired order of p's cone points and fix orientations accordingly */
867:       /* Vaclav's note: Here we only compare first 2 points of the cone. Full cone size would lead to stronger self-checking. */
868:       DMPlexOrientCell(dm, p, 2, masterCone);
869:     } else if (debug) {PetscSynchronizedPrintf(comm, " ==\n");}
870:   }
871:   if (debug) {
872:     PetscSynchronizedFlush(comm, NULL);
873:     MPI_Barrier(comm);
874:   }
875:   DMViewFromOptions(dm, NULL, "-after_fix_dm_view");
876:   PetscFree4(roots, leaves, rootsRanks, leavesRanks);
877:   PetscFree2(rmine1, rremote1);
878:   return(0);
879: }

881: static PetscErrorCode IntArrayViewFromOptions(MPI_Comm comm, const char opt[], const char name[], const char idxname[], const char valname[], PetscInt n, const PetscInt a[])
882: {
883:   PetscInt       idx;
884:   PetscMPIInt    rank;
885:   PetscBool      flg;

889:   PetscOptionsHasName(NULL, NULL, opt, &flg);
890:   if (!flg) return(0);
891:   MPI_Comm_rank(comm, &rank);
892:   PetscSynchronizedPrintf(comm, "[%d]%s:\n", rank, name);
893:   for (idx = 0; idx < n; ++idx) {PetscSynchronizedPrintf(comm, "[%d]%s %D %s %D\n", rank, idxname, idx, valname, a[idx]);}
894:   PetscSynchronizedFlush(comm, NULL);
895:   return(0);
896: }

898: static PetscErrorCode SFNodeArrayViewFromOptions(MPI_Comm comm, const char opt[], const char name[], const char idxname[], PetscInt n, const PetscSFNode a[])
899: {
900:   PetscInt       idx;
901:   PetscMPIInt    rank;
902:   PetscBool      flg;

906:   PetscOptionsHasName(NULL, NULL, opt, &flg);
907:   if (!flg) return(0);
908:   MPI_Comm_rank(comm, &rank);
909:   PetscSynchronizedPrintf(comm, "[%d]%s:\n", rank, name);
910:   if (idxname) {
911:     for (idx = 0; idx < n; ++idx) {PetscSynchronizedPrintf(comm, "[%d]%s %D rank %D index %D\n", rank, idxname, idx, a[idx].rank, a[idx].index);}
912:   } else {
913:     for (idx = 0; idx < n; ++idx) {PetscSynchronizedPrintf(comm, "[%d]rank %D index %D\n", rank, a[idx].rank, a[idx].index);}
914:   }
915:   PetscSynchronizedFlush(comm, NULL);
916:   return(0);
917: }

919: static PetscErrorCode DMPlexMapToLocalPoint(DM dm, PetscHMapIJ remotehash, PetscSFNode remotePoint, PetscInt *localPoint)
920: {
921:   PetscSF         sf;
922:   const PetscInt *locals;
923:   PetscMPIInt     rank;
924:   PetscErrorCode  ierr;

927:   MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);
928:   DMGetPointSF(dm, &sf);
929:   PetscSFGetGraph(sf, NULL, NULL, &locals, NULL);
930:   if (remotePoint.rank == rank) {
931:     *localPoint = remotePoint.index;
932:   } else {
933:     PetscHashIJKey key;
934:     PetscInt       l;

936:     key.i = remotePoint.index;
937:     key.j = remotePoint.rank;
938:     PetscHMapIJGet(remotehash, key, &l);
939:     if (l >= 0) {
940:       *localPoint = locals[l];
941:     } else PetscFunctionReturn(1);
942:   }
943:   return(0);
944: }

946: static PetscErrorCode DMPlexMapToGlobalPoint(DM dm, PetscInt localPoint, PetscSFNode *remotePoint)
947: {
948:   PetscSF            sf;
949:   const PetscInt    *locals, *rootdegree;
950:   const PetscSFNode *remotes;
951:   PetscInt           Nl, l;
952:   PetscMPIInt        rank;
953:   PetscErrorCode     ierr;

956:   MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);
957:   DMGetPointSF(dm, &sf);
958:   PetscSFGetGraph(sf, NULL, &Nl, &locals, &remotes);
959:   if (Nl < 0) goto owned;
960:   PetscSFComputeDegreeBegin(sf, &rootdegree);
961:   PetscSFComputeDegreeEnd(sf, &rootdegree);
962:   if (rootdegree[localPoint]) goto owned;
963:   PetscFindInt(localPoint, Nl, locals, &l);
964:   if (l < 0) PetscFunctionReturn(1);
965:   *remotePoint = remotes[l];
966:   return(0);
967:   owned:
968:   remotePoint->rank  = rank;
969:   remotePoint->index = localPoint;
970:   return(0);
971: }


974: static PetscErrorCode DMPlexPointIsShared(DM dm, PetscInt p, PetscBool *isShared)
975: {
976:   PetscSF         sf;
977:   const PetscInt *locals, *rootdegree;
978:   PetscInt        Nl, idx;
979:   PetscErrorCode  ierr;

982:   *isShared = PETSC_FALSE;
983:   DMGetPointSF(dm, &sf);
984:   PetscSFGetGraph(sf, NULL, &Nl, &locals, NULL);
985:   if (Nl < 0) return(0);
986:   PetscFindInt(p, Nl, locals, &idx);
987:   if (idx >= 0) {*isShared = PETSC_TRUE; return(0);}
988:   PetscSFComputeDegreeBegin(sf, &rootdegree);
989:   PetscSFComputeDegreeEnd(sf, &rootdegree);
990:   if (rootdegree[p] > 0) *isShared = PETSC_TRUE;
991:   return(0);
992: }

994: static PetscErrorCode DMPlexConeIsShared(DM dm, PetscInt p, PetscBool *isShared)
995: {
996:   const PetscInt *cone;
997:   PetscInt        coneSize, c;
998:   PetscBool       cShared = PETSC_TRUE;
999:   PetscErrorCode  ierr;

1002:   DMPlexGetConeSize(dm, p, &coneSize);
1003:   DMPlexGetCone(dm, p, &cone);
1004:   for (c = 0; c < coneSize; ++c) {
1005:     PetscBool pointShared;

1007:     DMPlexPointIsShared(dm, cone[c], &pointShared);
1008:     cShared = (PetscBool) (cShared && pointShared);
1009:   }
1010:   *isShared = coneSize ? cShared : PETSC_FALSE;
1011:   return(0);
1012: }

1014: static PetscErrorCode DMPlexGetConeMinimum(DM dm, PetscInt p, PetscSFNode *cpmin)
1015: {
1016:   const PetscInt *cone;
1017:   PetscInt        coneSize, c;
1018:   PetscSFNode     cmin = {PETSC_MAX_INT, PETSC_MAX_INT}, missing = {-1, -1};
1019:   PetscErrorCode  ierr;

1022:   DMPlexGetConeSize(dm, p, &coneSize);
1023:   DMPlexGetCone(dm, p, &cone);
1024:   for (c = 0; c < coneSize; ++c) {
1025:     PetscSFNode rcp;

1027:     DMPlexMapToGlobalPoint(dm, cone[c], &rcp);
1028:     if (ierr) {
1029:       cmin = missing;
1030:     } else {
1031:       cmin = (rcp.rank < cmin.rank) || (rcp.rank == cmin.rank && rcp.index < cmin.index) ? rcp : cmin;
1032:     }
1033:   }
1034:   *cpmin = coneSize ? cmin : missing;
1035:   return(0);
1036: }

1038: /*
1039:   Each shared face has an entry in the candidates array:
1040:     (-1, coneSize-1), {(global cone point)}
1041:   where the set is missing the point p which we use as the key for the face
1042: */
1043: static PetscErrorCode DMPlexAddSharedFace_Private(DM dm, PetscSection candidateSection, PetscSFNode candidates[], PetscHMapIJ faceHash, PetscInt p, PetscBool debug)
1044: {
1045:   MPI_Comm        comm;
1046:   const PetscInt *support;
1047:   PetscInt        supportSize, s, off = 0, idx = 0;
1048:   PetscMPIInt     rank;
1049:   PetscErrorCode  ierr;

1052:   PetscObjectGetComm((PetscObject) dm, &comm);
1053:   MPI_Comm_rank(comm, &rank);
1054:   DMPlexGetSupportSize(dm, p, &supportSize);
1055:   DMPlexGetSupport(dm, p, &support);
1056:   if (candidates) {PetscSectionGetOffset(candidateSection, p, &off);}
1057:   for (s = 0; s < supportSize; ++s) {
1058:     const PetscInt  face = support[s];
1059:     const PetscInt *cone;
1060:     PetscSFNode     cpmin={-1,-1}, rp={-1,-1};
1061:     PetscInt        coneSize, c, f;
1062:     PetscBool       isShared = PETSC_FALSE;
1063:     PetscHashIJKey  key;

1065:     /* Only add point once */
1066:     if (debug) {PetscSynchronizedPrintf(comm, "[%d]    Support face %D\n", rank, face);}
1067:     key.i = p;
1068:     key.j = face;
1069:     PetscHMapIJGet(faceHash, key, &f);
1070:     if (f >= 0) continue;
1071:     DMPlexConeIsShared(dm, face, &isShared);
1072:     DMPlexGetConeMinimum(dm, face, &cpmin);
1073:     DMPlexMapToGlobalPoint(dm, p, &rp);
1074:     if (debug) {
1075:       PetscSynchronizedPrintf(comm, "[%d]      Face point %D is shared: %d\n", rank, face, (int) isShared);
1076:       PetscSynchronizedPrintf(comm, "[%d]      Global point (%D, %D) Min Cone Point (%D, %D)\n", rank, rp.rank, rp.index, cpmin.rank, cpmin.index);
1077:     }
1078:     if (isShared && (rp.rank == cpmin.rank && rp.index == cpmin.index)) {
1079:       PetscHMapIJSet(faceHash, key, p);
1080:       if (candidates) {
1081:         if (debug) {PetscSynchronizedPrintf(comm, "[%d]    Adding shared face %D at idx %D\n[%d]     ", rank, face, idx, rank);}
1082:         DMPlexGetConeSize(dm, face, &coneSize);
1083:         DMPlexGetCone(dm, face, &cone);
1084:         candidates[off+idx].rank    = -1;
1085:         candidates[off+idx++].index = coneSize-1;
1086:         candidates[off+idx].rank    = rank;
1087:         candidates[off+idx++].index = face;
1088:         for (c = 0; c < coneSize; ++c) {
1089:           const PetscInt cp = cone[c];

1091:           if (cp == p) continue;
1092:           DMPlexMapToGlobalPoint(dm, cp, &candidates[off+idx]);
1093:           if (debug) {PetscSynchronizedPrintf(comm, " (%D,%D)", candidates[off+idx].rank, candidates[off+idx].index);}
1094:           ++idx;
1095:         }
1096:         if (debug) {PetscSynchronizedPrintf(comm, "\n");}
1097:       } else {
1098:         /* Add cone size to section */
1099:         if (debug) {PetscSynchronizedPrintf(comm, "[%d]    Scheduling shared face %D\n", rank, face);}
1100:         DMPlexGetConeSize(dm, face, &coneSize);
1101:         PetscHMapIJSet(faceHash, key, p);
1102:         PetscSectionAddDof(candidateSection, p, coneSize+1);
1103:       }
1104:     }
1105:   }
1106:   return(0);
1107: }

1109: /*@
1110:   DMPlexInterpolatePointSF - Insert interpolated points in the overlap into the PointSF in parallel, following local interpolation

1112:   Collective on dm

1114:   Input Parameters:
1115: + dm      - The interpolated DM
1116: - pointSF - The initial SF without interpolated points

1118:   Output Parameter:
1119: . pointSF - The SF including interpolated points

1121:   Level: developer

1123:    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

1125: .seealso: DMPlexInterpolate(), DMPlexUninterpolate()
1126: @*/
1127: PetscErrorCode DMPlexInterpolatePointSF(DM dm, PetscSF pointSF)
1128: {
1129:   MPI_Comm           comm;
1130:   PetscHMapIJ        remoteHash;
1131:   PetscHMapI         claimshash;
1132:   PetscSection       candidateSection, candidateRemoteSection, claimSection;
1133:   PetscSFNode       *candidates, *candidatesRemote, *claims;
1134:   const PetscInt    *localPoints, *rootdegree;
1135:   const PetscSFNode *remotePoints;
1136:   PetscInt           ov, Nr, r, Nl, l;
1137:   PetscInt           candidatesSize, candidatesRemoteSize, claimsSize;
1138:   PetscBool          flg, debug = PETSC_FALSE;
1139:   PetscMPIInt        rank;
1140:   PetscErrorCode     ierr;

1145:   DMPlexIsDistributed(dm, &flg);
1146:   if (!flg) return(0);
1147:   /* Set initial SF so that lower level queries work */
1148:   DMSetPointSF(dm, pointSF);
1149:   PetscObjectGetComm((PetscObject) dm, &comm);
1150:   MPI_Comm_rank(comm, &rank);
1151:   DMPlexGetOverlap(dm, &ov);
1152:   if (ov) SETERRQ(comm, PETSC_ERR_SUP, "Interpolation of overlapped DMPlex not implemented yet");
1153:   PetscOptionsHasName(NULL, ((PetscObject) dm)->prefix, "-dmplex_interp_debug", &debug);
1154:   PetscObjectViewFromOptions((PetscObject) dm, NULL, "-dm_interp_pre_view");
1155:   PetscObjectViewFromOptions((PetscObject) pointSF, NULL, "-petscsf_interp_pre_view");
1156:   PetscLogEventBegin(DMPLEX_InterpolateSF,dm,0,0,0);
1157:   /* Step 0: Precalculations */
1158:   PetscSFGetGraph(pointSF, &Nr, &Nl, &localPoints, &remotePoints);
1159:   if (Nr < 0) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "This DMPlex is distributed but input PointSF has no graph set");
1160:   PetscHMapIJCreate(&remoteHash);
1161:   for (l = 0; l < Nl; ++l) {
1162:     PetscHashIJKey key;
1163:     key.i = remotePoints[l].index;
1164:     key.j = remotePoints[l].rank;
1165:     PetscHMapIJSet(remoteHash, key, l);
1166:   }
1167:   /*   Compute root degree to identify shared points */
1168:   PetscSFComputeDegreeBegin(pointSF, &rootdegree);
1169:   PetscSFComputeDegreeEnd(pointSF, &rootdegree);
1170:   IntArrayViewFromOptions(comm, "-interp_root_degree_view", "Root degree", "point", "degree", Nr, rootdegree);
1171:   /*
1172:   1) Loop over each leaf point $p$ at depth $d$ in the SF
1173:   \item Get set $F(p)$ of faces $f$ in the support of $p$ for which
1174:   \begin{itemize}
1175:     \item all cone points of $f$ are shared
1176:     \item $p$ is the cone point with smallest canonical number
1177:   \end{itemize}
1178:   \item Send $F(p)$ and the cone of each face to the active root point $r(p)$
1179:   \item At the root, if at least two faces with a given cone are present, including a local face, mark the face as shared \label{alg:rootStep} and choose the root face
1180:   \item Send the root face from the root back to all leaf process
1181:   \item Leaf processes add the shared face to the SF
1182:   */
1183:   /* Step 1: Construct section+SFNode array
1184:        The section has entries for all shared faces for which we have a leaf point in the cone
1185:        The array holds candidate shared faces, each face is refered to by the leaf point */
1186:   PetscSectionCreate(comm, &candidateSection);
1187:   PetscSectionSetChart(candidateSection, 0, Nr);
1188:   {
1189:     PetscHMapIJ faceHash;

1191:     PetscHMapIJCreate(&faceHash);
1192:     for (l = 0; l < Nl; ++l) {
1193:       const PetscInt p = localPoints[l];

1195:       if (debug) {PetscSynchronizedPrintf(comm, "[%d]  First pass leaf point %D\n", rank, p);}
1196:       DMPlexAddSharedFace_Private(dm, candidateSection, NULL, faceHash, p, debug);
1197:     }
1198:     PetscHMapIJClear(faceHash);
1199:     PetscSectionSetUp(candidateSection);
1200:     PetscSectionGetStorageSize(candidateSection, &candidatesSize);
1201:     PetscMalloc1(candidatesSize, &candidates);
1202:     for (l = 0; l < Nl; ++l) {
1203:       const PetscInt p = localPoints[l];

1205:       if (debug) {PetscSynchronizedPrintf(comm, "[%d]  Second pass leaf point %D\n", rank, p);}
1206:       DMPlexAddSharedFace_Private(dm, candidateSection, candidates, faceHash, p, debug);
1207:     }
1208:     PetscHMapIJDestroy(&faceHash);
1209:     if (debug) {PetscSynchronizedFlush(comm, NULL);}
1210:   }
1211:   PetscObjectSetName((PetscObject) candidateSection, "Candidate Section");
1212:   PetscObjectViewFromOptions((PetscObject) candidateSection, NULL, "-petscsection_interp_candidate_view");
1213:   SFNodeArrayViewFromOptions(comm, "-petscsection_interp_candidate_view", "Candidates", NULL, candidatesSize, candidates);
1214:   /* Step 2: Gather candidate section / array pair into the root partition via inverse(multi(pointSF)). */
1215:   /*   Note that this section is indexed by offsets into leaves, not by point number */
1216:   {
1217:     PetscSF   sfMulti, sfInverse, sfCandidates;
1218:     PetscInt *remoteOffsets;

1220:     PetscSFGetMultiSF(pointSF, &sfMulti);
1221:     PetscSFCreateInverseSF(sfMulti, &sfInverse);
1222:     PetscSectionCreate(comm, &candidateRemoteSection);
1223:     PetscSFDistributeSection(sfInverse, candidateSection, &remoteOffsets, candidateRemoteSection);
1224:     PetscSFCreateSectionSF(sfInverse, candidateSection, remoteOffsets, candidateRemoteSection, &sfCandidates);
1225:     PetscSectionGetStorageSize(candidateRemoteSection, &candidatesRemoteSize);
1226:     PetscMalloc1(candidatesRemoteSize, &candidatesRemote);
1227:     PetscSFBcastBegin(sfCandidates, MPIU_2INT, candidates, candidatesRemote);
1228:     PetscSFBcastEnd(sfCandidates, MPIU_2INT, candidates, candidatesRemote);
1229:     PetscSFDestroy(&sfInverse);
1230:     PetscSFDestroy(&sfCandidates);
1231:     PetscFree(remoteOffsets);

1233:     PetscObjectSetName((PetscObject) candidateRemoteSection, "Remote Candidate Section");
1234:     PetscObjectViewFromOptions((PetscObject) candidateRemoteSection, NULL, "-petscsection_interp_candidate_remote_view");
1235:     SFNodeArrayViewFromOptions(comm, "-petscsection_interp_candidate_remote_view", "Remote Candidates", NULL, candidatesRemoteSize, candidatesRemote);
1236:   }
1237:   /* Step 3: At the root, if at least two faces with a given cone are present, including a local face, mark the face as shared and choose the root face */
1238:   {
1239:     PetscHashIJKLRemote faceTable;
1240:     PetscInt            idx, idx2;

1242:     PetscHashIJKLRemoteCreate(&faceTable);
1243:     /* There is a section point for every leaf attached to a given root point */
1244:     for (r = 0, idx = 0, idx2 = 0; r < Nr; ++r) {
1245:       PetscInt deg;

1247:       for (deg = 0; deg < rootdegree[r]; ++deg, ++idx) {
1248:         PetscInt offset, dof, d;

1250:         PetscSectionGetDof(candidateRemoteSection, idx, &dof);
1251:         PetscSectionGetOffset(candidateRemoteSection, idx, &offset);
1252:         /* dof may include many faces from the remote process */
1253:         for (d = 0; d < dof; ++d) {
1254:           const PetscInt         hidx  = offset+d;
1255:           const PetscInt         Np    = candidatesRemote[hidx].index+1;
1256:           const PetscSFNode      rface = candidatesRemote[hidx+1];
1257:           const PetscSFNode     *fcone = &candidatesRemote[hidx+2];
1258:           PetscSFNode            fcp0;
1259:           const PetscSFNode      pmax  = {PETSC_MAX_INT, PETSC_MAX_INT};
1260:           const PetscInt        *join  = NULL;
1261:           PetscHashIJKLRemoteKey key;
1262:           PetscHashIter          iter;
1263:           PetscBool              missing;
1264:           PetscInt               points[1024], p, joinSize;

1266:           if (debug) {PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]  Checking face (%D, %D) at (%D, %D, %D) with cone size %D\n", rank, rface.index, rface.rank, r, idx, d, Np);}
1267:           if (Np > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot handle faces with %D cone points", Np);
1268:           fcp0.rank  = rank;
1269:           fcp0.index = r;
1270:           d += Np;
1271:           /* Put remote face in hash table */
1272:           key.i = fcp0;
1273:           key.j = fcone[0];
1274:           key.k = Np > 2 ? fcone[1] : pmax;
1275:           key.l = Np > 3 ? fcone[2] : pmax;
1276:           PetscSortSFNode(Np, (PetscSFNode *) &key);
1277:           PetscHashIJKLRemotePut(faceTable, key, &iter, &missing);
1278:           if (missing) {
1279:             if (debug) {PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]  Setting remote face (%D, %D)\n", rank, rface.index, rface.rank);}
1280:             PetscHashIJKLRemoteIterSet(faceTable, iter, rface);
1281:           } else {
1282:             PetscSFNode oface;

1284:             PetscHashIJKLRemoteIterGet(faceTable, iter, &oface);
1285:             if ((rface.rank < oface.rank) || (rface.rank == oface.rank && rface.index < oface.index)) {
1286:               if (debug) {PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]  Replacing with remote face (%D, %D)\n", rank, rface.index, rface.rank);}
1287:               PetscHashIJKLRemoteIterSet(faceTable, iter, rface);
1288:             }
1289:           }
1290:           /* Check for local face */
1291:           points[0] = r;
1292:           for (p = 1; p < Np; ++p) {
1293:             DMPlexMapToLocalPoint(dm, remoteHash, fcone[p-1], &points[p]);
1294:             if (ierr) break; /* We got a point not in our overlap */
1295:             if (debug) {PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]  Checking local candidate %D\n", rank, points[p]);}
1296:           }
1297:           if (ierr) continue;
1298:           DMPlexGetJoin(dm, Np, points, &joinSize, &join);
1299:           if (joinSize == 1) {
1300:             PetscSFNode lface;
1301:             PetscSFNode oface;

1303:             /* Always replace with local face */
1304:             lface.rank  = rank;
1305:             lface.index = join[0];
1306:             PetscHashIJKLRemoteIterGet(faceTable, iter, &oface);
1307:             if (debug) {PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]  Replacing (%D, %D) with local face (%D, %D)\n", rank, oface.index, oface.rank, lface.index, lface.rank);}
1308:             PetscHashIJKLRemoteIterSet(faceTable, iter, lface);
1309:           }
1310:           DMPlexRestoreJoin(dm, Np, points, &joinSize, &join);
1311:         }
1312:       }
1313:       /* Put back faces for this root */
1314:       for (deg = 0; deg < rootdegree[r]; ++deg, ++idx2) {
1315:         PetscInt offset, dof, d;

1317:         PetscSectionGetDof(candidateRemoteSection, idx2, &dof);
1318:         PetscSectionGetOffset(candidateRemoteSection, idx2, &offset);
1319:         /* dof may include many faces from the remote process */
1320:         for (d = 0; d < dof; ++d) {
1321:           const PetscInt         hidx  = offset+d;
1322:           const PetscInt         Np    = candidatesRemote[hidx].index+1;
1323:           const PetscSFNode     *fcone = &candidatesRemote[hidx+2];
1324:           PetscSFNode            fcp0;
1325:           const PetscSFNode      pmax  = {PETSC_MAX_INT, PETSC_MAX_INT};
1326:           PetscHashIJKLRemoteKey key;
1327:           PetscHashIter          iter;
1328:           PetscBool              missing;

1330:           if (debug) {PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]  Entering face at (%D, %D)\n", rank, r, idx);}
1331:           if (Np > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot handle faces with %D cone points", Np);
1332:           fcp0.rank  = rank;
1333:           fcp0.index = r;
1334:           d += Np;
1335:           /* Find remote face in hash table */
1336:           key.i = fcp0;
1337:           key.j = fcone[0];
1338:           key.k = Np > 2 ? fcone[1] : pmax;
1339:           key.l = Np > 3 ? fcone[2] : pmax;
1340:           PetscSortSFNode(Np, (PetscSFNode *) &key);
1341:           if (debug) {PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d]    key (%D, %D) (%D, %D) (%D, %D) (%D, %D)\n", rank, key.i.rank, key.i.index, key.j.rank, key.j.index, key.k.rank, key.k.index, key.l.rank, key.l.index);}
1342:           PetscHashIJKLRemotePut(faceTable, key, &iter, &missing);
1343:           if (missing) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Root %D Idx %D ought to have an assoicated face", r, idx2);
1344:           else        {PetscHashIJKLRemoteIterGet(faceTable, iter, &candidatesRemote[hidx]);}
1345:         }
1346:       }
1347:     }
1348:     if (debug) {PetscSynchronizedFlush(PetscObjectComm((PetscObject) dm), NULL);}
1349:     PetscHashIJKLRemoteDestroy(&faceTable);
1350:   }
1351:   /* Step 4: Push back owned faces */
1352:   {
1353:     PetscSF      sfMulti, sfClaims, sfPointNew;
1354:     PetscSFNode *remotePointsNew;
1355:     PetscInt    *remoteOffsets, *localPointsNew;
1356:     PetscInt     pStart, pEnd, r, NlNew, p;

1358:     /* 4) Push claims back to receiver via the MultiSF and derive new pointSF mapping on receiver */
1359:     PetscSFGetMultiSF(pointSF, &sfMulti);
1360:     PetscSectionCreate(comm, &claimSection);
1361:     PetscSFDistributeSection(sfMulti, candidateRemoteSection, &remoteOffsets, claimSection);
1362:     PetscSFCreateSectionSF(sfMulti, candidateRemoteSection, remoteOffsets, claimSection, &sfClaims);
1363:     PetscSectionGetStorageSize(claimSection, &claimsSize);
1364:     PetscMalloc1(claimsSize, &claims);
1365:     for (p = 0; p < claimsSize; ++p) claims[p].rank = -1;
1366:     PetscSFBcastBegin(sfClaims, MPIU_2INT, candidatesRemote, claims);
1367:     PetscSFBcastEnd(sfClaims, MPIU_2INT, candidatesRemote, claims);
1368:     PetscSFDestroy(&sfClaims);
1369:     PetscFree(remoteOffsets);
1370:     PetscObjectSetName((PetscObject) claimSection, "Claim Section");
1371:     PetscObjectViewFromOptions((PetscObject) claimSection, NULL, "-petscsection_interp_claim_view");
1372:     SFNodeArrayViewFromOptions(comm, "-petscsection_interp_claim_view", "Claims", NULL, claimsSize, claims);
1373:     /* Step 5) Walk the original section of local supports and add an SF entry for each updated item */
1374:     /* TODO I should not have to do a join here since I already put the face and its cone in the candidate section */
1375:     PetscHMapICreate(&claimshash);
1376:     for (r = 0; r < Nr; ++r) {
1377:       PetscInt dof, off, d;

1379:       if (debug) {PetscSynchronizedPrintf(comm, "[%d]  Checking root for claims %D\n", rank, r);}
1380:       PetscSectionGetDof(candidateSection, r, &dof);
1381:       PetscSectionGetOffset(candidateSection, r, &off);
1382:       for (d = 0; d < dof;) {
1383:         if (claims[off+d].rank >= 0) {
1384:           const PetscInt  faceInd = off+d;
1385:           const PetscInt  Np      = candidates[off+d].index;
1386:           const PetscInt *join    = NULL;
1387:           PetscInt        joinSize, points[1024], c;

1389:           if (debug) {PetscSynchronizedPrintf(comm, "[%d]    Found claim for remote point (%D, %D)\n", rank, claims[faceInd].rank, claims[faceInd].index);}
1390:           points[0] = r;
1391:           if (debug) {PetscSynchronizedPrintf(comm, "[%d]      point %D\n", rank, points[0]);}
1392:           for (c = 0, d += 2; c < Np; ++c, ++d) {
1393:             DMPlexMapToLocalPoint(dm, remoteHash, candidates[off+d], &points[c+1]);
1394:             if (debug) {PetscSynchronizedPrintf(comm, "[%d]      point %D\n", rank, points[c+1]);}
1395:           }
1396:           DMPlexGetJoin(dm, Np+1, points, &joinSize, &join);
1397:           if (joinSize == 1) {
1398:             if (claims[faceInd].rank == rank) {
1399:               if (debug) {PetscSynchronizedPrintf(comm, "[%d]    Ignoring local face %D for non-remote partner\n", rank, join[0]);}
1400:             } else {
1401:               if (debug) {PetscSynchronizedPrintf(comm, "[%d]    Found local face %D\n", rank, join[0]);}
1402:               PetscHMapISet(claimshash, join[0], faceInd);
1403:             }
1404:           } else {
1405:             if (debug) {PetscSynchronizedPrintf(comm, "[%d]    Failed to find face\n", rank);}
1406:           }
1407:           DMPlexRestoreJoin(dm, Np+1, points, &joinSize, &join);
1408:         } else {
1409:           if (debug) {PetscSynchronizedPrintf(comm, "[%d]    No claim for point %D\n", rank, r);}
1410:           d += claims[off+d].index+1;
1411:         }
1412:       }
1413:     }
1414:     if (debug) {PetscSynchronizedFlush(comm, NULL);}
1415:     /* Step 6) Create new pointSF from hashed claims */
1416:     PetscHMapIGetSize(claimshash, &NlNew);
1417:     DMPlexGetChart(dm, &pStart, &pEnd);
1418:     PetscMalloc1(Nl + NlNew, &localPointsNew);
1419:     PetscMalloc1(Nl + NlNew, &remotePointsNew);
1420:     for (l = 0; l < Nl; ++l) {
1421:       localPointsNew[l] = localPoints[l];
1422:       remotePointsNew[l].index = remotePoints[l].index;
1423:       remotePointsNew[l].rank  = remotePoints[l].rank;
1424:     }
1425:     p = Nl;
1426:     PetscHMapIGetKeys(claimshash, &p, localPointsNew);
1427:     /* We sort new points, and assume they are numbered after all existing points */
1428:     PetscSortInt(NlNew, &localPointsNew[Nl]);
1429:     for (p = Nl; p < Nl + NlNew; ++p) {
1430:       PetscInt off;
1431:       PetscHMapIGet(claimshash, localPointsNew[p], &off);
1432:       if (claims[off].rank < 0 || claims[off].index < 0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid claim for local point %D, (%D, %D)", localPointsNew[p], claims[off].rank, claims[off].index);
1433:       remotePointsNew[p] = claims[off];
1434:     }
1435:     PetscSFCreate(comm, &sfPointNew);
1436:     PetscSFSetGraph(sfPointNew, pEnd-pStart, Nl+NlNew, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);
1437:     PetscSFSetUp(sfPointNew);
1438:     DMSetPointSF(dm, sfPointNew);
1439:     PetscObjectViewFromOptions((PetscObject) sfPointNew, NULL, "-petscsf_interp_view");
1440:     PetscSFDestroy(&sfPointNew);
1441:     PetscHMapIDestroy(&claimshash);
1442:   }
1443:   PetscHMapIJDestroy(&remoteHash);
1444:   PetscSectionDestroy(&candidateSection);
1445:   PetscSectionDestroy(&candidateRemoteSection);
1446:   PetscSectionDestroy(&claimSection);
1447:   PetscFree(candidates);
1448:   PetscFree(candidatesRemote);
1449:   PetscFree(claims);
1450:   PetscLogEventEnd(DMPLEX_InterpolateSF,dm,0,0,0);
1451:   return(0);
1452: }

1454: /*@
1455:   DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc.

1457:   Collective on dm

1459:   Input Parameters:
1460: + dm - The DMPlex object with only cells and vertices
1461: - dmInt - The interpolated DM

1463:   Output Parameter:
1464: . dmInt - The complete DMPlex object

1466:   Level: intermediate

1468:   Notes:
1469:     It does not copy over the coordinates.

1471:   Developer Notes:
1472:     It sets plex->interpolated = DMPLEX_INTERPOLATED_FULL.

1474: .seealso: DMPlexUninterpolate(), DMPlexCreateFromCellList(), DMPlexCopyCoordinates()
1475: @*/
1476: PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt)
1477: {
1478:   DMPlexInterpolatedFlag interpolated;
1479:   DM             idm, odm = dm;
1480:   PetscSF        sfPoint;
1481:   PetscInt       depth, dim, d;
1482:   const char    *name;
1483:   PetscBool      flg=PETSC_TRUE;

1489:   PetscLogEventBegin(DMPLEX_Interpolate,dm,0,0,0);
1490:   DMPlexGetDepth(dm, &depth);
1491:   DMGetDimension(dm, &dim);
1492:   DMPlexIsInterpolated(dm, &interpolated);
1493:   if (interpolated == DMPLEX_INTERPOLATED_PARTIAL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Not for partially interpolated meshes");
1494:   if (interpolated == DMPLEX_INTERPOLATED_FULL) {
1495:     PetscObjectReference((PetscObject) dm);
1496:     idm  = dm;
1497:   } else {
1498:     for (d = 1; d < dim; ++d) {
1499:       /* Create interpolated mesh */
1500:       DMCreate(PetscObjectComm((PetscObject)dm), &idm);
1501:       DMSetType(idm, DMPLEX);
1502:       DMSetDimension(idm, dim);
1503:       if (depth > 0) {
1504:         DMPlexInterpolateFaces_Internal(odm, 1, idm);
1505:         DMGetPointSF(odm, &sfPoint);
1506:         {
1507:           /* TODO: We need to systematically fix cases of distributed Plexes with no graph set */
1508:           PetscInt nroots;
1509:           PetscSFGetGraph(sfPoint, &nroots, NULL, NULL, NULL);
1510:           if (nroots >= 0) {DMPlexInterpolatePointSF(idm, sfPoint);}
1511:         }
1512:       }
1513:       if (odm != dm) {DMDestroy(&odm);}
1514:       odm = idm;
1515:     }
1516:     PetscObjectGetName((PetscObject) dm,  &name);
1517:     PetscObjectSetName((PetscObject) idm,  name);
1518:     DMPlexCopyCoordinates(dm, idm);
1519:     DMCopyLabels(dm, idm, PETSC_COPY_VALUES, PETSC_FALSE);
1520:     PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_interpolate_orient_interfaces", &flg, NULL);
1521:     if (flg) {DMPlexOrientInterface_Internal(idm);}
1522:   }
1523:   {
1524:     PetscBool            isper;
1525:     const PetscReal      *maxCell, *L;
1526:     const DMBoundaryType *bd;

1528:     DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);
1529:     DMSetPeriodicity(idm,isper,maxCell,L,bd);
1530:   }
1531:   /* This function makes the mesh fully interpolated on all ranks */
1532:   {
1533:     DM_Plex *plex = (DM_Plex *) idm->data;
1534:     plex->interpolated = plex->interpolatedCollective = DMPLEX_INTERPOLATED_FULL;
1535:   }
1536:   *dmInt = idm;
1537:   PetscLogEventEnd(DMPLEX_Interpolate,dm,0,0,0);
1538:   return(0);
1539: }

1541: /*@
1542:   DMPlexCopyCoordinates - Copy coordinates from one mesh to another with the same vertices

1544:   Collective on dmA

1546:   Input Parameter:
1547: . dmA - The DMPlex object with initial coordinates

1549:   Output Parameter:
1550: . dmB - The DMPlex object with copied coordinates

1552:   Level: intermediate

1554:   Note: This is typically used when adding pieces other than vertices to a mesh

1556: .seealso: DMCopyLabels(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMGetCoordinateSection()
1557: @*/
1558: PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB)
1559: {
1560:   Vec            coordinatesA, coordinatesB;
1561:   VecType        vtype;
1562:   PetscSection   coordSectionA, coordSectionB;
1563:   PetscScalar   *coordsA, *coordsB;
1564:   PetscInt       spaceDim, Nf, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d;
1565:   PetscInt       cStartA, cEndA, cStartB, cEndB, cS, cE;
1566:   PetscBool      lc = PETSC_FALSE;

1572:   if (dmA == dmB) return(0);
1573:   DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);
1574:   DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);
1575:   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);
1576:   DMPlexGetHeightStratum(dmA, 0, &cStartA, &cEndA);
1577:   DMPlexGetHeightStratum(dmB, 0, &cStartB, &cEndB);
1578:   DMGetCoordinateSection(dmA, &coordSectionA);
1579:   DMGetCoordinateSection(dmB, &coordSectionB);
1580:   if (coordSectionA == coordSectionB) return(0);
1581:   PetscSectionGetNumFields(coordSectionA, &Nf);
1582:   if (!Nf) return(0);
1583:   if (Nf > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of coordinate fields must be 1, not %D", Nf);
1584:   if (!coordSectionB) {
1585:     PetscInt dim;

1587:     PetscSectionCreate(PetscObjectComm((PetscObject) coordSectionA), &coordSectionB);
1588:     DMGetCoordinateDim(dmA, &dim);
1589:     DMSetCoordinateSection(dmB, dim, coordSectionB);
1590:     PetscObjectDereference((PetscObject) coordSectionB);
1591:   }
1592:   PetscSectionSetNumFields(coordSectionB, 1);
1593:   PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);
1594:   PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);
1595:   PetscSectionGetChart(coordSectionA, &cS, &cE);
1596:   if (cStartA <= cS && cS < cEndA) { /* localized coordinates */
1597:     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);
1598:     cS = cS - cStartA + cStartB;
1599:     cE = vEndB;
1600:     lc = PETSC_TRUE;
1601:   } else {
1602:     cS = vStartB;
1603:     cE = vEndB;
1604:   }
1605:   PetscSectionSetChart(coordSectionB, cS, cE);
1606:   for (v = vStartB; v < vEndB; ++v) {
1607:     PetscSectionSetDof(coordSectionB, v, spaceDim);
1608:     PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);
1609:   }
1610:   if (lc) { /* localized coordinates */
1611:     PetscInt c;

1613:     for (c = cS-cStartB; c < cEndB-cStartB; c++) {
1614:       PetscInt dof;

1616:       PetscSectionGetDof(coordSectionA, c + cStartA, &dof);
1617:       PetscSectionSetDof(coordSectionB, c + cStartB, dof);
1618:       PetscSectionSetFieldDof(coordSectionB, c + cStartB, 0, dof);
1619:     }
1620:   }
1621:   PetscSectionSetUp(coordSectionB);
1622:   PetscSectionGetStorageSize(coordSectionB, &coordSizeB);
1623:   DMGetCoordinatesLocal(dmA, &coordinatesA);
1624:   VecCreate(PETSC_COMM_SELF, &coordinatesB);
1625:   PetscObjectSetName((PetscObject) coordinatesB, "coordinates");
1626:   VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);
1627:   VecGetBlockSize(coordinatesA, &d);
1628:   VecSetBlockSize(coordinatesB, d);
1629:   VecGetType(coordinatesA, &vtype);
1630:   VecSetType(coordinatesB, vtype);
1631:   VecGetArray(coordinatesA, &coordsA);
1632:   VecGetArray(coordinatesB, &coordsB);
1633:   for (v = 0; v < vEndB-vStartB; ++v) {
1634:     PetscInt offA, offB;

1636:     PetscSectionGetOffset(coordSectionA, v + vStartA, &offA);
1637:     PetscSectionGetOffset(coordSectionB, v + vStartB, &offB);
1638:     for (d = 0; d < spaceDim; ++d) {
1639:       coordsB[offB+d] = coordsA[offA+d];
1640:     }
1641:   }
1642:   if (lc) { /* localized coordinates */
1643:     PetscInt c;

1645:     for (c = cS-cStartB; c < cEndB-cStartB; c++) {
1646:       PetscInt dof, offA, offB;

1648:       PetscSectionGetOffset(coordSectionA, c + cStartA, &offA);
1649:       PetscSectionGetOffset(coordSectionB, c + cStartB, &offB);
1650:       PetscSectionGetDof(coordSectionA, c + cStartA, &dof);
1651:       PetscArraycpy(coordsB + offB,coordsA + offA,dof);
1652:     }
1653:   }
1654:   VecRestoreArray(coordinatesA, &coordsA);
1655:   VecRestoreArray(coordinatesB, &coordsB);
1656:   DMSetCoordinatesLocal(dmB, coordinatesB);
1657:   VecDestroy(&coordinatesB);
1658:   return(0);
1659: }

1661: /*@
1662:   DMPlexUninterpolate - Take in a mesh with all intermediate faces, edges, etc. and return a cell-vertex mesh

1664:   Collective on dm

1666:   Input Parameter:
1667: . dm - The complete DMPlex object

1669:   Output Parameter:
1670: . dmUnint - The DMPlex object with only cells and vertices

1672:   Level: intermediate

1674:   Notes:
1675:     It does not copy over the coordinates.

1677:   Developer Notes:
1678:     It sets plex->interpolated = DMPLEX_INTERPOLATED_NONE.

1680: .seealso: DMPlexInterpolate(), DMPlexCreateFromCellList(), DMPlexCopyCoordinates()
1681: @*/
1682: PetscErrorCode DMPlexUninterpolate(DM dm, DM *dmUnint)
1683: {
1684:   DMPlexInterpolatedFlag interpolated;
1685:   DM             udm;
1686:   PetscInt       dim, vStart, vEnd, cStart, cEnd, cMax, c, maxConeSize = 0, *cone;

1692:   DMGetDimension(dm, &dim);
1693:   DMPlexIsInterpolated(dm, &interpolated);
1694:   if (interpolated == DMPLEX_INTERPOLATED_PARTIAL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Not for partially interpolated meshes");
1695:   if (interpolated == DMPLEX_INTERPOLATED_NONE || dim <= 1) {
1696:     /* in case dim <= 1 just keep the DMPLEX_INTERPOLATED_FULL flag */
1697:     PetscObjectReference((PetscObject) dm);
1698:     *dmUnint = dm;
1699:     return(0);
1700:   }
1701:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
1702:   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
1703:   DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);
1704:   DMCreate(PetscObjectComm((PetscObject) dm), &udm);
1705:   DMSetType(udm, DMPLEX);
1706:   DMSetDimension(udm, dim);
1707:   DMPlexSetChart(udm, cStart, vEnd);
1708:   for (c = cStart; c < cEnd; ++c) {
1709:     PetscInt *closure = NULL, closureSize, cl, coneSize = 0;

1711:     DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
1712:     for (cl = 0; cl < closureSize*2; cl += 2) {
1713:       const PetscInt p = closure[cl];

1715:       if ((p >= vStart) && (p < vEnd)) ++coneSize;
1716:     }
1717:     DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
1718:     DMPlexSetConeSize(udm, c, coneSize);
1719:     maxConeSize = PetscMax(maxConeSize, coneSize);
1720:   }
1721:   DMSetUp(udm);
1722:   PetscMalloc1(maxConeSize, &cone);
1723:   for (c = cStart; c < cEnd; ++c) {
1724:     PetscInt *closure = NULL, closureSize, cl, coneSize = 0;

1726:     DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
1727:     for (cl = 0; cl < closureSize*2; cl += 2) {
1728:       const PetscInt p = closure[cl];

1730:       if ((p >= vStart) && (p < vEnd)) cone[coneSize++] = p;
1731:     }
1732:     DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
1733:     DMPlexSetCone(udm, c, cone);
1734:   }
1735:   PetscFree(cone);
1736:   DMPlexSetHybridBounds(udm, cMax, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);
1737:   DMPlexSymmetrize(udm);
1738:   DMPlexStratify(udm);
1739:   /* Reduce SF */
1740:   {
1741:     PetscSF            sfPoint, sfPointUn;
1742:     const PetscSFNode *remotePoints;
1743:     const PetscInt    *localPoints;
1744:     PetscSFNode       *remotePointsUn;
1745:     PetscInt          *localPointsUn;
1746:     PetscInt           vEnd, numRoots, numLeaves, l;
1747:     PetscInt           numLeavesUn = 0, n = 0;
1748:     PetscErrorCode     ierr;

1750:     /* Get original SF information */
1751:     DMGetPointSF(dm, &sfPoint);
1752:     DMGetPointSF(udm, &sfPointUn);
1753:     DMPlexGetDepthStratum(dm, 0, NULL, &vEnd);
1754:     PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);
1755:     /* Allocate space for cells and vertices */
1756:     for (l = 0; l < numLeaves; ++l) if (localPoints[l] < vEnd) numLeavesUn++;
1757:     /* Fill in leaves */
1758:     if (vEnd >= 0) {
1759:       PetscMalloc1(numLeavesUn, &remotePointsUn);
1760:       PetscMalloc1(numLeavesUn, &localPointsUn);
1761:       for (l = 0; l < numLeaves; l++) {
1762:         if (localPoints[l] < vEnd) {
1763:           localPointsUn[n]        = localPoints[l];
1764:           remotePointsUn[n].rank  = remotePoints[l].rank;
1765:           remotePointsUn[n].index = remotePoints[l].index;
1766:           ++n;
1767:         }
1768:       }
1769:       if (n != numLeavesUn) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent number of leaves %d != %d", n, numLeavesUn);
1770:       PetscSFSetGraph(sfPointUn, vEnd, numLeavesUn, localPointsUn, PETSC_OWN_POINTER, remotePointsUn, PETSC_OWN_POINTER);
1771:     }
1772:   }
1773:   {
1774:     PetscBool            isper;
1775:     const PetscReal      *maxCell, *L;
1776:     const DMBoundaryType *bd;

1778:     DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);
1779:     DMSetPeriodicity(udm,isper,maxCell,L,bd);
1780:   }
1781:   /* This function makes the mesh fully uninterpolated on all ranks */
1782:   {
1783:     DM_Plex *plex = (DM_Plex *) udm->data;
1784:     plex->interpolated = plex->interpolatedCollective = DMPLEX_INTERPOLATED_NONE;
1785:   }
1786:   *dmUnint = udm;
1787:   return(0);
1788: }

1790: static PetscErrorCode DMPlexIsInterpolated_Internal(DM dm, DMPlexInterpolatedFlag *interpolated)
1791: {
1792:   PetscInt       coneSize, depth, dim, h, p, pStart, pEnd;
1793:   MPI_Comm       comm;

1797:   PetscObjectGetComm((PetscObject)dm, &comm);
1798:   DMPlexGetDepth(dm, &depth);
1799:   DMGetDimension(dm, &dim);

1801:   if (depth == dim) {
1802:     *interpolated = DMPLEX_INTERPOLATED_FULL;
1803:     if (!dim) goto finish;

1805:     /* Check points at height = dim are vertices (have no cones) */
1806:     DMPlexGetHeightStratum(dm, dim, &pStart, &pEnd);
1807:     for (p=pStart; p<pEnd; p++) {
1808:       DMPlexGetConeSize(dm, p, &coneSize);
1809:       if (coneSize) {
1810:         *interpolated = DMPLEX_INTERPOLATED_PARTIAL;
1811:         goto finish;
1812:       }
1813:     }

1815:     /* Check points at height < dim have cones */
1816:     for (h=0; h<dim; h++) {
1817:       DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);
1818:       for (p=pStart; p<pEnd; p++) {
1819:         DMPlexGetConeSize(dm, p, &coneSize);
1820:         if (!coneSize) {
1821:           *interpolated = DMPLEX_INTERPOLATED_PARTIAL;
1822:           goto finish;
1823:         }
1824:       }
1825:     }
1826:   } else if (depth == 1) {
1827:     *interpolated = DMPLEX_INTERPOLATED_NONE;
1828:   } else {
1829:     *interpolated = DMPLEX_INTERPOLATED_PARTIAL;
1830:   }
1831: finish:
1832:   return(0);
1833: }

1835: /*@
1836:   DMPlexIsInterpolated - Find out to what extent the DMPlex is topologically interpolated.

1838:   Not Collective

1840:   Input Parameter:
1841: . dm      - The DM object

1843:   Output Parameter:
1844: . interpolated - Flag whether the DM is interpolated

1846:   Level: intermediate

1848:   Notes:
1849:   Unlike DMPlexIsInterpolatedCollective(), this is NOT collective
1850:   so the results can be different on different ranks in special cases.
1851:   However, DMPlexInterpolate() guarantees the result is the same on all.

1853:   Unlike DMPlexIsInterpolatedCollective(), this cannot return DMPLEX_INTERPOLATED_MIXED.

1855:   Developer Notes:
1856:   Initially, plex->interpolated = DMPLEX_INTERPOLATED_INVALID.

1858:   If plex->interpolated == DMPLEX_INTERPOLATED_INVALID, DMPlexIsInterpolated_Internal() is called.
1859:   It checks the actual topology and sets plex->interpolated on each rank separately to one of
1860:   DMPLEX_INTERPOLATED_NONE, DMPLEX_INTERPOLATED_PARTIAL or DMPLEX_INTERPOLATED_FULL.

1862:   If plex->interpolated != DMPLEX_INTERPOLATED_INVALID, this function just returns plex->interpolated.

1864:   DMPlexInterpolate() sets plex->interpolated = DMPLEX_INTERPOLATED_FULL,
1865:   and DMPlexUninterpolate() sets plex->interpolated = DMPLEX_INTERPOLATED_NONE.

1867: .seealso: DMPlexInterpolate(), DMPlexIsInterpolatedCollective()
1868: @*/
1869: PetscErrorCode DMPlexIsInterpolated(DM dm, DMPlexInterpolatedFlag *interpolated)
1870: {
1871:   DM_Plex        *plex = (DM_Plex *) dm->data;
1872:   PetscErrorCode  ierr;

1877:   if (plex->interpolated < 0) {
1878:     DMPlexIsInterpolated_Internal(dm, &plex->interpolated);
1879:   } else {
1880: #if defined (PETSC_USE_DEBUG)
1881:     DMPlexInterpolatedFlag flg;

1883:     DMPlexIsInterpolated_Internal(dm, &flg);
1884:     if (flg != plex->interpolated) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Stashed DMPlexInterpolatedFlag %s is inconsistent with current %s", DMPlexInterpolatedFlags[plex->interpolated], DMPlexInterpolatedFlags[flg]);
1885: #endif
1886:   }
1887:   *interpolated = plex->interpolated;
1888:   return(0);
1889: }

1891: /*@
1892:   DMPlexIsInterpolatedCollective - Find out to what extent the DMPlex is topologically interpolated (in collective manner).

1894:   Collective

1896:   Input Parameter:
1897: . dm      - The DM object

1899:   Output Parameter:
1900: . interpolated - Flag whether the DM is interpolated

1902:   Level: intermediate

1904:   Notes:
1905:   Unlike DMPlexIsInterpolated(), this is collective so the results are guaranteed to be the same on all ranks.

1907:   This function will return DMPLEX_INTERPOLATED_MIXED if the results of DMPlexIsInterpolated() are different on different ranks.

1909:   Developer Notes:
1910:   Initially, plex->interpolatedCollective = DMPLEX_INTERPOLATED_INVALID.

1912:   If plex->interpolatedCollective == DMPLEX_INTERPOLATED_INVALID, this function calls DMPlexIsInterpolated() which sets plex->interpolated.
1913:   MPI_Allreduce() is then called and collectively consistent flag plex->interpolatedCollective is set and returned;
1914:   if plex->interpolated varies on different ranks, plex->interpolatedCollective = DMPLEX_INTERPOLATED_MIXED,
1915:   otherwise sets plex->interpolatedCollective = plex->interpolated.

1917:   If plex->interpolatedCollective != DMPLEX_INTERPOLATED_INVALID, this function just returns plex->interpolatedCollective.

1919: .seealso: DMPlexInterpolate(), DMPlexIsInterpolated()
1920: @*/
1921: PetscErrorCode DMPlexIsInterpolatedCollective(DM dm, DMPlexInterpolatedFlag *interpolated)
1922: {
1923:   DM_Plex        *plex = (DM_Plex *) dm->data;
1924:   PetscBool       debug=PETSC_FALSE;
1925:   PetscErrorCode  ierr;

1930:   PetscOptionsGetBool(((PetscObject) dm)->options, ((PetscObject) dm)->prefix, "-dm_plex_is_interpolated_collective_debug", &debug, NULL);
1931:   if (plex->interpolatedCollective < 0) {
1932:     DMPlexInterpolatedFlag  min, max;
1933:     MPI_Comm                comm;

1935:     PetscObjectGetComm((PetscObject)dm, &comm);
1936:     DMPlexIsInterpolated(dm, &plex->interpolatedCollective);
1937:     MPI_Allreduce(&plex->interpolatedCollective, &min, 1, MPIU_ENUM, MPI_MIN, comm);
1938:     MPI_Allreduce(&plex->interpolatedCollective, &max, 1, MPIU_ENUM, MPI_MAX, comm);
1939:     if (min != max) plex->interpolatedCollective = DMPLEX_INTERPOLATED_MIXED;
1940:     if (debug) {
1941:       PetscMPIInt rank;

1943:       MPI_Comm_rank(comm, &rank);
1944:       PetscSynchronizedPrintf(comm, "[%d] interpolated=%s interpolatedCollective=%s\n", rank, DMPlexInterpolatedFlags[plex->interpolated], DMPlexInterpolatedFlags[plex->interpolatedCollective]);
1945:       PetscSynchronizedFlush(comm, PETSC_STDOUT);
1946:     }
1947:   }
1948:   *interpolated = plex->interpolatedCollective;
1949:   return(0);
1950: }