xref: /petsc/src/dm/impls/plex/plexorient.c (revision 7b310ce2c7012fba036708a25c3ea2bb81f4d0b0)
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 face, 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        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   PetscBT        seenCells, flippedCells, seenFaces;
156   PetscInt      *faceFIFO, fTop, fBottom;
157   PetscInt       dim, h, cStart, cEnd, c, fStart, fEnd, face;
158   PetscBool      flg;
159   PetscErrorCode ierr;
160 
161   PetscFunctionBegin;
162   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
163   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-orientation_view", &flg);CHKERRQ(ierr);
164   /* Truth Table
165      mismatch    flips   do action   mismatch   flipA ^ flipB   action
166          F       0 flips     no         F             F           F
167          F       1 flip      yes        F             T           T
168          F       2 flips     no         T             F           T
169          T       0 flips     yes        T             T           F
170          T       1 flip      no
171          T       2 flips     yes
172   */
173   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
174   ierr = DMPlexGetVTKCellHeight(dm, &h);CHKERRQ(ierr);
175   ierr = DMPlexGetHeightStratum(dm, h,   &cStart, &cEnd);CHKERRQ(ierr);
176   ierr = DMPlexGetHeightStratum(dm, h+1, &fStart, &fEnd);CHKERRQ(ierr);
177   ierr = PetscBTCreate(cEnd - cStart, &seenCells);CHKERRQ(ierr);
178   ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr);
179   ierr = PetscBTCreate(cEnd - cStart, &flippedCells);CHKERRQ(ierr);
180   ierr = PetscBTMemzero(cEnd - cStart, flippedCells);CHKERRQ(ierr);
181   ierr = PetscBTCreate(fEnd - fStart, &seenFaces);CHKERRQ(ierr);
182   ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr);
183   ierr = PetscMalloc1((fEnd - fStart), &faceFIFO);CHKERRQ(ierr);
184   fTop = fBottom = 0;
185   /* Initialize FIFO with first cell */
186   if (cEnd > cStart) {
187     const PetscInt *cone;
188     PetscInt        coneSize;
189 
190     ierr = DMPlexGetConeSize(dm, cStart, &coneSize);CHKERRQ(ierr);
191     ierr = DMPlexGetCone(dm, cStart, &cone);CHKERRQ(ierr);
192     for (c = 0; c < coneSize; ++c) {
193       faceFIFO[fBottom++] = cone[c];
194       ierr = PetscBTSet(seenFaces, cone[c]-fStart);CHKERRQ(ierr);
195     }
196   }
197   /* Consider each face in FIFO */
198   while (fTop < fBottom) {
199     ierr = DMPlexCheckFace_Internal(dm, face, faceFIFO, &fTop, &fBottom, cStart, fStart, fEnd, seenCells, flippedCells, seenFaces);CHKERRQ(ierr);
200   }
201   if (flg) {
202     PetscViewer v;
203     PetscMPIInt rank;
204 
205     ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
206     ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr);
207     ierr = PetscViewerASCIISynchronizedAllow(v, PETSC_TRUE);CHKERRQ(ierr);
208     ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial seen faces:\n", rank);CHKERRQ(ierr);
209     ierr = PetscBTView(fEnd-fStart, seenFaces,    v);CHKERRQ(ierr);
210     ierr = PetscViewerASCIISynchronizedAllow(v, PETSC_TRUE);CHKERRQ(ierr);
211     ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank);CHKERRQ(ierr);
212     ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr);
213   }
214   /* Now all subdomains are oriented, but we need a consistent parallel orientation */
215   {
216     /* Find a representative face (edge) separating pairs of procs */
217     PetscSF            sf;
218     const PetscInt    *lpoints;
219     const PetscSFNode *rpoints;
220     PetscInt          *neighbors, *nranks;
221     PetscInt           numLeaves, numRoots, numNeighbors = 0, l, n;
222 
223     ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
224     ierr = PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints);CHKERRQ(ierr);
225     if (numLeaves >= 0) {
226       const PetscInt *cone, *ornt, *support;
227       PetscInt        coneSize, supportSize;
228       int            *rornt, *lornt; /* PetscSF cannot handle smaller than int */
229       PetscBool      *match, flipped = PETSC_FALSE;
230 
231       ierr = PetscMalloc1(numLeaves,&neighbors);CHKERRQ(ierr);
232       /* I know this is p^2 time in general, but for bounded degree its alright */
233       for (l = 0; l < numLeaves; ++l) {
234         const PetscInt face = lpoints[l];
235         if ((face >= fStart) && (face < fEnd)) {
236           const PetscInt rank = rpoints[l].rank;
237           for (n = 0; n < numNeighbors; ++n) if (rank == rpoints[neighbors[n]].rank) break;
238           if (n >= numNeighbors) {
239             PetscInt supportSize;
240             ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
241             if (supportSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %d", supportSize);
242             neighbors[numNeighbors++] = l;
243           }
244         }
245       }
246       ierr = PetscCalloc4(numNeighbors,&match,numNeighbors,&nranks,numRoots,&rornt,numRoots,&lornt);CHKERRQ(ierr);
247       for (face = fStart; face < fEnd; ++face) {
248         ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
249         if (supportSize != 1) continue;
250         ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr);
251 
252         ierr = DMPlexGetCone(dm, support[0], &cone);CHKERRQ(ierr);
253         ierr = DMPlexGetConeSize(dm, support[0], &coneSize);CHKERRQ(ierr);
254         ierr = DMPlexGetConeOrientation(dm, support[0], &ornt);CHKERRQ(ierr);
255         for (c = 0; c < coneSize; ++c) if (cone[c] == face) break;
256         if (dim == 1) {
257           /* Use cone position instead, shifted to -1 or 1 */
258           rornt[face] = c*2-1;
259         } else {
260           if (PetscBTLookup(flippedCells, support[0]-cStart)) rornt[face] = ornt[c] < 0 ? -1 :  1;
261           else                                                rornt[face] = ornt[c] < 0 ?  1 : -1;
262         }
263       }
264       /* Mark each edge with match or nomatch */
265       ierr = PetscSFBcastBegin(sf, MPI_INT, rornt, lornt);CHKERRQ(ierr);
266       ierr = PetscSFBcastEnd(sf, MPI_INT, rornt, lornt);CHKERRQ(ierr);
267       for (n = 0; n < numNeighbors; ++n) {
268         const PetscInt face = lpoints[neighbors[n]];
269 
270         if (rornt[face]*lornt[face] < 0) match[n] = PETSC_TRUE;
271         else                             match[n] = PETSC_FALSE;
272         nranks[n] = rpoints[neighbors[n]].rank;
273       }
274       /* Collect the graph on 0 */
275       {
276         Mat          G;
277         PetscBT      seenProcs, flippedProcs;
278         PetscInt    *procFIFO, pTop, pBottom;
279         PetscInt    *adj = NULL;
280         PetscBool   *val = NULL;
281         PetscMPIInt *recvcounts = NULL, *displs = NULL, p;
282         PetscMPIInt  N = numNeighbors, numProcs = 0, rank;
283         PetscInt     debug = 0;
284 
285         ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
286         if (!rank) {ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);}
287         ierr = PetscCalloc2(numProcs,&recvcounts,numProcs+1,&displs);CHKERRQ(ierr);
288         ierr = MPI_Gather(&N, 1, MPI_INT, recvcounts, 1, MPI_INT, 0, comm);CHKERRQ(ierr);
289         for (p = 0; p < numProcs; ++p) {
290           displs[p+1] = displs[p] + recvcounts[p];
291         }
292         if (!rank) {ierr = PetscMalloc2(displs[numProcs],&adj,displs[numProcs],&val);CHKERRQ(ierr);}
293         ierr = MPI_Gatherv(nranks, numNeighbors, MPIU_INT, adj, recvcounts, displs, MPIU_INT, 0, comm);CHKERRQ(ierr);
294         ierr = MPI_Gatherv(match, numNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
295         if (debug) {
296           for (p = 0; p < numProcs; ++p) {
297             ierr = PetscPrintf(comm, "Proc %d:\n", p);
298             for (n = 0; n < recvcounts[p]; ++n) {
299               ierr = PetscPrintf(comm, "  edge %d (%d):\n", adj[displs[p]+n], val[displs[p]+n]);
300             }
301           }
302         }
303         /* Symmetrize the graph */
304         ierr = MatCreate(PETSC_COMM_SELF, &G);CHKERRQ(ierr);
305         ierr = MatSetSizes(G, numProcs, numProcs, numProcs, numProcs);CHKERRQ(ierr);
306         ierr = MatSetUp(G);CHKERRQ(ierr);
307         for (p = 0; p < numProcs; ++p) {
308           for (n = 0; n < recvcounts[p]; ++n) {
309             const PetscInt    r = p;
310             const PetscInt    q = adj[displs[p]+n];
311             const PetscScalar o = val[displs[p]+n] ? 1.0 : 0.0;
312 
313             ierr = MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES);CHKERRQ(ierr);
314             ierr = MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES);CHKERRQ(ierr);
315           }
316         }
317         ierr = MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
318         ierr = MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
319 
320         ierr = PetscBTCreate(numProcs, &seenProcs);CHKERRQ(ierr);
321         ierr = PetscBTMemzero(numProcs, seenProcs);CHKERRQ(ierr);
322         ierr = PetscBTCreate(numProcs, &flippedProcs);CHKERRQ(ierr);
323         ierr = PetscBTMemzero(numProcs, flippedProcs);CHKERRQ(ierr);
324         ierr = PetscMalloc1(numProcs,&procFIFO);CHKERRQ(ierr);
325         pTop = pBottom = 0;
326         for (p = 0; p < numProcs; ++p) {
327           if (PetscBTLookup(seenProcs, p)) continue;
328           /* Initialize FIFO with next proc */
329           procFIFO[pBottom++] = p;
330           ierr = PetscBTSet(seenProcs, p);CHKERRQ(ierr);
331           /* Consider each proc in FIFO */
332           while (pTop < pBottom) {
333             const PetscScalar *ornt;
334             const PetscInt    *neighbors;
335             PetscInt           proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors;
336 
337             proc     = procFIFO[pTop++];
338             flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0;
339             ierr = MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt);CHKERRQ(ierr);
340             /* Loop over neighboring procs */
341             for (n = 0; n < numNeighbors; ++n) {
342               nproc    = neighbors[n];
343               mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1;
344               seen     = PetscBTLookup(seenProcs, nproc);
345               flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0;
346 
347               if (mismatch ^ (flippedA ^ flippedB)) {
348                 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);
349                 if (!flippedB) {
350                   ierr = PetscBTSet(flippedProcs, nproc);CHKERRQ(ierr);
351               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
352               } else if (mismatch && flippedA && flippedB) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
353               if (!seen) {
354                 procFIFO[pBottom++] = nproc;
355                 ierr = PetscBTSet(seenProcs, nproc);CHKERRQ(ierr);
356               }
357             }
358           }
359         }
360         ierr = PetscFree(procFIFO);CHKERRQ(ierr);
361         ierr = MatDestroy(&G);CHKERRQ(ierr);
362 
363         ierr = PetscFree2(recvcounts,displs);CHKERRQ(ierr);
364         ierr = PetscFree2(adj,val);CHKERRQ(ierr);
365         {
366           PetscBool *flips;
367 
368           ierr = PetscMalloc1(numProcs,&flips);CHKERRQ(ierr);
369           for (p = 0; p < numProcs; ++p) {
370             flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE;
371             if (debug && flips[p]) {ierr = PetscPrintf(comm, "Flipping Proc %d:\n", p);}
372           }
373           ierr = MPI_Scatter(flips, 1, MPIU_BOOL, &flipped, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
374           ierr = PetscFree(flips);CHKERRQ(ierr);
375         }
376         ierr = PetscBTDestroy(&seenProcs);CHKERRQ(ierr);
377         ierr = PetscBTDestroy(&flippedProcs);CHKERRQ(ierr);
378       }
379       ierr = PetscFree4(match,nranks,rornt,lornt);CHKERRQ(ierr);
380       ierr = PetscFree(neighbors);CHKERRQ(ierr);
381       if (flipped) {for (c = cStart; c < cEnd; ++c) {ierr = PetscBTNegate(flippedCells, c-cStart);CHKERRQ(ierr);}}
382     }
383   }
384   if (flg) {
385     PetscViewer v;
386     PetscMPIInt rank;
387 
388     ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
389     ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr);
390     ierr = PetscViewerASCIISynchronizedAllow(v, PETSC_TRUE);CHKERRQ(ierr);
391     ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank);CHKERRQ(ierr);
392     ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr);
393   }
394   /* Reverse flipped cells in the mesh */
395   for (c = cStart; c < cEnd; ++c) {
396     if (PetscBTLookup(flippedCells, c-cStart)) {ierr = DMPlexReverseCell(dm, c);CHKERRQ(ierr);}
397   }
398   ierr = PetscBTDestroy(&seenCells);CHKERRQ(ierr);
399   ierr = PetscBTDestroy(&flippedCells);CHKERRQ(ierr);
400   ierr = PetscBTDestroy(&seenFaces);CHKERRQ(ierr);
401   ierr = PetscFree(faceFIFO);CHKERRQ(ierr);
402   PetscFunctionReturn(0);
403 }
404