xref: /petsc/src/dm/impls/plex/plexorient.c (revision 3952672843d878b0f648d746d6debcd9b6bd8fcc)
1 #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 #include <petscsf.h>
3 
4 /*@
5   DMPlexReverseCell - Give a mesh cell the opposite orientation
6 
7   Input Parameters:
8 + dm   - The DM
9 - cell - The cell number
10 
11   Note: The modification of the DM is done in-place.
12 
13   Level: advanced
14 
15 .seealso: DMPlexOrient(), DMCreate(), DMPLEX
16 @*/
17 PetscErrorCode DMPlexReverseCell(DM dm, PetscInt cell)
18 {
19   /* Note that the reverse orientation ro of a face with orientation o is:
20 
21        ro = o >= 0 ? -(faceSize - o) : faceSize + o
22 
23      where faceSize is the size of the cone for the face.
24   */
25   const PetscInt *cone,    *coneO, *support;
26   PetscInt       *revcone, *revconeO;
27   PetscInt        maxConeSize, coneSize, supportSize, faceSize, cp, sp;
28   PetscErrorCode  ierr;
29 
30   PetscFunctionBegin;
31   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
32   ierr = DMGetWorkArray(dm, maxConeSize, PETSC_INT, &revcone);CHKERRQ(ierr);
33   ierr = DMGetWorkArray(dm, maxConeSize, PETSC_INT, &revconeO);CHKERRQ(ierr);
34   /* Reverse cone, and reverse orientations of faces */
35   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
36   ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
37   ierr = DMPlexGetConeOrientation(dm, cell, &coneO);CHKERRQ(ierr);
38   for (cp = 0; cp < coneSize; ++cp) {
39     const PetscInt rcp = coneSize-cp-1;
40 
41     ierr = DMPlexGetConeSize(dm, cone[rcp], &faceSize);CHKERRQ(ierr);
42     revcone[cp]  = cone[rcp];
43     revconeO[cp] = coneO[rcp] >= 0 ? -(faceSize-coneO[rcp]) : faceSize+coneO[rcp];
44   }
45   ierr = DMPlexSetCone(dm, cell, revcone);CHKERRQ(ierr);
46   ierr = DMPlexSetConeOrientation(dm, cell, revconeO);CHKERRQ(ierr);
47   /* Reverse orientation of this cell in the support hypercells */
48   faceSize = coneSize;
49   ierr = DMPlexGetSupportSize(dm, cell, &supportSize);CHKERRQ(ierr);
50   ierr = DMPlexGetSupport(dm, cell, &support);CHKERRQ(ierr);
51   for (sp = 0; sp < supportSize; ++sp) {
52     ierr = DMPlexGetConeSize(dm, support[sp], &coneSize);CHKERRQ(ierr);
53     ierr = DMPlexGetCone(dm, support[sp], &cone);CHKERRQ(ierr);
54     ierr = DMPlexGetConeOrientation(dm, support[sp], &coneO);CHKERRQ(ierr);
55     for (cp = 0; cp < coneSize; ++cp) {
56       if (cone[cp] != cell) continue;
57       ierr = DMPlexInsertConeOrientation(dm, support[sp], cp, coneO[cp] >= 0 ? -(faceSize-coneO[cp]) : faceSize+coneO[cp]);CHKERRQ(ierr);
58     }
59   }
60   ierr = DMRestoreWorkArray(dm, maxConeSize, PETSC_INT, &revcone);CHKERRQ(ierr);
61   ierr = DMRestoreWorkArray(dm, maxConeSize, PETSC_INT, &revconeO);CHKERRQ(ierr);
62   PetscFunctionReturn(0);
63 }
64 
65 /*
66   - Checks face match
67     - Flips non-matching
68   - Inserts faces of support cells in FIFO
69 */
70 static PetscErrorCode DMPlexCheckFace_Internal(DM dm, PetscInt *faceFIFO, PetscInt *fTop, PetscInt *fBottom, PetscInt cStart, PetscInt fStart, PetscInt fEnd, PetscBT seenCells, PetscBT flippedCells, PetscBT seenFaces)
71 {
72   const PetscInt *support, *coneA, *coneB, *coneOA, *coneOB;
73   PetscInt        supportSize, coneSizeA, coneSizeB, posA = -1, posB = -1;
74   PetscInt        face, dim, seenA, flippedA, seenB, flippedB, mismatch, c;
75   PetscErrorCode  ierr;
76 
77   PetscFunctionBegin;
78   face = faceFIFO[(*fTop)++];
79   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
80   ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
81   ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr);
82   if (supportSize < 2) PetscFunctionReturn(0);
83   if (supportSize != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %d", supportSize);
84   seenA    = PetscBTLookup(seenCells,    support[0]-cStart);
85   flippedA = PetscBTLookup(flippedCells, support[0]-cStart) ? 1 : 0;
86   seenB    = PetscBTLookup(seenCells,    support[1]-cStart);
87   flippedB = PetscBTLookup(flippedCells, support[1]-cStart) ? 1 : 0;
88 
89   ierr = DMPlexGetConeSize(dm, support[0], &coneSizeA);CHKERRQ(ierr);
90   ierr = DMPlexGetConeSize(dm, support[1], &coneSizeB);CHKERRQ(ierr);
91   ierr = DMPlexGetCone(dm, support[0], &coneA);CHKERRQ(ierr);
92   ierr = DMPlexGetCone(dm, support[1], &coneB);CHKERRQ(ierr);
93   ierr = DMPlexGetConeOrientation(dm, support[0], &coneOA);CHKERRQ(ierr);
94   ierr = DMPlexGetConeOrientation(dm, support[1], &coneOB);CHKERRQ(ierr);
95   for (c = 0; c < coneSizeA; ++c) {
96     if (!PetscBTLookup(seenFaces, coneA[c]-fStart)) {
97       faceFIFO[(*fBottom)++] = coneA[c];
98       ierr = PetscBTSet(seenFaces, coneA[c]-fStart);CHKERRQ(ierr);
99     }
100     if (coneA[c] == face) posA = c;
101     if (*fBottom > fEnd-fStart) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %d was pushed exceeding capacity %d > %d", coneA[c], *fBottom, fEnd-fStart);
102   }
103   if (posA < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d could not be located in cell %d", face, support[0]);
104   for (c = 0; c < coneSizeB; ++c) {
105     if (!PetscBTLookup(seenFaces, coneB[c]-fStart)) {
106       faceFIFO[(*fBottom)++] = coneB[c];
107       ierr = PetscBTSet(seenFaces, coneB[c]-fStart);CHKERRQ(ierr);
108     }
109     if (coneB[c] == face) posB = c;
110     if (*fBottom > fEnd-fStart) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %d was pushed exceeding capacity %d > %d", coneA[c], *fBottom, fEnd-fStart);
111   }
112   if (posB < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d could not be located in cell %d", face, support[1]);
113 
114   if (dim == 1) {
115     mismatch = posA == posB;
116   } else {
117     mismatch = coneOA[posA] == coneOB[posB];
118   }
119 
120   if (mismatch ^ (flippedA ^ flippedB)) {
121     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]);
122     if (!seenA && !flippedA) {
123       ierr = PetscBTSet(flippedCells, support[0]-cStart);CHKERRQ(ierr);
124     } else if (!seenB && !flippedB) {
125       ierr = PetscBTSet(flippedCells, support[1]-cStart);CHKERRQ(ierr);
126     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
127   } else if (mismatch && flippedA && flippedB) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
128   ierr = PetscBTSet(seenCells, support[0]-cStart);CHKERRQ(ierr);
129   ierr = PetscBTSet(seenCells, support[1]-cStart);CHKERRQ(ierr);
130   PetscFunctionReturn(0);
131 }
132 
133 /*@
134   DMPlexOrient - Give a consistent orientation to the input mesh
135 
136   Input Parameters:
137 . dm - The DM
138 
139   Note: The orientation data for the DM are change in-place.
140 $ This routine will fail for non-orientable surfaces, such as the Moebius strip.
141 
142   Level: advanced
143 
144 .seealso: DMCreate(), DMPLEX
145 @*/
146 PetscErrorCode DMPlexOrient(DM dm)
147 {
148   MPI_Comm           comm;
149   PetscSF            sf;
150   const PetscInt    *lpoints;
151   const PetscSFNode *rpoints;
152   PetscSFNode       *rorntComp = NULL, *lorntComp = NULL;
153   PetscInt          *numNeighbors, **neighbors;
154   PetscSFNode       *nrankComp;
155   PetscBool         *match, *flipped;
156   PetscBT            seenCells, flippedCells, seenFaces;
157   PetscInt          *faceFIFO, fTop, fBottom, *cellComp, *faceComp;
158   PetscInt           numLeaves, numRoots, dim, h, cStart, cEnd, c, cell, fStart, fEnd, face, off, totNeighbors = 0;
159   PetscMPIInt        rank, size, numComponents, comp = 0;
160   PetscBool          flg, flg2;
161   PetscErrorCode     ierr;
162 
163   PetscFunctionBegin;
164   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
165   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
166   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
167   ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view", &flg);CHKERRQ(ierr);
168   ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view_synchronized", &flg2);CHKERRQ(ierr);
169   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
170   ierr = PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints);CHKERRQ(ierr);
171   /* Truth Table
172      mismatch    flips   do action   mismatch   flipA ^ flipB   action
173          F       0 flips     no         F             F           F
174          F       1 flip      yes        F             T           T
175          F       2 flips     no         T             F           T
176          T       0 flips     yes        T             T           F
177          T       1 flip      no
178          T       2 flips     yes
179   */
180   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
181   ierr = DMPlexGetVTKCellHeight(dm, &h);CHKERRQ(ierr);
182   ierr = DMPlexGetHeightStratum(dm, h,   &cStart, &cEnd);CHKERRQ(ierr);
183   ierr = DMPlexGetHeightStratum(dm, h+1, &fStart, &fEnd);CHKERRQ(ierr);
184   ierr = PetscBTCreate(cEnd - cStart, &seenCells);CHKERRQ(ierr);
185   ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr);
186   ierr = PetscBTCreate(cEnd - cStart, &flippedCells);CHKERRQ(ierr);
187   ierr = PetscBTMemzero(cEnd - cStart, flippedCells);CHKERRQ(ierr);
188   ierr = PetscBTCreate(fEnd - fStart, &seenFaces);CHKERRQ(ierr);
189   ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr);
190   ierr = PetscCalloc3(fEnd - fStart, &faceFIFO, cEnd-cStart, &cellComp, fEnd-fStart, &faceComp);CHKERRQ(ierr);
191   /*
192    OLD STYLE
193    - Add an integer array over cells and faces (component) for connected component number
194    Foreach component
195      - Mark the initial cell as seen
196      - Process component as usual
197      - Set component for all seenCells
198      - Wipe seenCells and seenFaces (flippedCells can stay)
199    - Generate parallel adjacency for component using SF and seenFaces
200    - Collect numComponents adj data from each proc to 0
201    - Build same serial graph
202    - Use same solver
203    - Use Scatterv to to send back flipped flags for each component
204    - Negate flippedCells by component
205 
206    NEW STYLE
207    - Create the adj on each process
208    - Bootstrap to complete graph on proc 0
209   */
210   /* Loop over components */
211   for (cell = cStart; cell < cEnd; ++cell) cellComp[cell-cStart] = -1;
212   do {
213     /* Look for first unmarked cell */
214     for (cell = cStart; cell < cEnd; ++cell) if (cellComp[cell-cStart] < 0) break;
215     if (cell >= cEnd) break;
216     /* Initialize FIFO with first cell in component */
217     {
218       const PetscInt *cone;
219       PetscInt        coneSize;
220 
221       fTop = fBottom = 0;
222       ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
223       ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
224       for (c = 0; c < coneSize; ++c) {
225         faceFIFO[fBottom++] = cone[c];
226         ierr = PetscBTSet(seenFaces, cone[c]-fStart);CHKERRQ(ierr);
227       }
228       ierr = PetscBTSet(seenCells, cell-cStart);CHKERRQ(ierr);
229     }
230     /* Consider each face in FIFO */
231     while (fTop < fBottom) {
232       ierr = DMPlexCheckFace_Internal(dm, faceFIFO, &fTop, &fBottom, cStart, fStart, fEnd, seenCells, flippedCells, seenFaces);CHKERRQ(ierr);
233     }
234     /* Set component for cells and faces */
235     for (cell = 0; cell < cEnd-cStart; ++cell) {
236       if (PetscBTLookup(seenCells, cell)) cellComp[cell] = comp;
237     }
238     for (face = 0; face < fEnd-fStart; ++face) {
239       if (PetscBTLookup(seenFaces, face)) faceComp[face] = comp;
240     }
241     /* Wipe seenCells and seenFaces for next component */
242     ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr);
243     ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr);
244     ++comp;
245   } while (1);
246   numComponents = comp;
247   if (flg) {
248     PetscViewer v;
249 
250     ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr);
251     ierr = PetscViewerASCIIPushSynchronized(v);CHKERRQ(ierr);
252     ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank);CHKERRQ(ierr);
253     ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr);
254     ierr = PetscViewerFlush(v);CHKERRQ(ierr);
255     ierr = PetscViewerASCIIPopSynchronized(v);CHKERRQ(ierr);
256   }
257   /* Now all subdomains are oriented, but we need a consistent parallel orientation */
258   if (numLeaves >= 0) {
259     /* Store orientations of boundary faces*/
260     ierr = PetscCalloc2(numRoots,&rorntComp,numRoots,&lorntComp);CHKERRQ(ierr);
261     for (face = fStart; face < fEnd; ++face) {
262       const PetscInt *cone, *support, *ornt;
263       PetscInt        coneSize, supportSize;
264 
265       ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
266       if (supportSize != 1) continue;
267       ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr);
268 
269       ierr = DMPlexGetCone(dm, support[0], &cone);CHKERRQ(ierr);
270       ierr = DMPlexGetConeSize(dm, support[0], &coneSize);CHKERRQ(ierr);
271       ierr = DMPlexGetConeOrientation(dm, support[0], &ornt);CHKERRQ(ierr);
272       for (c = 0; c < coneSize; ++c) if (cone[c] == face) break;
273       if (dim == 1) {
274         /* Use cone position instead, shifted to -1 or 1 */
275         if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = 1-c*2;
276         else                                                rorntComp[face].rank = c*2-1;
277       } else {
278         if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 :  1;
279         else                                                rorntComp[face].rank = ornt[c] < 0 ?  1 : -1;
280       }
281       rorntComp[face].index = faceComp[face-fStart];
282     }
283     /* Communicate boundary edge orientations */
284     ierr = PetscSFBcastBegin(sf, MPIU_2INT, rorntComp, lorntComp);CHKERRQ(ierr);
285     ierr = PetscSFBcastEnd(sf, MPIU_2INT, rorntComp, lorntComp);CHKERRQ(ierr);
286   }
287   /* Get process adjacency */
288   ierr = PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors);CHKERRQ(ierr);
289   for (comp = 0; comp < numComponents; ++comp) {
290     PetscInt  l, n;
291 
292     numNeighbors[comp] = 0;
293     ierr = PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp]);CHKERRQ(ierr);
294     /* I know this is p^2 time in general, but for bounded degree its alright */
295     for (l = 0; l < numLeaves; ++l) {
296       const PetscInt face = lpoints[l];
297 
298       /* Find a representative face (edge) separating pairs of procs */
299       if ((face >= fStart) && (face < fEnd) && (faceComp[face-fStart] == comp)) {
300         const PetscInt rrank = rpoints[l].rank;
301         const PetscInt rcomp = lorntComp[face].index;
302 
303         for (n = 0; n < numNeighbors[comp]; ++n) if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break;
304         if (n >= numNeighbors[comp]) {
305           PetscInt supportSize;
306 
307           ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
308           if (supportSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %d", supportSize);
309           if (flg && !flg2) {ierr = PetscPrintf(PETSC_COMM_SELF, "[%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);CHKERRQ(ierr);}
310           neighbors[comp][numNeighbors[comp]++] = l;
311         }
312       }
313     }
314     totNeighbors += numNeighbors[comp];
315   }
316   if (flg2) {
317     PetscInt p;
318     for (p = 0; p < size; p++) {
319       if (p == rank) {
320         for (comp = 0; comp < numComponents; ++comp) {
321           PetscInt  n;
322 
323           for (n = 0; n < numNeighbors[comp]; n++) {
324             PetscInt l    = neighbors[comp][n];
325             PetscInt face = lpoints[l];
326             PetscInt rrank = rpoints[l].rank;
327             PetscInt rcomp = lorntComp[face].index;
328 
329             ierr = PetscPrintf(PETSC_COMM_SELF, "[%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);CHKERRQ(ierr);
330           }
331         }
332       }
333       ierr = PetscBarrier((PetscObject)dm);CHKERRQ(ierr);
334     }
335   }
336   ierr = PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match);CHKERRQ(ierr);
337   for (comp = 0, off = 0; comp < numComponents; ++comp) {
338     PetscInt n;
339 
340     for (n = 0; n < numNeighbors[comp]; ++n, ++off) {
341       const PetscInt face = lpoints[neighbors[comp][n]];
342       const PetscInt o    = rorntComp[face].rank*lorntComp[face].rank;
343 
344       if      (o < 0) match[off] = PETSC_TRUE;
345       else if (o > 0) match[off] = PETSC_FALSE;
346       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);
347       nrankComp[off].rank  = rpoints[neighbors[comp][n]].rank;
348       nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index;
349     }
350     ierr = PetscFree(neighbors[comp]);CHKERRQ(ierr);
351   }
352   /* Collect the graph on 0 */
353   if (numLeaves >= 0) {
354     Mat          G;
355     PetscBT      seenProcs, flippedProcs;
356     PetscInt    *procFIFO, pTop, pBottom;
357     PetscInt    *N   = NULL, *Noff;
358     PetscSFNode *adj = NULL;
359     PetscBool   *val = NULL;
360     PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc, p, o;
361     PetscMPIInt  size = 0;
362 
363     ierr = PetscCalloc1(numComponents, &flipped);CHKERRQ(ierr);
364     if (!rank) {ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);}
365     ierr = PetscCalloc4(size, &recvcounts, size+1, &displs, size, &Nc, size+1, &Noff);CHKERRQ(ierr);
366     ierr = MPI_Gather(&numComponents, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm);CHKERRQ(ierr);
367     for (p = 0; p < size; ++p) {
368       displs[p+1] = displs[p] + Nc[p];
369     }
370     if (!rank) {ierr = PetscMalloc1(displs[size],&N);CHKERRQ(ierr);}
371     ierr = MPI_Gatherv(numNeighbors, numComponents, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm);CHKERRQ(ierr);
372     for (p = 0, o = 0; p < size; ++p) {
373       recvcounts[p] = 0;
374       for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o];
375       displs[p+1] = displs[p] + recvcounts[p];
376     }
377     if (!rank) {ierr = PetscMalloc2(displs[size], &adj, displs[size], &val);CHKERRQ(ierr);}
378     ierr = MPI_Gatherv(nrankComp, totNeighbors, MPIU_2INT, adj, recvcounts, displs, MPIU_2INT, 0, comm);CHKERRQ(ierr);
379     ierr = MPI_Gatherv(match, totNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
380     ierr = PetscFree2(numNeighbors, neighbors);CHKERRQ(ierr);
381     if (!rank) {
382       for (p = 1; p <= size; ++p) {Noff[p] = Noff[p-1] + Nc[p-1];}
383       if (flg) {
384         PetscInt n;
385 
386         for (p = 0, off = 0; p < size; ++p) {
387           for (c = 0; c < Nc[p]; ++c) {
388             ierr = PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %d:\n", p, c);CHKERRQ(ierr);
389             for (n = 0; n < N[Noff[p]+c]; ++n, ++off) {
390               ierr = PetscPrintf(PETSC_COMM_SELF, "  edge (%d, %d) (%d):\n", adj[off].rank, adj[off].index, val[off]);CHKERRQ(ierr);
391             }
392           }
393         }
394       }
395       /* Symmetrize the graph */
396       ierr = MatCreate(PETSC_COMM_SELF, &G);CHKERRQ(ierr);
397       ierr = MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size]);CHKERRQ(ierr);
398       ierr = MatSetUp(G);CHKERRQ(ierr);
399       for (p = 0, off = 0; p < size; ++p) {
400         for (c = 0; c < Nc[p]; ++c) {
401           const PetscInt r = Noff[p]+c;
402           PetscInt       n;
403 
404           for (n = 0; n < N[r]; ++n, ++off) {
405             const PetscInt    q = Noff[adj[off].rank] + adj[off].index;
406             const PetscScalar o = val[off] ? 1.0 : 0.0;
407 
408             ierr = MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES);CHKERRQ(ierr);
409             ierr = MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES);CHKERRQ(ierr);
410           }
411         }
412       }
413       ierr = MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
414       ierr = MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
415 
416       ierr = PetscBTCreate(Noff[size], &seenProcs);CHKERRQ(ierr);
417       ierr = PetscBTMemzero(Noff[size], seenProcs);CHKERRQ(ierr);
418       ierr = PetscBTCreate(Noff[size], &flippedProcs);CHKERRQ(ierr);
419       ierr = PetscBTMemzero(Noff[size], flippedProcs);CHKERRQ(ierr);
420       ierr = PetscMalloc1(Noff[size], &procFIFO);CHKERRQ(ierr);
421       pTop = pBottom = 0;
422       for (p = 0; p < Noff[size]; ++p) {
423         if (PetscBTLookup(seenProcs, p)) continue;
424         /* Initialize FIFO with next proc */
425         procFIFO[pBottom++] = p;
426         ierr = PetscBTSet(seenProcs, p);CHKERRQ(ierr);
427         /* Consider each proc in FIFO */
428         while (pTop < pBottom) {
429           const PetscScalar *ornt;
430           const PetscInt    *neighbors;
431           PetscInt           proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n;
432 
433           proc     = procFIFO[pTop++];
434           flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0;
435           ierr = MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt);CHKERRQ(ierr);
436           /* Loop over neighboring procs */
437           for (n = 0; n < numNeighbors; ++n) {
438             nproc    = neighbors[n];
439             mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1;
440             seen     = PetscBTLookup(seenProcs, nproc);
441             flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0;
442 
443             if (mismatch ^ (flippedA ^ flippedB)) {
444               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);
445               if (!flippedB) {
446                 ierr = PetscBTSet(flippedProcs, nproc);CHKERRQ(ierr);
447               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
448             } else if (mismatch && flippedA && flippedB) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
449             if (!seen) {
450               procFIFO[pBottom++] = nproc;
451               ierr = PetscBTSet(seenProcs, nproc);CHKERRQ(ierr);
452             }
453           }
454         }
455       }
456       ierr = PetscFree(procFIFO);CHKERRQ(ierr);
457       ierr = MatDestroy(&G);CHKERRQ(ierr);
458       ierr = PetscFree2(adj, val);CHKERRQ(ierr);
459       ierr = PetscBTDestroy(&seenProcs);CHKERRQ(ierr);
460     }
461     /* Scatter flip flags */
462     {
463       PetscBool *flips = NULL;
464 
465       if (!rank) {
466         ierr = PetscMalloc1(Noff[size], &flips);CHKERRQ(ierr);
467         for (p = 0; p < Noff[size]; ++p) {
468           flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE;
469           if (flg && flips[p]) {ierr = PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p);CHKERRQ(ierr);}
470         }
471         for (p = 0; p < size; ++p) {
472           displs[p+1] = displs[p] + Nc[p];
473         }
474       }
475       ierr = MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, numComponents, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
476       ierr = PetscFree(flips);CHKERRQ(ierr);
477     }
478     if (!rank) {ierr = PetscBTDestroy(&flippedProcs);CHKERRQ(ierr);}
479     ierr = PetscFree(N);CHKERRQ(ierr);
480     ierr = PetscFree4(recvcounts, displs, Nc, Noff);CHKERRQ(ierr);
481     ierr = PetscFree2(nrankComp, match);CHKERRQ(ierr);
482 
483     /* Decide whether to flip cells in each component */
484     for (c = 0; c < cEnd-cStart; ++c) {if (flipped[cellComp[c]]) {ierr = PetscBTNegate(flippedCells, c);CHKERRQ(ierr);}}
485     ierr = PetscFree(flipped);CHKERRQ(ierr);
486   }
487   if (flg) {
488     PetscViewer v;
489 
490     ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr);
491     ierr = PetscViewerASCIIPushSynchronized(v);CHKERRQ(ierr);
492     ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank);CHKERRQ(ierr);
493     ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr);
494     ierr = PetscViewerFlush(v);CHKERRQ(ierr);
495     ierr = PetscViewerASCIIPopSynchronized(v);CHKERRQ(ierr);
496   }
497   /* Reverse flipped cells in the mesh */
498   for (c = cStart; c < cEnd; ++c) {
499     if (PetscBTLookup(flippedCells, c-cStart)) {ierr = DMPlexReverseCell(dm, c);CHKERRQ(ierr);}
500   }
501   ierr = PetscBTDestroy(&seenCells);CHKERRQ(ierr);
502   ierr = PetscBTDestroy(&flippedCells);CHKERRQ(ierr);
503   ierr = PetscBTDestroy(&seenFaces);CHKERRQ(ierr);
504   ierr = PetscFree2(numNeighbors, neighbors);CHKERRQ(ierr);
505   ierr = PetscFree2(rorntComp, lorntComp);CHKERRQ(ierr);
506   ierr = PetscFree3(faceFIFO, cellComp, faceComp);CHKERRQ(ierr);
507   PetscFunctionReturn(0);
508 }
509