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