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