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