xref: /petsc/src/dm/impls/plex/plexorient.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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 %" PetscInt_FMT " not in [%" PetscInt_FMT ",%" PetscInt_FMT ") for %s %" PetscInt_FMT, 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 %" PetscInt_FMT, 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 %" PetscInt_FMT " was pushed exceeding capacity %" PetscInt_FMT " > %" PetscInt_FMT, coneA[c], *fBottom, fEnd - fStart);
103   }
104   PetscCheck(posA >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " could not be located in cell %" PetscInt_FMT, 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 %" PetscInt_FMT " was pushed exceeding capacity %" PetscInt_FMT " > %" PetscInt_FMT, coneA[c], *fBottom, fEnd - fStart);
112   }
113   PetscCheck(posB >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " could not be located in cell %" PetscInt_FMT, 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 %" PetscInt_FMT " and %" PetscInt_FMT " 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 PetscCheck(!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, *locSupport = NULL;
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)
216       if (cellComp[cell - cStart] < 0) break;
217     if (cell >= cEnd) break;
218     /* Initialize FIFO with first cell in component */
219     {
220       const PetscInt *cone;
221       PetscInt        coneSize;
222 
223       fTop = fBottom = 0;
224       PetscCall(DMPlexGetConeSize(dm, cell, &coneSize));
225       PetscCall(DMPlexGetCone(dm, cell, &cone));
226       for (c = 0; c < coneSize; ++c) {
227         faceFIFO[fBottom++] = cone[c];
228         PetscCall(PetscBTSet(seenFaces, cone[c] - fStart));
229       }
230       PetscCall(PetscBTSet(seenCells, cell - cStart));
231     }
232     /* Consider each face in FIFO */
233     while (fTop < fBottom) PetscCall(DMPlexCheckFace_Internal(dm, faceFIFO, &fTop, &fBottom, cStart, fStart, fEnd, seenCells, flippedCells, seenFaces));
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     PetscCall(PetscBTMemzero(fEnd - fStart, seenFaces));
243     PetscCall(PetscBTMemzero(cEnd - cStart, seenCells));
244     ++comp;
245   } while (1);
246   numComponents = comp;
247   if (flg) {
248     PetscViewer v;
249 
250     PetscCall(PetscViewerASCIIGetStdout(comm, &v));
251     PetscCall(PetscViewerASCIIPushSynchronized(v));
252     PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank));
253     PetscCall(PetscBTView(cEnd - cStart, flippedCells, v));
254     PetscCall(PetscViewerFlush(v));
255     PetscCall(PetscViewerASCIIPopSynchronized(v));
256   }
257   /* Now all subdomains are oriented, but we need a consistent parallel orientation */
258   if (numLeaves >= 0) {
259     PetscInt maxSupportSize, neighbor;
260 
261     /* Store orientations of boundary faces*/
262     PetscCall(DMPlexGetMaxSizes(dm, NULL, &maxSupportSize));
263     PetscCall(PetscCalloc3(numRoots, &rorntComp, numRoots, &lorntComp, maxSupportSize, &locSupport));
264     for (face = fStart; face < fEnd; ++face) {
265       const PetscInt *cone, *support, *ornt;
266       PetscInt        coneSize, supportSize, Ns = 0, s, l;
267 
268       PetscCall(DMPlexGetSupportSize(dm, face, &supportSize));
269       /* Ignore overlapping cells */
270       PetscCall(DMPlexGetSupport(dm, face, &support));
271       for (s = 0; s < supportSize; ++s) {
272         PetscCall(PetscFindInt(support[s], numLeaves, lpoints, &l));
273         if (l >= 0) continue;
274         locSupport[Ns++] = support[s];
275       }
276       if (Ns != 1) continue;
277       neighbor = locSupport[0];
278       PetscCall(DMPlexGetCone(dm, neighbor, &cone));
279       PetscCall(DMPlexGetConeSize(dm, neighbor, &coneSize));
280       PetscCall(DMPlexGetConeOrientation(dm, neighbor, &ornt));
281       for (c = 0; c < coneSize; ++c)
282         if (cone[c] == face) break;
283       if (dim == 1) {
284         /* Use cone position instead, shifted to -1 or 1 */
285         if (PetscBTLookup(flippedCells, neighbor - cStart)) rorntComp[face].rank = 1 - c * 2;
286         else rorntComp[face].rank = c * 2 - 1;
287       } else {
288         if (PetscBTLookup(flippedCells, neighbor - cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 : 1;
289         else rorntComp[face].rank = ornt[c] < 0 ? 1 : -1;
290       }
291       rorntComp[face].index = faceComp[face - fStart];
292     }
293     /* Communicate boundary edge orientations */
294     PetscCall(PetscSFBcastBegin(sf, MPIU_2INT, rorntComp, lorntComp, MPI_REPLACE));
295     PetscCall(PetscSFBcastEnd(sf, MPIU_2INT, rorntComp, lorntComp, MPI_REPLACE));
296   }
297   /* Get process adjacency */
298   PetscCall(PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors));
299   viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm));
300   if (flg2) PetscCall(PetscViewerASCIIPushSynchronized(viewer));
301   PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &selfviewer));
302   for (comp = 0; comp < numComponents; ++comp) {
303     PetscInt l, n;
304 
305     numNeighbors[comp] = 0;
306     PetscCall(PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp]));
307     /* I know this is p^2 time in general, but for bounded degree its alright */
308     for (l = 0; l < numLeaves; ++l) {
309       const PetscInt face = lpoints[l];
310 
311       /* Find a representative face (edge) separating pairs of procs */
312       if ((face >= fStart) && (face < fEnd) && (faceComp[face - fStart] == comp) && rorntComp[face].rank) {
313         const PetscInt rrank = rpoints[l].rank;
314         const PetscInt rcomp = lorntComp[face].index;
315 
316         for (n = 0; n < numNeighbors[comp]; ++n)
317           if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break;
318         if (n >= numNeighbors[comp]) {
319           PetscInt supportSize;
320 
321           PetscCall(DMPlexGetSupportSize(dm, face, &supportSize));
322           PetscCheck(supportSize == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %" PetscInt_FMT, supportSize);
323           if (flg)
324             PetscCall(PetscViewerASCIIPrintf(selfviewer, "[%d]: component %d, Found representative leaf %" PetscInt_FMT " (face %" PetscInt_FMT ") connecting to face %" PetscInt_FMT " on (%" PetscInt_FMT ", %" PetscInt_FMT ") with orientation %" PetscInt_FMT "\n", rank, comp, l, face,
325                                              rpoints[l].index, rrank, rcomp, lorntComp[face].rank));
326           neighbors[comp][numNeighbors[comp]++] = l;
327         }
328       }
329     }
330     totNeighbors += numNeighbors[comp];
331   }
332   PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &selfviewer));
333   PetscCall(PetscViewerFlush(viewer));
334   if (flg2) PetscCall(PetscViewerASCIIPopSynchronized(viewer));
335   PetscCall(PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match));
336   for (comp = 0, off = 0; comp < numComponents; ++comp) {
337     PetscInt n;
338 
339     for (n = 0; n < numNeighbors[comp]; ++n, ++off) {
340       const PetscInt face = lpoints[neighbors[comp][n]];
341       const PetscInt o    = rorntComp[face].rank * lorntComp[face].rank;
342 
343       if (o < 0) match[off] = PETSC_TRUE;
344       else if (o > 0) match[off] = PETSC_FALSE;
345       else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid face %" PetscInt_FMT " (%" PetscInt_FMT ", %" PetscInt_FMT ") neighbor: %" PetscInt_FMT " comp: %d", face, rorntComp[face].rank, lorntComp[face].rank, neighbors[comp][n], comp);
346       nrankComp[off].rank  = rpoints[neighbors[comp][n]].rank;
347       nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index;
348     }
349     PetscCall(PetscFree(neighbors[comp]));
350   }
351   /* Collect the graph on 0 */
352   if (numLeaves >= 0) {
353     Mat          G;
354     PetscBT      seenProcs, flippedProcs;
355     PetscInt    *procFIFO, pTop, pBottom;
356     PetscInt    *N          = NULL, *Noff;
357     PetscSFNode *adj        = NULL;
358     PetscBool   *val        = NULL;
359     PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc, p, o;
360     PetscMPIInt  size = 0;
361 
362     PetscCall(PetscCalloc1(numComponents, &flipped));
363     if (rank == 0) PetscCallMPI(MPI_Comm_size(comm, &size));
364     PetscCall(PetscCalloc4(size, &recvcounts, size + 1, &displs, size, &Nc, size + 1, &Noff));
365     PetscCallMPI(MPI_Gather(&numComponents, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm));
366     for (p = 0; p < size; ++p) displs[p + 1] = displs[p] + Nc[p];
367     if (rank == 0) PetscCall(PetscMalloc1(displs[size], &N));
368     PetscCallMPI(MPI_Gatherv(numNeighbors, numComponents, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm));
369     for (p = 0, o = 0; p < size; ++p) {
370       recvcounts[p] = 0;
371       for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o];
372       displs[p + 1] = displs[p] + recvcounts[p];
373     }
374     if (rank == 0) PetscCall(PetscMalloc2(displs[size], &adj, displs[size], &val));
375     PetscCallMPI(MPI_Gatherv(nrankComp, totNeighbors, MPIU_2INT, adj, recvcounts, displs, MPIU_2INT, 0, comm));
376     PetscCallMPI(MPI_Gatherv(match, totNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm));
377     PetscCall(PetscFree2(numNeighbors, neighbors));
378     if (rank == 0) {
379       for (p = 1; p <= size; ++p) Noff[p] = Noff[p - 1] + Nc[p - 1];
380       if (flg) {
381         PetscInt n;
382 
383         for (p = 0, off = 0; p < size; ++p) {
384           for (c = 0; c < Nc[p]; ++c) {
385             PetscCall(PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %" PetscInt_FMT ":\n", p, c));
386             for (n = 0; n < N[Noff[p] + c]; ++n, ++off) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  edge (%" PetscInt_FMT ", %" PetscInt_FMT ") (%s):\n", adj[off].rank, adj[off].index, PetscBools[val[off]]));
387           }
388         }
389       }
390       /* Symmetrize the graph */
391       PetscCall(MatCreate(PETSC_COMM_SELF, &G));
392       PetscCall(MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size]));
393       PetscCall(MatSetUp(G));
394       for (p = 0, off = 0; p < size; ++p) {
395         for (c = 0; c < Nc[p]; ++c) {
396           const PetscInt r = Noff[p] + c;
397           PetscInt       n;
398 
399           for (n = 0; n < N[r]; ++n, ++off) {
400             const PetscInt    q = Noff[adj[off].rank] + adj[off].index;
401             const PetscScalar o = val[off] ? 1.0 : 0.0;
402 
403             PetscCall(MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES));
404             PetscCall(MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES));
405           }
406         }
407       }
408       PetscCall(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY));
409       PetscCall(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY));
410 
411       PetscCall(PetscBTCreate(Noff[size], &seenProcs));
412       PetscCall(PetscBTMemzero(Noff[size], seenProcs));
413       PetscCall(PetscBTCreate(Noff[size], &flippedProcs));
414       PetscCall(PetscBTMemzero(Noff[size], flippedProcs));
415       PetscCall(PetscMalloc1(Noff[size], &procFIFO));
416       pTop = pBottom = 0;
417       for (p = 0; p < Noff[size]; ++p) {
418         if (PetscBTLookup(seenProcs, p)) continue;
419         /* Initialize FIFO with next proc */
420         procFIFO[pBottom++] = p;
421         PetscCall(PetscBTSet(seenProcs, p));
422         /* Consider each proc in FIFO */
423         while (pTop < pBottom) {
424           const PetscScalar *ornt;
425           const PetscInt    *neighbors;
426           PetscInt           proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n;
427 
428           proc     = procFIFO[pTop++];
429           flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0;
430           PetscCall(MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt));
431           /* Loop over neighboring procs */
432           for (n = 0; n < numNeighbors; ++n) {
433             nproc    = neighbors[n];
434             mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1;
435             seen     = PetscBTLookup(seenProcs, nproc);
436             flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0;
437 
438             if (mismatch ^ (flippedA ^ flippedB)) {
439               PetscCheck(!seen, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Previously seen procs %" PetscInt_FMT " and %" PetscInt_FMT " do not match: Fault mesh is non-orientable", proc, nproc);
440               if (!flippedB) {
441                 PetscCall(PetscBTSet(flippedProcs, nproc));
442               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
443             } else PetscCheck(!mismatch || !flippedA || !flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
444             if (!seen) {
445               procFIFO[pBottom++] = nproc;
446               PetscCall(PetscBTSet(seenProcs, nproc));
447             }
448           }
449         }
450       }
451       PetscCall(PetscFree(procFIFO));
452       PetscCall(MatDestroy(&G));
453       PetscCall(PetscFree2(adj, val));
454       PetscCall(PetscBTDestroy(&seenProcs));
455     }
456     /* Scatter flip flags */
457     {
458       PetscBool *flips = NULL;
459 
460       if (rank == 0) {
461         PetscCall(PetscMalloc1(Noff[size], &flips));
462         for (p = 0; p < Noff[size]; ++p) {
463           flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE;
464           if (flg && flips[p]) PetscCall(PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p));
465         }
466         for (p = 0; p < size; ++p) displs[p + 1] = displs[p] + Nc[p];
467       }
468       PetscCallMPI(MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, numComponents, MPIU_BOOL, 0, comm));
469       PetscCall(PetscFree(flips));
470     }
471     if (rank == 0) PetscCall(PetscBTDestroy(&flippedProcs));
472     PetscCall(PetscFree(N));
473     PetscCall(PetscFree4(recvcounts, displs, Nc, Noff));
474     PetscCall(PetscFree2(nrankComp, match));
475 
476     /* Decide whether to flip cells in each component */
477     for (c = 0; c < cEnd - cStart; ++c) {
478       if (flipped[cellComp[c]]) PetscCall(PetscBTNegate(flippedCells, c));
479     }
480     PetscCall(PetscFree(flipped));
481   }
482   if (flg) {
483     PetscViewer v;
484 
485     PetscCall(PetscViewerASCIIGetStdout(comm, &v));
486     PetscCall(PetscViewerASCIIPushSynchronized(v));
487     PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank));
488     PetscCall(PetscBTView(cEnd - cStart, flippedCells, v));
489     PetscCall(PetscViewerFlush(v));
490     PetscCall(PetscViewerASCIIPopSynchronized(v));
491   }
492   /* Reverse flipped cells in the mesh */
493   for (c = cStart; c < cEnd; ++c) {
494     if (PetscBTLookup(flippedCells, c - cStart)) PetscCall(DMPlexOrientPoint(dm, c, -1));
495   }
496   PetscCall(PetscBTDestroy(&seenCells));
497   PetscCall(PetscBTDestroy(&flippedCells));
498   PetscCall(PetscBTDestroy(&seenFaces));
499   PetscCall(PetscFree2(numNeighbors, neighbors));
500   PetscCall(PetscFree3(rorntComp, lorntComp, locSupport));
501   PetscCall(PetscFree3(faceFIFO, cellComp, faceComp));
502   PetscFunctionReturn(0);
503 }
504