xref: /petsc/src/dm/impls/plex/plexorient.c (revision 7e09831bc92b7f5435999991a3b282b8ca6d1e95)
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__ "DMPlexOrient"
69 /*@
70   DMPlexOrient - Give a consistent orientation to the input mesh
71 
72   Input Parameters:
73 . dm - The DM
74 
75   Note: The orientation data for the DM are change in-place.
76 $ This routine will fail for non-orientable surfaces, such as the Moebius strip.
77 
78   Level: advanced
79 
80 .seealso: DMCreate(), DMPLEX
81 @*/
82 PetscErrorCode DMPlexOrient(DM dm)
83 {
84   MPI_Comm       comm;
85   PetscBT        seenCells, flippedCells, seenFaces;
86   PetscInt      *faceFIFO, fTop, fBottom;
87   PetscInt       dim, h, cStart, cEnd, c, fStart, fEnd, face;
88   PetscBool      flg;
89   PetscErrorCode ierr;
90 
91   PetscFunctionBegin;
92   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
93   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-orientation_view", &flg);CHKERRQ(ierr);
94   /* Truth Table
95      mismatch    flips   do action   mismatch   flipA ^ flipB   action
96          F       0 flips     no         F             F           F
97          F       1 flip      yes        F             T           T
98          F       2 flips     no         T             F           T
99          T       0 flips     yes        T             T           F
100          T       1 flip      no
101          T       2 flips     yes
102   */
103   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
104   ierr = DMPlexGetVTKCellHeight(dm, &h);CHKERRQ(ierr);
105   ierr = DMPlexGetHeightStratum(dm, h,   &cStart, &cEnd);CHKERRQ(ierr);
106   ierr = DMPlexGetHeightStratum(dm, h+1, &fStart, &fEnd);CHKERRQ(ierr);
107   ierr = PetscBTCreate(cEnd - cStart, &seenCells);CHKERRQ(ierr);
108   ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr);
109   ierr = PetscBTCreate(cEnd - cStart, &flippedCells);CHKERRQ(ierr);
110   ierr = PetscBTMemzero(cEnd - cStart, flippedCells);CHKERRQ(ierr);
111   ierr = PetscBTCreate(fEnd - fStart, &seenFaces);CHKERRQ(ierr);
112   ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr);
113   ierr = PetscMalloc1((fEnd - fStart), &faceFIFO);CHKERRQ(ierr);
114   fTop = fBottom = 0;
115   /* Initialize FIFO with first cell */
116   if (cEnd > cStart) {
117     const PetscInt *cone;
118     PetscInt        coneSize;
119 
120     ierr = DMPlexGetConeSize(dm, cStart, &coneSize);CHKERRQ(ierr);
121     ierr = DMPlexGetCone(dm, cStart, &cone);CHKERRQ(ierr);
122     for (c = 0; c < coneSize; ++c) {
123       faceFIFO[fBottom++] = cone[c];
124       ierr = PetscBTSet(seenFaces, cone[c]-fStart);CHKERRQ(ierr);
125     }
126   }
127   /* Consider each face in FIFO */
128   while (fTop < fBottom) {
129     const PetscInt *support, *coneA, *coneB, *coneOA, *coneOB;
130     PetscInt        supportSize, coneSizeA, coneSizeB, posA = -1, posB = -1;
131     PetscInt        seenA, flippedA, seenB, flippedB, mismatch;
132 
133     face = faceFIFO[fTop++];
134     ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
135     ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr);
136     if (supportSize < 2) continue;
137     if (supportSize != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %d", supportSize);
138     seenA    = PetscBTLookup(seenCells,    support[0]-cStart);
139     flippedA = PetscBTLookup(flippedCells, support[0]-cStart) ? 1 : 0;
140     seenB    = PetscBTLookup(seenCells,    support[1]-cStart);
141     flippedB = PetscBTLookup(flippedCells, support[1]-cStart) ? 1 : 0;
142 
143     ierr = DMPlexGetConeSize(dm, support[0], &coneSizeA);CHKERRQ(ierr);
144     ierr = DMPlexGetConeSize(dm, support[1], &coneSizeB);CHKERRQ(ierr);
145     ierr = DMPlexGetCone(dm, support[0], &coneA);CHKERRQ(ierr);
146     ierr = DMPlexGetCone(dm, support[1], &coneB);CHKERRQ(ierr);
147     ierr = DMPlexGetConeOrientation(dm, support[0], &coneOA);CHKERRQ(ierr);
148     ierr = DMPlexGetConeOrientation(dm, support[1], &coneOB);CHKERRQ(ierr);
149     for (c = 0; c < coneSizeA; ++c) {
150       if (!PetscBTLookup(seenFaces, coneA[c]-fStart)) {
151         faceFIFO[fBottom++] = coneA[c];
152         ierr = PetscBTSet(seenFaces, coneA[c]-fStart);CHKERRQ(ierr);
153       }
154       if (coneA[c] == face) posA = c;
155       if (fBottom > fEnd-fStart) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %d was pushed exceeding capacity %d > %d", coneA[c], fBottom, fEnd-fStart);
156     }
157     if (posA < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d could not be located in cell %d", face, support[0]);
158     for (c = 0; c < coneSizeB; ++c) {
159       if (!PetscBTLookup(seenFaces, coneB[c]-fStart)) {
160         faceFIFO[fBottom++] = coneB[c];
161         ierr = PetscBTSet(seenFaces, coneB[c]-fStart);CHKERRQ(ierr);
162       }
163       if (coneB[c] == face) posB = c;
164       if (fBottom > fEnd-fStart) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %d was pushed exceeding capacity %d > %d", coneA[c], fBottom, fEnd-fStart);
165     }
166     if (posB < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d could not be located in cell %d", face, support[1]);
167 
168     if (dim == 1) {
169       mismatch = posA == posB;
170     } else {
171       mismatch = coneOA[posA] == coneOB[posB];
172     }
173 
174     if (mismatch ^ (flippedA ^ flippedB)) {
175       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]);
176       if (!seenA && !flippedA) {
177         ierr = PetscBTSet(flippedCells, support[0]-cStart);CHKERRQ(ierr);
178       } else if (!seenB && !flippedB) {
179         ierr = PetscBTSet(flippedCells, support[1]-cStart);CHKERRQ(ierr);
180       } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
181     } else if (mismatch && flippedA && flippedB) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
182     ierr = PetscBTSet(seenCells, support[0]-cStart);CHKERRQ(ierr);
183     ierr = PetscBTSet(seenCells, support[1]-cStart);CHKERRQ(ierr);
184   }
185   if (flg) {
186     PetscViewer v;
187     PetscMPIInt rank;
188 
189     ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
190     ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr);
191     ierr = PetscViewerASCIISynchronizedAllow(v, PETSC_TRUE);CHKERRQ(ierr);
192     ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial seen faces:\n", rank);CHKERRQ(ierr);
193     ierr = PetscBTView(fEnd-fStart, seenFaces,    v);CHKERRQ(ierr);
194     ierr = PetscViewerASCIISynchronizedAllow(v, PETSC_TRUE);CHKERRQ(ierr);
195     ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank);CHKERRQ(ierr);
196     ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr);
197   }
198   /* Now all subdomains are oriented, but we need a consistent parallel orientation */
199   {
200     /* Find a representative face (edge) separating pairs of procs */
201     PetscSF            sf;
202     const PetscInt    *lpoints;
203     const PetscSFNode *rpoints;
204     PetscInt          *neighbors, *nranks;
205     PetscInt           numLeaves, numRoots, numNeighbors = 0, l, n;
206 
207     ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
208     ierr = PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints);CHKERRQ(ierr);
209     if (numLeaves >= 0) {
210       const PetscInt *cone, *ornt, *support;
211       PetscInt        coneSize, supportSize;
212       int            *rornt, *lornt; /* PetscSF cannot handle smaller than int */
213       PetscBool      *match, flipped = PETSC_FALSE;
214 
215       ierr = PetscMalloc1(numLeaves,&neighbors);CHKERRQ(ierr);
216       /* I know this is p^2 time in general, but for bounded degree its alright */
217       for (l = 0; l < numLeaves; ++l) {
218         const PetscInt face = lpoints[l];
219         if ((face >= fStart) && (face < fEnd)) {
220           const PetscInt rank = rpoints[l].rank;
221           for (n = 0; n < numNeighbors; ++n) if (rank == rpoints[neighbors[n]].rank) break;
222           if (n >= numNeighbors) {
223             PetscInt supportSize;
224             ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
225             if (supportSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %d", supportSize);
226             neighbors[numNeighbors++] = l;
227           }
228         }
229       }
230       ierr = PetscCalloc4(numNeighbors,&match,numNeighbors,&nranks,numRoots,&rornt,numRoots,&lornt);CHKERRQ(ierr);
231       for (face = fStart; face < fEnd; ++face) {
232         ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
233         if (supportSize != 1) continue;
234         ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr);
235 
236         ierr = DMPlexGetCone(dm, support[0], &cone);CHKERRQ(ierr);
237         ierr = DMPlexGetConeSize(dm, support[0], &coneSize);CHKERRQ(ierr);
238         ierr = DMPlexGetConeOrientation(dm, support[0], &ornt);CHKERRQ(ierr);
239         for (c = 0; c < coneSize; ++c) if (cone[c] == face) break;
240         if (dim == 1) {
241           /* Use cone position instead, shifted to -1 or 1 */
242           rornt[face] = c*2-1;
243         } else {
244           if (PetscBTLookup(flippedCells, support[0]-cStart)) rornt[face] = ornt[c] < 0 ? -1 :  1;
245           else                                                rornt[face] = ornt[c] < 0 ?  1 : -1;
246         }
247       }
248       /* Mark each edge with match or nomatch */
249       ierr = PetscSFBcastBegin(sf, MPI_INT, rornt, lornt);CHKERRQ(ierr);
250       ierr = PetscSFBcastEnd(sf, MPI_INT, rornt, lornt);CHKERRQ(ierr);
251       for (n = 0; n < numNeighbors; ++n) {
252         const PetscInt face = lpoints[neighbors[n]];
253 
254         if (rornt[face]*lornt[face] < 0) match[n] = PETSC_TRUE;
255         else                             match[n] = PETSC_FALSE;
256         nranks[n] = rpoints[neighbors[n]].rank;
257       }
258       /* Collect the graph on 0 */
259       {
260         Mat          G;
261         PetscBT      seenProcs, flippedProcs;
262         PetscInt    *procFIFO, pTop, pBottom;
263         PetscInt    *adj = NULL;
264         PetscBool   *val = NULL;
265         PetscMPIInt *recvcounts = NULL, *displs = NULL, p;
266         PetscMPIInt  N = numNeighbors, numProcs = 0, rank;
267         PetscInt     debug = 0;
268 
269         ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
270         if (!rank) {ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);}
271         ierr = PetscCalloc2(numProcs,&recvcounts,numProcs+1,&displs);CHKERRQ(ierr);
272         ierr = MPI_Gather(&N, 1, MPI_INT, recvcounts, 1, MPI_INT, 0, comm);CHKERRQ(ierr);
273         for (p = 0; p < numProcs; ++p) {
274           displs[p+1] = displs[p] + recvcounts[p];
275         }
276         if (!rank) {ierr = PetscMalloc2(displs[numProcs],&adj,displs[numProcs],&val);CHKERRQ(ierr);}
277         ierr = MPI_Gatherv(nranks, numNeighbors, MPIU_INT, adj, recvcounts, displs, MPIU_INT, 0, comm);CHKERRQ(ierr);
278         ierr = MPI_Gatherv(match, numNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
279         if (debug) {
280           for (p = 0; p < numProcs; ++p) {
281             ierr = PetscPrintf(comm, "Proc %d:\n", p);
282             for (n = 0; n < recvcounts[p]; ++n) {
283               ierr = PetscPrintf(comm, "  edge %d (%d):\n", adj[displs[p]+n], val[displs[p]+n]);
284             }
285           }
286         }
287         /* Symmetrize the graph */
288         ierr = MatCreate(PETSC_COMM_SELF, &G);CHKERRQ(ierr);
289         ierr = MatSetSizes(G, numProcs, numProcs, numProcs, numProcs);CHKERRQ(ierr);
290         ierr = MatSetUp(G);CHKERRQ(ierr);
291         for (p = 0; p < numProcs; ++p) {
292           for (n = 0; n < recvcounts[p]; ++n) {
293             const PetscInt    r = p;
294             const PetscInt    q = adj[displs[p]+n];
295             const PetscScalar o = val[displs[p]+n] ? 1.0 : 0.0;
296 
297             ierr = MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES);CHKERRQ(ierr);
298             ierr = MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES);CHKERRQ(ierr);
299           }
300         }
301         ierr = MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
302         ierr = MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
303 
304         ierr = PetscBTCreate(numProcs, &seenProcs);CHKERRQ(ierr);
305         ierr = PetscBTMemzero(numProcs, seenProcs);CHKERRQ(ierr);
306         ierr = PetscBTCreate(numProcs, &flippedProcs);CHKERRQ(ierr);
307         ierr = PetscBTMemzero(numProcs, flippedProcs);CHKERRQ(ierr);
308         ierr = PetscMalloc1(numProcs,&procFIFO);CHKERRQ(ierr);
309         pTop = pBottom = 0;
310         for (p = 0; p < numProcs; ++p) {
311           if (PetscBTLookup(seenProcs, p)) continue;
312           /* Initialize FIFO with next proc */
313           procFIFO[pBottom++] = p;
314           ierr = PetscBTSet(seenProcs, p);CHKERRQ(ierr);
315           /* Consider each proc in FIFO */
316           while (pTop < pBottom) {
317             const PetscScalar *ornt;
318             const PetscInt    *neighbors;
319             PetscInt           proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors;
320 
321             proc     = procFIFO[pTop++];
322             flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0;
323             ierr = MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt);CHKERRQ(ierr);
324             /* Loop over neighboring procs */
325             for (n = 0; n < numNeighbors; ++n) {
326               nproc    = neighbors[n];
327               mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1;
328               seen     = PetscBTLookup(seenProcs, nproc);
329               flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0;
330 
331               if (mismatch ^ (flippedA ^ flippedB)) {
332                 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);
333                 if (!flippedB) {
334                   ierr = PetscBTSet(flippedProcs, nproc);CHKERRQ(ierr);
335               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
336               } else if (mismatch && flippedA && flippedB) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
337               if (!seen) {
338                 procFIFO[pBottom++] = nproc;
339                 ierr = PetscBTSet(seenProcs, nproc);CHKERRQ(ierr);
340               }
341             }
342           }
343         }
344         ierr = PetscFree(procFIFO);CHKERRQ(ierr);
345         ierr = MatDestroy(&G);CHKERRQ(ierr);
346 
347         ierr = PetscFree2(recvcounts,displs);CHKERRQ(ierr);
348         ierr = PetscFree2(adj,val);CHKERRQ(ierr);
349         {
350           PetscBool *flips;
351 
352           ierr = PetscMalloc1(numProcs,&flips);CHKERRQ(ierr);
353           for (p = 0; p < numProcs; ++p) {
354             flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE;
355             if (debug && flips[p]) {ierr = PetscPrintf(comm, "Flipping Proc %d:\n", p);}
356           }
357           ierr = MPI_Scatter(flips, 1, MPIU_BOOL, &flipped, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
358           ierr = PetscFree(flips);CHKERRQ(ierr);
359         }
360         ierr = PetscBTDestroy(&seenProcs);CHKERRQ(ierr);
361         ierr = PetscBTDestroy(&flippedProcs);CHKERRQ(ierr);
362       }
363       ierr = PetscFree4(match,nranks,rornt,lornt);CHKERRQ(ierr);
364       ierr = PetscFree(neighbors);CHKERRQ(ierr);
365       if (flipped) {for (c = cStart; c < cEnd; ++c) {ierr = PetscBTNegate(flippedCells, c-cStart);CHKERRQ(ierr);}}
366     }
367   }
368   if (flg) {
369     PetscViewer v;
370     PetscMPIInt rank;
371 
372     ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
373     ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr);
374     ierr = PetscViewerASCIISynchronizedAllow(v, PETSC_TRUE);CHKERRQ(ierr);
375     ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank);CHKERRQ(ierr);
376     ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr);
377   }
378   /* Reverse flipped cells in the mesh */
379   for (c = cStart; c < cEnd; ++c) {
380     if (PetscBTLookup(flippedCells, c-cStart)) {ierr = DMPlexReverseCell(dm, c);CHKERRQ(ierr);}
381   }
382   ierr = PetscBTDestroy(&seenCells);CHKERRQ(ierr);
383   ierr = PetscBTDestroy(&flippedCells);CHKERRQ(ierr);
384   ierr = PetscBTDestroy(&seenFaces);CHKERRQ(ierr);
385   ierr = PetscFree(faceFIFO);CHKERRQ(ierr);
386   PetscFunctionReturn(0);
387 }
388