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