Actual source code: plexorient.c

petsc-3.13.6 2020-09-29
Report Typos and Errors
  1:  #include <petsc/private/dmpleximpl.h>
  2:  #include <petscsf.h>

  4: /*@
  5:   DMPlexCompareOrientations - Compare the cone of the given DAG point (cell) with the given reference cone (with the same cone points modulo order), and return relative orientation.

  7:   Not Collective

  9:   Input Parameters:
 10: + dm              - The DM (DMPLEX)
 11: . p               - The DAG point whose cone is compared
 12: . masterConeSize  - Number of the reference cone points passed (at least 2 and at most size of the cone of p)
 13: - masterCone      - The reference cone points

 15:   Output Parameters:
 16: + start           - The new starting point within the cone of p to make it conforming with the reference cone
 17: - reverse         - The flag whether the order of the cone points should be reversed

 19:   Level: advanced

 21: .seealso: DMPlexOrient(), DMPlexOrientCell()
 22: @*/
 23: PetscErrorCode DMPlexCompareOrientations(DM dm, PetscInt p, PetscInt masterConeSize, const PetscInt masterCone[], PetscInt *start, PetscBool *reverse)
 24: {
 25:   PetscInt        coneSize;
 26:   const PetscInt *cone;
 27:   PetscInt        i, start_;
 28:   PetscBool       reverse_;
 29:   PetscErrorCode  ierr;

 33:   DMPlexGetConeSize(dm, p, &coneSize);
 34:   if (coneSize < 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %D has no cone", p);
 35:   DMPlexGetCone(dm, p, &cone);
 36:   if (masterConeSize < 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %D: masterConeSize must be at least 2", p);
 37:   if (masterConeSize > coneSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %D: masterConeSize must be at most coneSize", p);
 38:   start_ = 0;
 39:   for (i=0; i<coneSize; i++) {
 40:     if (cone[i] == masterCone[0]) {
 41:       start_ = i;
 42:       break;
 43:     }
 44:   }
 45:   if (PetscUnlikely(i==coneSize)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Point %D: starting point of reference cone not found in slave cone", p);
 46:   reverse_ = PETSC_FALSE;
 47:   for (i=0; i<masterConeSize; i++) {if (cone[(start_+i)%coneSize] != masterCone[i]) break;}
 48:   if (i != masterConeSize) {
 49:     reverse_ = PETSC_TRUE;
 50:     for (i=0; i<masterConeSize; i++) {if (cone[(coneSize+start_-i)%coneSize] != masterCone[i]) break;}
 51:     if (i < masterConeSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Point %D: cone has non-conforming order of points with respect to reference cone", p);
 52:   }
 53:   if (start) *start = start_;
 54:   if (reverse) *reverse = reverse_;
 55:   if (PetscUnlikely(cone[start_] != masterCone[0])) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D: cone[%d] = %d != %d = masterCone[0]", p, start_, cone[start_], masterCone[0]);
 56:   return(0);
 57: }

 59: /*@
 60:   DMPlexOrientCell - Set the desired order of cone points of this DAG point, and fix orientations accordingly.

 62:   Not Collective

 64:   Input Parameters:
 65: + dm              - The DM
 66: . p               - The DAG point (from interval given by DMPlexGetChart())
 67: . masterConeSize  - Number of specified cone points (at least 2)
 68: - masterCone      - Specified cone points, i.e. ordered subset of current cone in DAG numbering (not cone-local numbering)

 70:   Level: intermediate

 72: .seealso: DMPlexOrient(), DMPlexGetCone(), DMPlexGetConeOrientation(), DMPlexInterpolate(), DMPlexGetChart()
 73: @*/
 74: PetscErrorCode DMPlexOrientCell(DM dm, PetscInt p, PetscInt masterConeSize, const PetscInt masterCone[])
 75: {
 76:   PetscInt        coneSize;
 77:   PetscInt        start1=0;
 78:   PetscBool       reverse1=PETSC_FALSE;
 79:   PetscErrorCode  ierr;

 84:   if (masterConeSize == 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "masterConeSize cannot be 1");
 85:   DMPlexGetConeSize(dm, p, &coneSize);
 86:   if (!coneSize) return(0); /* do nothing for points with no cone */
 87:   DMPlexCompareOrientations(dm, p, masterConeSize, masterCone, &start1, &reverse1);
 88:   DMPlexOrientCell_Internal(dm, p, start1, reverse1);
 89: #if defined(PETSC_USE_DEBUG)
 90:   {
 91:     PetscInt        c;
 92:     const PetscInt *cone;
 93:     DMPlexGetCone(dm, p, &cone);
 94:     for (c = 0; c < masterConeSize; c++) {
 95:       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);
 96:     }
 97:   }
 98: #endif
 99:   return(0);
100: }

102: PetscErrorCode DMPlexOrientCell_Internal(DM dm, PetscInt p, PetscInt start1, PetscBool reverse1)
103: {
104:   PetscInt        i, j, k, maxConeSize, coneSize, coneConeSize, supportSize, supportConeSize;
105:   PetscInt        start0, start;
106:   PetscBool       reverse0, reverse;
107:   PetscInt        newornt;
108:   const PetscInt *cone=NULL, *support=NULL, *supportCone=NULL, *ornts=NULL;
109:   PetscInt       *newcone=NULL, *newornts=NULL;
110:   PetscErrorCode  ierr;

113:   if (!start1 && !reverse1) return(0);
114:   DMPlexGetConeSize(dm, p, &coneSize);
115:   if (!coneSize) return(0); /* do nothing for points with no cone */
116:   DMPlexGetCone(dm, p, &cone);
117:   DMPlexGetMaxSizes(dm, &maxConeSize, NULL);
118:   /* permute p's cone and orientations */
119:   DMPlexGetConeOrientation(dm, p, &ornts);
120:   DMGetWorkArray(dm, maxConeSize, MPIU_INT, &newcone);
121:   DMGetWorkArray(dm, maxConeSize, MPIU_INT, &newornts);
122:   DMPlexFixFaceOrientations_Permute_Private(coneSize, cone, start1, reverse1, newcone);
123:   DMPlexFixFaceOrientations_Permute_Private(coneSize, ornts, start1, reverse1, newornts);
124:   /* if direction of p (face) is flipped, flip also p's cone points (edges) */
125:   if (reverse1) {
126:     for (i=0; i<coneSize; i++) {
127:       DMPlexGetConeSize(dm, cone[i], &coneConeSize);
128:       DMPlexFixFaceOrientations_Translate_Private(newornts[i], &start0, &reverse0);
129:       DMPlexFixFaceOrientations_Combine_Private(coneConeSize, start0, reverse0, 1, PETSC_FALSE, &start, &reverse);
130:       DMPlexFixFaceOrientations_TranslateBack_Private(coneConeSize, start, reverse, &newornts[i]);
131:     }
132:   }
133:   DMPlexSetConeOrientation(dm, p, newornts);
134:   /* fix oriention of p within cones of p's support points */
135:   DMPlexGetSupport(dm, p, &support);
136:   DMPlexGetSupportSize(dm, p, &supportSize);
137:   for (j=0; j<supportSize; j++) {
138:     DMPlexGetCone(dm, support[j], &supportCone);
139:     DMPlexGetConeSize(dm, support[j], &supportConeSize);
140:     DMPlexGetConeOrientation(dm, support[j], &ornts);
141:     for (k=0; k<supportConeSize; k++) {
142:       if (supportCone[k] != p) continue;
143:       DMPlexFixFaceOrientations_Translate_Private(ornts[k], &start0, &reverse0);
144:       DMPlexFixFaceOrientations_Combine_Private(coneSize, start0, reverse0, start1, reverse1, &start, &reverse);
145:       DMPlexFixFaceOrientations_TranslateBack_Private(coneSize, start, reverse, &newornt);
146:       DMPlexInsertConeOrientation(dm, support[j], k, newornt);
147:     }
148:   }
149:   /* rewrite cone */
150:   DMPlexSetCone(dm, p, newcone);
151:   DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &newcone);
152:   DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &newornts);
153:   return(0);
154: }

156: /*@
157:   DMPlexReverseCell - Give a mesh cell the opposite orientation

159:   Input Parameters:
160: + dm   - The DM
161: - cell - The cell number

163:   Note: The modification of the DM is done in-place.

165:   Level: advanced

167: .seealso: DMPlexOrient(), DMCreate(), DMPLEX
168: @*/
169: PetscErrorCode DMPlexReverseCell(DM dm, PetscInt cell)
170: {
171:   /* Note that the reverse orientation ro of a face with orientation o is:

173:        ro = o >= 0 ? -(faceSize - o) : faceSize + o

175:      where faceSize is the size of the cone for the face.
176:   */
177:   const PetscInt *cone,    *coneO, *support;
178:   PetscInt       *revcone, *revconeO;
179:   PetscInt        maxConeSize, coneSize, supportSize, faceSize, cp, sp;
180:   PetscErrorCode  ierr;

183:   DMPlexGetMaxSizes(dm, &maxConeSize, NULL);
184:   DMGetWorkArray(dm, maxConeSize, MPIU_INT, &revcone);
185:   DMGetWorkArray(dm, maxConeSize, MPIU_INT, &revconeO);
186:   /* Reverse cone, and reverse orientations of faces */
187:   DMPlexGetConeSize(dm, cell, &coneSize);
188:   DMPlexGetCone(dm, cell, &cone);
189:   DMPlexGetConeOrientation(dm, cell, &coneO);
190:   for (cp = 0; cp < coneSize; ++cp) {
191:     const PetscInt rcp = coneSize-cp-1;

193:     DMPlexGetConeSize(dm, cone[rcp], &faceSize);
194:     revcone[cp]  = cone[rcp];
195:     revconeO[cp] = coneO[rcp] >= 0 ? -(faceSize-coneO[rcp]) : faceSize+coneO[rcp];
196:   }
197:   DMPlexSetCone(dm, cell, revcone);
198:   DMPlexSetConeOrientation(dm, cell, revconeO);
199:   /* Reverse orientation of this cell in the support hypercells */
200:   faceSize = coneSize;
201:   DMPlexGetSupportSize(dm, cell, &supportSize);
202:   DMPlexGetSupport(dm, cell, &support);
203:   for (sp = 0; sp < supportSize; ++sp) {
204:     DMPlexGetConeSize(dm, support[sp], &coneSize);
205:     DMPlexGetCone(dm, support[sp], &cone);
206:     DMPlexGetConeOrientation(dm, support[sp], &coneO);
207:     for (cp = 0; cp < coneSize; ++cp) {
208:       if (cone[cp] != cell) continue;
209:       DMPlexInsertConeOrientation(dm, support[sp], cp, coneO[cp] >= 0 ? -(faceSize-coneO[cp]) : faceSize+coneO[cp]);
210:     }
211:   }
212:   DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &revcone);
213:   DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &revconeO);
214:   return(0);
215: }

217: /*
218:   - Checks face match
219:     - Flips non-matching
220:   - Inserts faces of support cells in FIFO
221: */
222: static PetscErrorCode DMPlexCheckFace_Internal(DM dm, PetscInt *faceFIFO, PetscInt *fTop, PetscInt *fBottom, PetscInt cStart, PetscInt fStart, PetscInt fEnd, PetscBT seenCells, PetscBT flippedCells, PetscBT seenFaces)
223: {
224:   const PetscInt *support, *coneA, *coneB, *coneOA, *coneOB;
225:   PetscInt        supportSize, coneSizeA, coneSizeB, posA = -1, posB = -1;
226:   PetscInt        face, dim, seenA, flippedA, seenB, flippedB, mismatch, c;
227:   PetscErrorCode  ierr;

230:   face = faceFIFO[(*fTop)++];
231:   DMGetDimension(dm, &dim);
232:   DMPlexGetSupportSize(dm, face, &supportSize);
233:   DMPlexGetSupport(dm, face, &support);
234:   if (supportSize < 2) return(0);
235:   if (supportSize != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %d", supportSize);
236:   seenA    = PetscBTLookup(seenCells,    support[0]-cStart);
237:   flippedA = PetscBTLookup(flippedCells, support[0]-cStart) ? 1 : 0;
238:   seenB    = PetscBTLookup(seenCells,    support[1]-cStart);
239:   flippedB = PetscBTLookup(flippedCells, support[1]-cStart) ? 1 : 0;

241:   DMPlexGetConeSize(dm, support[0], &coneSizeA);
242:   DMPlexGetConeSize(dm, support[1], &coneSizeB);
243:   DMPlexGetCone(dm, support[0], &coneA);
244:   DMPlexGetCone(dm, support[1], &coneB);
245:   DMPlexGetConeOrientation(dm, support[0], &coneOA);
246:   DMPlexGetConeOrientation(dm, support[1], &coneOB);
247:   for (c = 0; c < coneSizeA; ++c) {
248:     if (!PetscBTLookup(seenFaces, coneA[c]-fStart)) {
249:       faceFIFO[(*fBottom)++] = coneA[c];
250:       PetscBTSet(seenFaces, coneA[c]-fStart);
251:     }
252:     if (coneA[c] == face) posA = c;
253:     if (*fBottom > fEnd-fStart) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %d was pushed exceeding capacity %d > %d", coneA[c], *fBottom, fEnd-fStart);
254:   }
255:   if (posA < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d could not be located in cell %d", face, support[0]);
256:   for (c = 0; c < coneSizeB; ++c) {
257:     if (!PetscBTLookup(seenFaces, coneB[c]-fStart)) {
258:       faceFIFO[(*fBottom)++] = coneB[c];
259:       PetscBTSet(seenFaces, coneB[c]-fStart);
260:     }
261:     if (coneB[c] == face) posB = c;
262:     if (*fBottom > fEnd-fStart) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %d was pushed exceeding capacity %d > %d", coneA[c], *fBottom, fEnd-fStart);
263:   }
264:   if (posB < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d could not be located in cell %d", face, support[1]);

266:   if (dim == 1) {
267:     mismatch = posA == posB;
268:   } else {
269:     mismatch = coneOA[posA] == coneOB[posB];
270:   }

272:   if (mismatch ^ (flippedA ^ flippedB)) {
273:     if (seenA && seenB) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Previously seen cells %d and %d do not match: Fault mesh is non-orientable", support[0], support[1]);
274:     if (!seenA && !flippedA) {
275:       PetscBTSet(flippedCells, support[0]-cStart);
276:     } else if (!seenB && !flippedB) {
277:       PetscBTSet(flippedCells, support[1]-cStart);
278:     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
279:   } else if (mismatch && flippedA && flippedB) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
280:   PetscBTSet(seenCells, support[0]-cStart);
281:   PetscBTSet(seenCells, support[1]-cStart);
282:   return(0);
283: }

285: /*@
286:   DMPlexOrient - Give a consistent orientation to the input mesh

288:   Input Parameters:
289: . dm - The DM

291:   Note: The orientation data for the DM are change in-place.
292: $ This routine will fail for non-orientable surfaces, such as the Moebius strip.

294:   Level: advanced

296: .seealso: DMCreate(), DMPLEX
297: @*/
298: PetscErrorCode DMPlexOrient(DM dm)
299: {
300:   MPI_Comm           comm;
301:   PetscSF            sf;
302:   const PetscInt    *lpoints;
303:   const PetscSFNode *rpoints;
304:   PetscSFNode       *rorntComp = NULL, *lorntComp = NULL;
305:   PetscInt          *numNeighbors, **neighbors;
306:   PetscSFNode       *nrankComp;
307:   PetscBool         *match, *flipped;
308:   PetscBT            seenCells, flippedCells, seenFaces;
309:   PetscInt          *faceFIFO, fTop, fBottom, *cellComp, *faceComp;
310:   PetscInt           numLeaves, numRoots, dim, h, cStart, cEnd, c, cell, fStart, fEnd, face, off, totNeighbors = 0;
311:   PetscMPIInt        rank, size, numComponents, comp = 0;
312:   PetscBool          flg, flg2;
313:   PetscViewer        viewer = NULL, selfviewer = NULL;
314:   PetscErrorCode     ierr;

317:   PetscObjectGetComm((PetscObject) dm, &comm);
318:   MPI_Comm_rank(comm, &rank);
319:   MPI_Comm_size(comm, &size);
320:   PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view", &flg);
321:   PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view_synchronized", &flg2);
322:   DMGetPointSF(dm, &sf);
323:   PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints);
324:   /* Truth Table
325:      mismatch    flips   do action   mismatch   flipA ^ flipB   action
326:          F       0 flips     no         F             F           F
327:          F       1 flip      yes        F             T           T
328:          F       2 flips     no         T             F           T
329:          T       0 flips     yes        T             T           F
330:          T       1 flip      no
331:          T       2 flips     yes
332:   */
333:   DMGetDimension(dm, &dim);
334:   DMPlexGetVTKCellHeight(dm, &h);
335:   DMPlexGetHeightStratum(dm, h,   &cStart, &cEnd);
336:   DMPlexGetHeightStratum(dm, h+1, &fStart, &fEnd);
337:   PetscBTCreate(cEnd - cStart, &seenCells);
338:   PetscBTMemzero(cEnd - cStart, seenCells);
339:   PetscBTCreate(cEnd - cStart, &flippedCells);
340:   PetscBTMemzero(cEnd - cStart, flippedCells);
341:   PetscBTCreate(fEnd - fStart, &seenFaces);
342:   PetscBTMemzero(fEnd - fStart, seenFaces);
343:   PetscCalloc3(fEnd - fStart, &faceFIFO, cEnd-cStart, &cellComp, fEnd-fStart, &faceComp);
344:   /*
345:    OLD STYLE
346:    - Add an integer array over cells and faces (component) for connected component number
347:    Foreach component
348:      - Mark the initial cell as seen
349:      - Process component as usual
350:      - Set component for all seenCells
351:      - Wipe seenCells and seenFaces (flippedCells can stay)
352:    - Generate parallel adjacency for component using SF and seenFaces
353:    - Collect numComponents adj data from each proc to 0
354:    - Build same serial graph
355:    - Use same solver
356:    - Use Scatterv to to send back flipped flags for each component
357:    - Negate flippedCells by component

359:    NEW STYLE
360:    - Create the adj on each process
361:    - Bootstrap to complete graph on proc 0
362:   */
363:   /* Loop over components */
364:   for (cell = cStart; cell < cEnd; ++cell) cellComp[cell-cStart] = -1;
365:   do {
366:     /* Look for first unmarked cell */
367:     for (cell = cStart; cell < cEnd; ++cell) if (cellComp[cell-cStart] < 0) break;
368:     if (cell >= cEnd) break;
369:     /* Initialize FIFO with first cell in component */
370:     {
371:       const PetscInt *cone;
372:       PetscInt        coneSize;

374:       fTop = fBottom = 0;
375:       DMPlexGetConeSize(dm, cell, &coneSize);
376:       DMPlexGetCone(dm, cell, &cone);
377:       for (c = 0; c < coneSize; ++c) {
378:         faceFIFO[fBottom++] = cone[c];
379:         PetscBTSet(seenFaces, cone[c]-fStart);
380:       }
381:       PetscBTSet(seenCells, cell-cStart);
382:     }
383:     /* Consider each face in FIFO */
384:     while (fTop < fBottom) {
385:       DMPlexCheckFace_Internal(dm, faceFIFO, &fTop, &fBottom, cStart, fStart, fEnd, seenCells, flippedCells, seenFaces);
386:     }
387:     /* Set component for cells and faces */
388:     for (cell = 0; cell < cEnd-cStart; ++cell) {
389:       if (PetscBTLookup(seenCells, cell)) cellComp[cell] = comp;
390:     }
391:     for (face = 0; face < fEnd-fStart; ++face) {
392:       if (PetscBTLookup(seenFaces, face)) faceComp[face] = comp;
393:     }
394:     /* Wipe seenCells and seenFaces for next component */
395:     PetscBTMemzero(fEnd - fStart, seenFaces);
396:     PetscBTMemzero(cEnd - cStart, seenCells);
397:     ++comp;
398:   } while (1);
399:   numComponents = comp;
400:   if (flg) {
401:     PetscViewer v;

403:     PetscViewerASCIIGetStdout(comm, &v);
404:     PetscViewerASCIIPushSynchronized(v);
405:     PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank);
406:     PetscBTView(cEnd-cStart, flippedCells, v);
407:     PetscViewerFlush(v);
408:     PetscViewerASCIIPopSynchronized(v);
409:   }
410:   /* Now all subdomains are oriented, but we need a consistent parallel orientation */
411:   if (numLeaves >= 0) {
412:     /* Store orientations of boundary faces*/
413:     PetscCalloc2(numRoots,&rorntComp,numRoots,&lorntComp);
414:     for (face = fStart; face < fEnd; ++face) {
415:       const PetscInt *cone, *support, *ornt;
416:       PetscInt        coneSize, supportSize;

418:       DMPlexGetSupportSize(dm, face, &supportSize);
419:       if (supportSize != 1) continue;
420:       DMPlexGetSupport(dm, face, &support);

422:       DMPlexGetCone(dm, support[0], &cone);
423:       DMPlexGetConeSize(dm, support[0], &coneSize);
424:       DMPlexGetConeOrientation(dm, support[0], &ornt);
425:       for (c = 0; c < coneSize; ++c) if (cone[c] == face) break;
426:       if (dim == 1) {
427:         /* Use cone position instead, shifted to -1 or 1 */
428:         if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = 1-c*2;
429:         else                                                rorntComp[face].rank = c*2-1;
430:       } else {
431:         if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 :  1;
432:         else                                                rorntComp[face].rank = ornt[c] < 0 ?  1 : -1;
433:       }
434:       rorntComp[face].index = faceComp[face-fStart];
435:     }
436:     /* Communicate boundary edge orientations */
437:     PetscSFBcastBegin(sf, MPIU_2INT, rorntComp, lorntComp);
438:     PetscSFBcastEnd(sf, MPIU_2INT, rorntComp, lorntComp);
439:   }
440:   /* Get process adjacency */
441:   PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors);
442:   viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm));
443:   if (flg2) {PetscViewerASCIIPushSynchronized(viewer);}
444:   PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&selfviewer);
445:   for (comp = 0; comp < numComponents; ++comp) {
446:     PetscInt  l, n;

448:     numNeighbors[comp] = 0;
449:     PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp]);
450:     /* I know this is p^2 time in general, but for bounded degree its alright */
451:     for (l = 0; l < numLeaves; ++l) {
452:       const PetscInt face = lpoints[l];

454:       /* Find a representative face (edge) separating pairs of procs */
455:       if ((face >= fStart) && (face < fEnd) && (faceComp[face-fStart] == comp)) {
456:         const PetscInt rrank = rpoints[l].rank;
457:         const PetscInt rcomp = lorntComp[face].index;

459:         for (n = 0; n < numNeighbors[comp]; ++n) if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break;
460:         if (n >= numNeighbors[comp]) {
461:           PetscInt supportSize;

463:           DMPlexGetSupportSize(dm, face, &supportSize);
464:           if (supportSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %d", supportSize);
465:           if (flg) {PetscViewerASCIIPrintf(selfviewer, "[%d]: component %d, Found representative leaf %d (face %d) connecting to face %d on (%d, %d) with orientation %d\n", rank, comp, l, face, rpoints[l].index, rrank, rcomp, lorntComp[face].rank);}
466:           neighbors[comp][numNeighbors[comp]++] = l;
467:         }
468:       }
469:     }
470:     totNeighbors += numNeighbors[comp];
471:   }
472:   PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&selfviewer);
473:   PetscViewerFlush(viewer);
474:   if (flg2) {PetscViewerASCIIPopSynchronized(viewer);}
475:   PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match);
476:   for (comp = 0, off = 0; comp < numComponents; ++comp) {
477:     PetscInt n;

479:     for (n = 0; n < numNeighbors[comp]; ++n, ++off) {
480:       const PetscInt face = lpoints[neighbors[comp][n]];
481:       const PetscInt o    = rorntComp[face].rank*lorntComp[face].rank;

483:       if      (o < 0) match[off] = PETSC_TRUE;
484:       else if (o > 0) match[off] = PETSC_FALSE;
485:       else SETERRQ5(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid face %d (%d, %d) neighbor: %d comp: %d", face, rorntComp[face], lorntComp[face], neighbors[comp][n], comp);
486:       nrankComp[off].rank  = rpoints[neighbors[comp][n]].rank;
487:       nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index;
488:     }
489:     PetscFree(neighbors[comp]);
490:   }
491:   /* Collect the graph on 0 */
492:   if (numLeaves >= 0) {
493:     Mat          G;
494:     PetscBT      seenProcs, flippedProcs;
495:     PetscInt    *procFIFO, pTop, pBottom;
496:     PetscInt    *N   = NULL, *Noff;
497:     PetscSFNode *adj = NULL;
498:     PetscBool   *val = NULL;
499:     PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc, p, o;
500:     PetscMPIInt  size = 0;

502:     PetscCalloc1(numComponents, &flipped);
503:     if (!rank) {MPI_Comm_size(comm, &size);}
504:     PetscCalloc4(size, &recvcounts, size+1, &displs, size, &Nc, size+1, &Noff);
505:     MPI_Gather(&numComponents, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm);
506:     for (p = 0; p < size; ++p) {
507:       displs[p+1] = displs[p] + Nc[p];
508:     }
509:     if (!rank) {PetscMalloc1(displs[size],&N);}
510:     MPI_Gatherv(numNeighbors, numComponents, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm);
511:     for (p = 0, o = 0; p < size; ++p) {
512:       recvcounts[p] = 0;
513:       for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o];
514:       displs[p+1] = displs[p] + recvcounts[p];
515:     }
516:     if (!rank) {PetscMalloc2(displs[size], &adj, displs[size], &val);}
517:     MPI_Gatherv(nrankComp, totNeighbors, MPIU_2INT, adj, recvcounts, displs, MPIU_2INT, 0, comm);
518:     MPI_Gatherv(match, totNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm);
519:     PetscFree2(numNeighbors, neighbors);
520:     if (!rank) {
521:       for (p = 1; p <= size; ++p) {Noff[p] = Noff[p-1] + Nc[p-1];}
522:       if (flg) {
523:         PetscInt n;

525:         for (p = 0, off = 0; p < size; ++p) {
526:           for (c = 0; c < Nc[p]; ++c) {
527:             PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %d:\n", p, c);
528:             for (n = 0; n < N[Noff[p]+c]; ++n, ++off) {
529:               PetscPrintf(PETSC_COMM_SELF, "  edge (%d, %d) (%d):\n", adj[off].rank, adj[off].index, val[off]);
530:             }
531:           }
532:         }
533:       }
534:       /* Symmetrize the graph */
535:       MatCreate(PETSC_COMM_SELF, &G);
536:       MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size]);
537:       MatSetUp(G);
538:       for (p = 0, off = 0; p < size; ++p) {
539:         for (c = 0; c < Nc[p]; ++c) {
540:           const PetscInt r = Noff[p]+c;
541:           PetscInt       n;

543:           for (n = 0; n < N[r]; ++n, ++off) {
544:             const PetscInt    q = Noff[adj[off].rank] + adj[off].index;
545:             const PetscScalar o = val[off] ? 1.0 : 0.0;

547:             MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES);
548:             MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES);
549:           }
550:         }
551:       }
552:       MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);
553:       MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);

555:       PetscBTCreate(Noff[size], &seenProcs);
556:       PetscBTMemzero(Noff[size], seenProcs);
557:       PetscBTCreate(Noff[size], &flippedProcs);
558:       PetscBTMemzero(Noff[size], flippedProcs);
559:       PetscMalloc1(Noff[size], &procFIFO);
560:       pTop = pBottom = 0;
561:       for (p = 0; p < Noff[size]; ++p) {
562:         if (PetscBTLookup(seenProcs, p)) continue;
563:         /* Initialize FIFO with next proc */
564:         procFIFO[pBottom++] = p;
565:         PetscBTSet(seenProcs, p);
566:         /* Consider each proc in FIFO */
567:         while (pTop < pBottom) {
568:           const PetscScalar *ornt;
569:           const PetscInt    *neighbors;
570:           PetscInt           proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n;

572:           proc     = procFIFO[pTop++];
573:           flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0;
574:           MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt);
575:           /* Loop over neighboring procs */
576:           for (n = 0; n < numNeighbors; ++n) {
577:             nproc    = neighbors[n];
578:             mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1;
579:             seen     = PetscBTLookup(seenProcs, nproc);
580:             flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0;

582:             if (mismatch ^ (flippedA ^ flippedB)) {
583:               if (seen) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Previously seen procs %d and %d do not match: Fault mesh is non-orientable", proc, nproc);
584:               if (!flippedB) {
585:                 PetscBTSet(flippedProcs, nproc);
586:               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
587:             } else if (mismatch && flippedA && flippedB) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
588:             if (!seen) {
589:               procFIFO[pBottom++] = nproc;
590:               PetscBTSet(seenProcs, nproc);
591:             }
592:           }
593:         }
594:       }
595:       PetscFree(procFIFO);
596:       MatDestroy(&G);
597:       PetscFree2(adj, val);
598:       PetscBTDestroy(&seenProcs);
599:     }
600:     /* Scatter flip flags */
601:     {
602:       PetscBool *flips = NULL;

604:       if (!rank) {
605:         PetscMalloc1(Noff[size], &flips);
606:         for (p = 0; p < Noff[size]; ++p) {
607:           flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE;
608:           if (flg && flips[p]) {PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p);}
609:         }
610:         for (p = 0; p < size; ++p) {
611:           displs[p+1] = displs[p] + Nc[p];
612:         }
613:       }
614:       MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, numComponents, MPIU_BOOL, 0, comm);
615:       PetscFree(flips);
616:     }
617:     if (!rank) {PetscBTDestroy(&flippedProcs);}
618:     PetscFree(N);
619:     PetscFree4(recvcounts, displs, Nc, Noff);
620:     PetscFree2(nrankComp, match);

622:     /* Decide whether to flip cells in each component */
623:     for (c = 0; c < cEnd-cStart; ++c) {if (flipped[cellComp[c]]) {PetscBTNegate(flippedCells, c);}}
624:     PetscFree(flipped);
625:   }
626:   if (flg) {
627:     PetscViewer v;

629:     PetscViewerASCIIGetStdout(comm, &v);
630:     PetscViewerASCIIPushSynchronized(v);
631:     PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank);
632:     PetscBTView(cEnd-cStart, flippedCells, v);
633:     PetscViewerFlush(v);
634:     PetscViewerASCIIPopSynchronized(v);
635:   }
636:   /* Reverse flipped cells in the mesh */
637:   for (c = cStart; c < cEnd; ++c) {
638:     if (PetscBTLookup(flippedCells, c-cStart)) {
639:       DMPlexReverseCell(dm, c);
640:     }
641:   }
642:   PetscBTDestroy(&seenCells);
643:   PetscBTDestroy(&flippedCells);
644:   PetscBTDestroy(&seenFaces);
645:   PetscFree2(numNeighbors, neighbors);
646:   PetscFree2(rorntComp, lorntComp);
647:   PetscFree3(faceFIFO, cellComp, faceComp);
648:   return(0);
649: }