xref: /petsc/src/dm/impls/plex/plexorient.c (revision e1d831098ca2bcddfd7593004fb80ef665cdae95)
184961fc4SMatthew G. Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
284961fc4SMatthew G. Knepley #include <petscsf.h>
384961fc4SMatthew G. Knepley 
484961fc4SMatthew G. Knepley #undef __FUNCT__
584961fc4SMatthew G. Knepley #define __FUNCT__ "DMPlexReverseCell"
684961fc4SMatthew G. Knepley /*@
784961fc4SMatthew G. Knepley   DMPlexReverseCell - Give a mesh cell the opposite orientation
884961fc4SMatthew G. Knepley 
984961fc4SMatthew G. Knepley   Input Parameters:
1084961fc4SMatthew G. Knepley + dm   - The DM
1184961fc4SMatthew G. Knepley - cell - The cell number
1284961fc4SMatthew G. Knepley 
1384961fc4SMatthew G. Knepley   Note: The modification of the DM is done in-place.
1484961fc4SMatthew G. Knepley 
1584961fc4SMatthew G. Knepley   Level: advanced
1684961fc4SMatthew G. Knepley 
1784961fc4SMatthew G. Knepley .seealso: DMPlexOrient(), DMCreate(), DMPLEX
1884961fc4SMatthew G. Knepley @*/
1984961fc4SMatthew G. Knepley PetscErrorCode DMPlexReverseCell(DM dm, PetscInt cell)
2084961fc4SMatthew G. Knepley {
2184961fc4SMatthew G. Knepley   /* Note that the reverse orientation ro of a face with orientation o is:
2284961fc4SMatthew G. Knepley 
2384961fc4SMatthew G. Knepley        ro = o >= 0 ? -(faceSize - o) : faceSize + o
2484961fc4SMatthew G. Knepley 
2584961fc4SMatthew G. Knepley      where faceSize is the size of the cone for the face.
2684961fc4SMatthew G. Knepley   */
2784961fc4SMatthew G. Knepley   const PetscInt *cone,    *coneO, *support;
2884961fc4SMatthew G. Knepley   PetscInt       *revcone, *revconeO;
2984961fc4SMatthew G. Knepley   PetscInt        maxConeSize, coneSize, supportSize, faceSize, cp, sp;
3084961fc4SMatthew G. Knepley   PetscErrorCode  ierr;
3184961fc4SMatthew G. Knepley 
3284961fc4SMatthew G. Knepley   PetscFunctionBegin;
3384961fc4SMatthew G. Knepley   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
3484961fc4SMatthew G. Knepley   ierr = DMGetWorkArray(dm, maxConeSize, PETSC_INT, &revcone);CHKERRQ(ierr);
3584961fc4SMatthew G. Knepley   ierr = DMGetWorkArray(dm, maxConeSize, PETSC_INT, &revconeO);CHKERRQ(ierr);
3684961fc4SMatthew G. Knepley   /* Reverse cone, and reverse orientations of faces */
3784961fc4SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
3884961fc4SMatthew G. Knepley   ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
3984961fc4SMatthew G. Knepley   ierr = DMPlexGetConeOrientation(dm, cell, &coneO);CHKERRQ(ierr);
4084961fc4SMatthew G. Knepley   for (cp = 0; cp < coneSize; ++cp) {
4184961fc4SMatthew G. Knepley     const PetscInt rcp = coneSize-cp-1;
4284961fc4SMatthew G. Knepley 
4384961fc4SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, cone[rcp], &faceSize);CHKERRQ(ierr);
4484961fc4SMatthew G. Knepley     revcone[cp]  = cone[rcp];
4584961fc4SMatthew G. Knepley     revconeO[cp] = coneO[rcp] >= 0 ? -(faceSize-coneO[rcp]) : faceSize+coneO[rcp];
4684961fc4SMatthew G. Knepley   }
4784961fc4SMatthew G. Knepley   ierr = DMPlexSetCone(dm, cell, revcone);CHKERRQ(ierr);
4884961fc4SMatthew G. Knepley   ierr = DMPlexSetConeOrientation(dm, cell, revconeO);CHKERRQ(ierr);
4984961fc4SMatthew G. Knepley   /* Reverse orientation of this cell in the support hypercells */
5084961fc4SMatthew G. Knepley   faceSize = coneSize;
5184961fc4SMatthew G. Knepley   ierr = DMPlexGetSupportSize(dm, cell, &supportSize);CHKERRQ(ierr);
5284961fc4SMatthew G. Knepley   ierr = DMPlexGetSupport(dm, cell, &support);CHKERRQ(ierr);
5384961fc4SMatthew G. Knepley   for (sp = 0; sp < supportSize; ++sp) {
5484961fc4SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, support[sp], &coneSize);CHKERRQ(ierr);
5584961fc4SMatthew G. Knepley     ierr = DMPlexGetCone(dm, support[sp], &cone);CHKERRQ(ierr);
5684961fc4SMatthew G. Knepley     ierr = DMPlexGetConeOrientation(dm, support[sp], &coneO);CHKERRQ(ierr);
5784961fc4SMatthew G. Knepley     for (cp = 0; cp < coneSize; ++cp) {
5884961fc4SMatthew G. Knepley       if (cone[cp] != cell) continue;
5984961fc4SMatthew G. Knepley       ierr = DMPlexInsertConeOrientation(dm, support[sp], cp, coneO[cp] >= 0 ? -(faceSize-coneO[cp]) : faceSize+coneO[cp]);CHKERRQ(ierr);
6084961fc4SMatthew G. Knepley     }
6184961fc4SMatthew G. Knepley   }
6284961fc4SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, maxConeSize, PETSC_INT, &revcone);CHKERRQ(ierr);
6384961fc4SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, maxConeSize, PETSC_INT, &revconeO);CHKERRQ(ierr);
6484961fc4SMatthew G. Knepley   PetscFunctionReturn(0);
6584961fc4SMatthew G. Knepley }
6684961fc4SMatthew G. Knepley 
6784961fc4SMatthew G. Knepley #undef __FUNCT__
687b310ce2SMatthew G. Knepley #define __FUNCT__ "DMPlexCheckFace_Internal"
697b310ce2SMatthew G. Knepley /*
707b310ce2SMatthew G. Knepley   - Checks face match
717b310ce2SMatthew G. Knepley     - Flips non-matching
727b310ce2SMatthew G. Knepley   - Inserts faces of support cells in FIFO
737b310ce2SMatthew G. Knepley */
747b310ce2SMatthew G. Knepley 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)
757b310ce2SMatthew G. Knepley {
767b310ce2SMatthew G. Knepley   const PetscInt *support, *coneA, *coneB, *coneOA, *coneOB;
777b310ce2SMatthew G. Knepley   PetscInt        supportSize, coneSizeA, coneSizeB, posA = -1, posB = -1;
787b310ce2SMatthew G. Knepley   PetscInt        dim, seenA, flippedA, seenB, flippedB, mismatch, c;
797b310ce2SMatthew G. Knepley   PetscErrorCode  ierr;
807b310ce2SMatthew G. Knepley 
817b310ce2SMatthew G. Knepley   PetscFunctionBegin;
827b310ce2SMatthew G. Knepley   face = faceFIFO[(*fTop)++];
837b310ce2SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
847b310ce2SMatthew G. Knepley   ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
857b310ce2SMatthew G. Knepley   ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr);
867b310ce2SMatthew G. Knepley   if (supportSize < 2) PetscFunctionReturn(0);
877b310ce2SMatthew G. Knepley   if (supportSize != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %d", supportSize);
887b310ce2SMatthew G. Knepley   seenA    = PetscBTLookup(seenCells,    support[0]-cStart);
897b310ce2SMatthew G. Knepley   flippedA = PetscBTLookup(flippedCells, support[0]-cStart) ? 1 : 0;
907b310ce2SMatthew G. Knepley   seenB    = PetscBTLookup(seenCells,    support[1]-cStart);
917b310ce2SMatthew G. Knepley   flippedB = PetscBTLookup(flippedCells, support[1]-cStart) ? 1 : 0;
927b310ce2SMatthew G. Knepley 
937b310ce2SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, support[0], &coneSizeA);CHKERRQ(ierr);
947b310ce2SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, support[1], &coneSizeB);CHKERRQ(ierr);
957b310ce2SMatthew G. Knepley   ierr = DMPlexGetCone(dm, support[0], &coneA);CHKERRQ(ierr);
967b310ce2SMatthew G. Knepley   ierr = DMPlexGetCone(dm, support[1], &coneB);CHKERRQ(ierr);
977b310ce2SMatthew G. Knepley   ierr = DMPlexGetConeOrientation(dm, support[0], &coneOA);CHKERRQ(ierr);
987b310ce2SMatthew G. Knepley   ierr = DMPlexGetConeOrientation(dm, support[1], &coneOB);CHKERRQ(ierr);
997b310ce2SMatthew G. Knepley   for (c = 0; c < coneSizeA; ++c) {
1007b310ce2SMatthew G. Knepley     if (!PetscBTLookup(seenFaces, coneA[c]-fStart)) {
1017b310ce2SMatthew G. Knepley       faceFIFO[(*fBottom)++] = coneA[c];
1027b310ce2SMatthew G. Knepley       ierr = PetscBTSet(seenFaces, coneA[c]-fStart);CHKERRQ(ierr);
1037b310ce2SMatthew G. Knepley     }
1047b310ce2SMatthew G. Knepley     if (coneA[c] == face) posA = c;
1057b310ce2SMatthew G. Knepley     if (*fBottom > fEnd-fStart) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %d was pushed exceeding capacity %d > %d", coneA[c], *fBottom, fEnd-fStart);
1067b310ce2SMatthew G. Knepley   }
1077b310ce2SMatthew G. Knepley   if (posA < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d could not be located in cell %d", face, support[0]);
1087b310ce2SMatthew G. Knepley   for (c = 0; c < coneSizeB; ++c) {
1097b310ce2SMatthew G. Knepley     if (!PetscBTLookup(seenFaces, coneB[c]-fStart)) {
1107b310ce2SMatthew G. Knepley       faceFIFO[(*fBottom)++] = coneB[c];
1117b310ce2SMatthew G. Knepley       ierr = PetscBTSet(seenFaces, coneB[c]-fStart);CHKERRQ(ierr);
1127b310ce2SMatthew G. Knepley     }
1137b310ce2SMatthew G. Knepley     if (coneB[c] == face) posB = c;
1147b310ce2SMatthew G. Knepley     if (*fBottom > fEnd-fStart) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %d was pushed exceeding capacity %d > %d", coneA[c], *fBottom, fEnd-fStart);
1157b310ce2SMatthew G. Knepley   }
1167b310ce2SMatthew G. Knepley   if (posB < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d could not be located in cell %d", face, support[1]);
1177b310ce2SMatthew G. Knepley 
1187b310ce2SMatthew G. Knepley   if (dim == 1) {
1197b310ce2SMatthew G. Knepley     mismatch = posA == posB;
1207b310ce2SMatthew G. Knepley   } else {
1217b310ce2SMatthew G. Knepley     mismatch = coneOA[posA] == coneOB[posB];
1227b310ce2SMatthew G. Knepley   }
1237b310ce2SMatthew G. Knepley 
1247b310ce2SMatthew G. Knepley   if (mismatch ^ (flippedA ^ flippedB)) {
1257b310ce2SMatthew G. Knepley     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]);
1267b310ce2SMatthew G. Knepley     if (!seenA && !flippedA) {
1277b310ce2SMatthew G. Knepley       ierr = PetscBTSet(flippedCells, support[0]-cStart);CHKERRQ(ierr);
1287b310ce2SMatthew G. Knepley     } else if (!seenB && !flippedB) {
1297b310ce2SMatthew G. Knepley       ierr = PetscBTSet(flippedCells, support[1]-cStart);CHKERRQ(ierr);
1307b310ce2SMatthew G. Knepley     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
1317b310ce2SMatthew G. Knepley   } else if (mismatch && flippedA && flippedB) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
1327b310ce2SMatthew G. Knepley   ierr = PetscBTSet(seenCells, support[0]-cStart);CHKERRQ(ierr);
1337b310ce2SMatthew G. Knepley   ierr = PetscBTSet(seenCells, support[1]-cStart);CHKERRQ(ierr);
1347b310ce2SMatthew G. Knepley   PetscFunctionReturn(0);
1357b310ce2SMatthew G. Knepley }
1367b310ce2SMatthew G. Knepley 
1377b310ce2SMatthew G. Knepley #undef __FUNCT__
13884961fc4SMatthew G. Knepley #define __FUNCT__ "DMPlexOrient"
13984961fc4SMatthew G. Knepley /*@
14084961fc4SMatthew G. Knepley   DMPlexOrient - Give a consistent orientation to the input mesh
14184961fc4SMatthew G. Knepley 
14284961fc4SMatthew G. Knepley   Input Parameters:
14384961fc4SMatthew G. Knepley . dm - The DM
14484961fc4SMatthew G. Knepley 
14584961fc4SMatthew G. Knepley   Note: The orientation data for the DM are change in-place.
14684961fc4SMatthew G. Knepley $ This routine will fail for non-orientable surfaces, such as the Moebius strip.
14784961fc4SMatthew G. Knepley 
14884961fc4SMatthew G. Knepley   Level: advanced
14984961fc4SMatthew G. Knepley 
15084961fc4SMatthew G. Knepley .seealso: DMCreate(), DMPLEX
15184961fc4SMatthew G. Knepley @*/
15284961fc4SMatthew G. Knepley PetscErrorCode DMPlexOrient(DM dm)
15384961fc4SMatthew G. Knepley {
15484961fc4SMatthew G. Knepley   MPI_Comm           comm;
155*e1d83109SMatthew G. Knepley   PetscSF            sf;
156*e1d83109SMatthew G. Knepley   const PetscInt    *lpoints;
157*e1d83109SMatthew G. Knepley   const PetscSFNode *rpoints;
158*e1d83109SMatthew G. Knepley   PetscSFNode       *rorntComp = NULL, *lorntComp = NULL;
159*e1d83109SMatthew G. Knepley   PetscInt          *numNeighbors, **neighbors;
160*e1d83109SMatthew G. Knepley   PetscSFNode       *nrankComp;
161*e1d83109SMatthew G. Knepley   PetscBool         *match, *flipped;
16284961fc4SMatthew G. Knepley   PetscBT            seenCells, flippedCells, seenFaces;
163*e1d83109SMatthew G. Knepley   PetscInt          *faceFIFO, fTop, fBottom, *cellComp, *faceComp;
164*e1d83109SMatthew G. Knepley   PetscInt           numLeaves, numRoots, dim, h, cStart, cEnd, c, cell, fStart, fEnd, face, numComponents, off, totNeighbors = 0, comp = 0;
16584961fc4SMatthew G. Knepley   PetscBool          flg;
16684961fc4SMatthew G. Knepley   PetscErrorCode     ierr;
16784961fc4SMatthew G. Knepley 
16884961fc4SMatthew G. Knepley   PetscFunctionBegin;
16984961fc4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
17084961fc4SMatthew G. Knepley   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-orientation_view", &flg);CHKERRQ(ierr);
171*e1d83109SMatthew G. Knepley   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
172*e1d83109SMatthew G. Knepley   ierr = PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints);CHKERRQ(ierr);
17384961fc4SMatthew G. Knepley   /* Truth Table
17484961fc4SMatthew G. Knepley      mismatch    flips   do action   mismatch   flipA ^ flipB   action
17584961fc4SMatthew G. Knepley          F       0 flips     no         F             F           F
17684961fc4SMatthew G. Knepley          F       1 flip      yes        F             T           T
17784961fc4SMatthew G. Knepley          F       2 flips     no         T             F           T
17884961fc4SMatthew G. Knepley          T       0 flips     yes        T             T           F
17984961fc4SMatthew G. Knepley          T       1 flip      no
18084961fc4SMatthew G. Knepley          T       2 flips     yes
18184961fc4SMatthew G. Knepley   */
18284961fc4SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
18384961fc4SMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(dm, &h);CHKERRQ(ierr);
18484961fc4SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, h,   &cStart, &cEnd);CHKERRQ(ierr);
18584961fc4SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, h+1, &fStart, &fEnd);CHKERRQ(ierr);
18684961fc4SMatthew G. Knepley   ierr = PetscBTCreate(cEnd - cStart, &seenCells);CHKERRQ(ierr);
18784961fc4SMatthew G. Knepley   ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr);
18884961fc4SMatthew G. Knepley   ierr = PetscBTCreate(cEnd - cStart, &flippedCells);CHKERRQ(ierr);
18984961fc4SMatthew G. Knepley   ierr = PetscBTMemzero(cEnd - cStart, flippedCells);CHKERRQ(ierr);
19084961fc4SMatthew G. Knepley   ierr = PetscBTCreate(fEnd - fStart, &seenFaces);CHKERRQ(ierr);
19184961fc4SMatthew G. Knepley   ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr);
192*e1d83109SMatthew G. Knepley   ierr = PetscCalloc3(fEnd - fStart, &faceFIFO, cEnd-cStart, &cellComp, fEnd-fStart, &faceComp);CHKERRQ(ierr);
193*e1d83109SMatthew G. Knepley   /*
194*e1d83109SMatthew G. Knepley    OLD STYLE
195*e1d83109SMatthew G. Knepley    - Add an integer array over cells and faces (component) for connected component number
196*e1d83109SMatthew G. Knepley    Foreach component
197*e1d83109SMatthew G. Knepley      - Mark the initial cell as seen
198*e1d83109SMatthew G. Knepley      - Process component as usual
199*e1d83109SMatthew G. Knepley      - Set component for all seenCells
200*e1d83109SMatthew G. Knepley      - Wipe seenCells and seenFaces (flippedCells can stay)
201*e1d83109SMatthew G. Knepley    - Generate parallel adjacency for component using SF and seenFaces
202*e1d83109SMatthew G. Knepley    - Collect numComponents adj data from each proc to 0
203*e1d83109SMatthew G. Knepley    - Build same serial graph
204*e1d83109SMatthew G. Knepley    - Use same solver
205*e1d83109SMatthew G. Knepley    - Use Scatterv to to send back flipped flags for each component
206*e1d83109SMatthew G. Knepley    - Negate flippedCells by component
207*e1d83109SMatthew G. Knepley 
208*e1d83109SMatthew G. Knepley    NEW STYLE
209*e1d83109SMatthew G. Knepley    - Create the adj on each process
210*e1d83109SMatthew G. Knepley    - Bootstrap to complete graph on proc 0
211*e1d83109SMatthew G. Knepley   */
212*e1d83109SMatthew G. Knepley   /* Loop over components */
213*e1d83109SMatthew G. Knepley   for (cell = cStart; cell < cEnd; ++cell) cellComp[cell-cStart] = -1;
214*e1d83109SMatthew G. Knepley   do {
215*e1d83109SMatthew G. Knepley     /* Look for first unmarked cell */
216*e1d83109SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) if (cellComp[cell-cStart] < 0) break;
217*e1d83109SMatthew G. Knepley     if (cell >= cEnd) break;
218*e1d83109SMatthew G. Knepley     /* Initialize FIFO with first cell in component */
219*e1d83109SMatthew G. Knepley     {
22084961fc4SMatthew G. Knepley       const PetscInt *cone;
22184961fc4SMatthew G. Knepley       PetscInt        coneSize;
22284961fc4SMatthew G. Knepley 
223*e1d83109SMatthew G. Knepley       fTop = fBottom = 0;
224*e1d83109SMatthew G. Knepley       ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
225*e1d83109SMatthew G. Knepley       ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
22684961fc4SMatthew G. Knepley       for (c = 0; c < coneSize; ++c) {
22784961fc4SMatthew G. Knepley         faceFIFO[fBottom++] = cone[c];
22884961fc4SMatthew G. Knepley         ierr = PetscBTSet(seenFaces, cone[c]-fStart);CHKERRQ(ierr);
22984961fc4SMatthew G. Knepley       }
230*e1d83109SMatthew G. Knepley       ierr = PetscBTSet(seenCells, cell);CHKERRQ(ierr);
23184961fc4SMatthew G. Knepley     }
23284961fc4SMatthew G. Knepley     /* Consider each face in FIFO */
23384961fc4SMatthew G. Knepley     while (fTop < fBottom) {
2347b310ce2SMatthew G. Knepley       ierr = DMPlexCheckFace_Internal(dm, face, faceFIFO, &fTop, &fBottom, cStart, fStart, fEnd, seenCells, flippedCells, seenFaces);CHKERRQ(ierr);
23584961fc4SMatthew G. Knepley     }
236*e1d83109SMatthew G. Knepley     /* Set component for cells and faces */
237*e1d83109SMatthew G. Knepley     for (cell = 0; cell < cEnd-cStart; ++cell) {
238*e1d83109SMatthew G. Knepley       if (PetscBTLookup(seenCells, cell)) cellComp[cell] = comp;
239*e1d83109SMatthew G. Knepley     }
240*e1d83109SMatthew G. Knepley     for (face = 0; face < fEnd-fStart; ++face) {
241*e1d83109SMatthew G. Knepley       if (PetscBTLookup(seenFaces, face)) faceComp[face] = comp;
242*e1d83109SMatthew G. Knepley     }
243*e1d83109SMatthew G. Knepley     /* Wipe seenCells and seenFaces for next component */
244*e1d83109SMatthew G. Knepley     ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr);
245*e1d83109SMatthew G. Knepley     ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr);
246*e1d83109SMatthew G. Knepley     ++comp;
247*e1d83109SMatthew G. Knepley   } while (1);
248*e1d83109SMatthew G. Knepley   numComponents = comp;
24984961fc4SMatthew G. Knepley   if (flg) {
25084961fc4SMatthew G. Knepley     PetscViewer v;
25184961fc4SMatthew G. Knepley     PetscMPIInt rank;
25284961fc4SMatthew G. Knepley 
25384961fc4SMatthew G. Knepley     ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
2547e09831bSMatthew G. Knepley     ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr);
2557e09831bSMatthew G. Knepley     ierr = PetscViewerASCIISynchronizedAllow(v, PETSC_TRUE);CHKERRQ(ierr);
25684961fc4SMatthew G. Knepley     ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank);CHKERRQ(ierr);
25784961fc4SMatthew G. Knepley     ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr);
25884961fc4SMatthew G. Knepley   }
25984961fc4SMatthew G. Knepley   /* Now all subdomains are oriented, but we need a consistent parallel orientation */
26084961fc4SMatthew G. Knepley   if (numLeaves >= 0) {
261*e1d83109SMatthew G. Knepley     /* Store orientations of boundary faces*/
262*e1d83109SMatthew G. Knepley     ierr = PetscCalloc2(numRoots,&rorntComp,numRoots,&lorntComp);CHKERRQ(ierr);
26384961fc4SMatthew G. Knepley     for (face = fStart; face < fEnd; ++face) {
264*e1d83109SMatthew G. Knepley       const PetscInt *cone, *support, *ornt;
265*e1d83109SMatthew G. Knepley       PetscInt        coneSize, supportSize;
266*e1d83109SMatthew G. Knepley 
26784961fc4SMatthew G. Knepley       ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
26884961fc4SMatthew G. Knepley       if (supportSize != 1) continue;
26984961fc4SMatthew G. Knepley       ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr);
27084961fc4SMatthew G. Knepley 
27184961fc4SMatthew G. Knepley       ierr = DMPlexGetCone(dm, support[0], &cone);CHKERRQ(ierr);
27284961fc4SMatthew G. Knepley       ierr = DMPlexGetConeSize(dm, support[0], &coneSize);CHKERRQ(ierr);
27384961fc4SMatthew G. Knepley       ierr = DMPlexGetConeOrientation(dm, support[0], &ornt);CHKERRQ(ierr);
27484961fc4SMatthew G. Knepley       for (c = 0; c < coneSize; ++c) if (cone[c] == face) break;
27584961fc4SMatthew G. Knepley       if (dim == 1) {
27684961fc4SMatthew G. Knepley         /* Use cone position instead, shifted to -1 or 1 */
277*e1d83109SMatthew G. Knepley         rorntComp[face].rank = c*2-1;
27884961fc4SMatthew G. Knepley       } else {
279*e1d83109SMatthew G. Knepley         if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 :  1;
280*e1d83109SMatthew G. Knepley         else                                                rorntComp[face].rank = ornt[c] < 0 ?  1 : -1;
28184961fc4SMatthew G. Knepley       }
282*e1d83109SMatthew G. Knepley       rorntComp[face].index = faceComp[face-fStart];
28384961fc4SMatthew G. Knepley     }
284*e1d83109SMatthew G. Knepley     /* Communicate boundary edge orientations */
285*e1d83109SMatthew G. Knepley     ierr = PetscSFBcastBegin(sf, MPIU_2INT, rorntComp, lorntComp);CHKERRQ(ierr);
286*e1d83109SMatthew G. Knepley     ierr = PetscSFBcastEnd(sf, MPIU_2INT, rorntComp, lorntComp);CHKERRQ(ierr);
287*e1d83109SMatthew G. Knepley   }
288*e1d83109SMatthew G. Knepley   /* Get process adjacency */
289*e1d83109SMatthew G. Knepley   ierr = PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors);CHKERRQ(ierr);
290*e1d83109SMatthew G. Knepley   for (comp = 0; comp < numComponents; ++comp) {
291*e1d83109SMatthew G. Knepley     PetscInt  l, n;
29284961fc4SMatthew G. Knepley 
293*e1d83109SMatthew G. Knepley     numNeighbors[comp] = 0;
294*e1d83109SMatthew G. Knepley     ierr = PetscMalloc1(numLeaves, &neighbors[comp]);CHKERRQ(ierr);
295*e1d83109SMatthew G. Knepley     /* I know this is p^2 time in general, but for bounded degree its alright */
296*e1d83109SMatthew G. Knepley     for (l = 0; l < numLeaves; ++l) {
297*e1d83109SMatthew G. Knepley       const PetscInt face = lpoints[l];
298*e1d83109SMatthew G. Knepley 
299*e1d83109SMatthew G. Knepley       /* Find a representative face (edge) separating pairs of procs */
300*e1d83109SMatthew G. Knepley       if ((face >= fStart) && (face < fEnd) && (faceComp[face-fStart] == comp)) {
301*e1d83109SMatthew G. Knepley         const PetscInt rank  = rpoints[l].rank;
302*e1d83109SMatthew G. Knepley         const PetscInt rcomp = lorntComp[face].index;
303*e1d83109SMatthew G. Knepley 
304*e1d83109SMatthew G. Knepley         for (n = 0; n < numNeighbors[comp]; ++n) if ((rank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break;
305*e1d83109SMatthew G. Knepley         if (n >= numNeighbors[comp]) {
306*e1d83109SMatthew G. Knepley           PetscInt supportSize;
307*e1d83109SMatthew G. Knepley 
308*e1d83109SMatthew G. Knepley           ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
309*e1d83109SMatthew G. Knepley           if (supportSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %d", supportSize);
310*e1d83109SMatthew G. Knepley           neighbors[comp][numNeighbors[comp]++] = l;
311*e1d83109SMatthew G. Knepley         }
312*e1d83109SMatthew G. Knepley       }
313*e1d83109SMatthew G. Knepley     }
314*e1d83109SMatthew G. Knepley     totNeighbors += numNeighbors[comp];
315*e1d83109SMatthew G. Knepley   }
316*e1d83109SMatthew G. Knepley   ierr = PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match);CHKERRQ(ierr);
317*e1d83109SMatthew G. Knepley   for (comp = 0, off = 0; comp < numComponents; ++comp) {
318*e1d83109SMatthew G. Knepley     PetscInt n;
319*e1d83109SMatthew G. Knepley 
320*e1d83109SMatthew G. Knepley     for (n = 0; n < numNeighbors[comp]; ++n, ++off) {
321*e1d83109SMatthew G. Knepley       const PetscInt face = lpoints[neighbors[comp][n]];
322*e1d83109SMatthew G. Knepley       const PetscInt o    = rorntComp[face].rank*lorntComp[face].rank;
323*e1d83109SMatthew G. Knepley 
324*e1d83109SMatthew G. Knepley       if      (o < 0) match[off] = PETSC_TRUE;
325*e1d83109SMatthew G. Knepley       else if (o > 0) match[off] = PETSC_FALSE;
326*e1d83109SMatthew G. Knepley       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);
327*e1d83109SMatthew G. Knepley       nrankComp[off].rank  = rpoints[neighbors[comp][n]].rank;
328*e1d83109SMatthew G. Knepley       nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index;
329*e1d83109SMatthew G. Knepley     }
330*e1d83109SMatthew G. Knepley     ierr = PetscFree(neighbors[comp]);CHKERRQ(ierr);
33184961fc4SMatthew G. Knepley   }
33284961fc4SMatthew G. Knepley   /* Collect the graph on 0 */
333*e1d83109SMatthew G. Knepley   if (numLeaves >= 0) {
33484961fc4SMatthew G. Knepley     Mat          G;
33584961fc4SMatthew G. Knepley     PetscBT      seenProcs, flippedProcs;
33684961fc4SMatthew G. Knepley     PetscInt    *procFIFO, pTop, pBottom;
337*e1d83109SMatthew G. Knepley     PetscInt    *N   = NULL, *Noff, *Nc;
338*e1d83109SMatthew G. Knepley     PetscSFNode *adj = NULL;
33984961fc4SMatthew G. Knepley     PetscBool   *val = NULL;
340*e1d83109SMatthew G. Knepley     PetscMPIInt *recvcounts = NULL, *displs = NULL, p, o;
341*e1d83109SMatthew G. Knepley     PetscMPIInt  numProcs = 0, rank;
34284961fc4SMatthew G. Knepley     PetscInt     debug = 0;
34384961fc4SMatthew G. Knepley 
344*e1d83109SMatthew G. Knepley     ierr = PetscCalloc1(numComponents, &flipped);CHKERRQ(ierr);
34584961fc4SMatthew G. Knepley     ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
34684961fc4SMatthew G. Knepley     if (!rank) {ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);}
347*e1d83109SMatthew G. Knepley     ierr = PetscCalloc4(numProcs, &recvcounts, numProcs+1, &displs, numProcs, &Nc, numProcs+1, &Noff);CHKERRQ(ierr);
348*e1d83109SMatthew G. Knepley     ierr = MPI_Gather(&numComponents, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm);CHKERRQ(ierr);
34984961fc4SMatthew G. Knepley     for (p = 0; p < numProcs; ++p) {
350*e1d83109SMatthew G. Knepley       displs[p+1] = displs[p] + Nc[p];
351*e1d83109SMatthew G. Knepley     }
352*e1d83109SMatthew G. Knepley     if (!rank) {ierr = PetscMalloc1(displs[numProcs],&N);CHKERRQ(ierr);}
353*e1d83109SMatthew G. Knepley     ierr = MPI_Gatherv(numNeighbors, numComponents, MPI_INT, N, Nc, displs, MPI_INT, 0, comm);CHKERRQ(ierr);
354*e1d83109SMatthew G. Knepley     for (p = 0, o = 0; p < numProcs; ++p) {
355*e1d83109SMatthew G. Knepley       recvcounts[p] = 0;
356*e1d83109SMatthew G. Knepley       for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o];
35784961fc4SMatthew G. Knepley       displs[p+1] = displs[p] + recvcounts[p];
35884961fc4SMatthew G. Knepley     }
35984961fc4SMatthew G. Knepley     if (!rank) {ierr = PetscMalloc2(displs[numProcs], &adj, displs[numProcs], &val);CHKERRQ(ierr);}
360*e1d83109SMatthew G. Knepley     ierr = MPI_Gatherv(nrankComp, totNeighbors, MPIU_2INT, adj, recvcounts, displs, MPIU_2INT, 0, comm);CHKERRQ(ierr);
361*e1d83109SMatthew G. Knepley     ierr = MPI_Gatherv(match, totNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
362*e1d83109SMatthew G. Knepley     ierr = PetscFree2(numNeighbors, neighbors);CHKERRQ(ierr);
363*e1d83109SMatthew G. Knepley     if (!rank) {
364*e1d83109SMatthew G. Knepley       for (p = 1; p <= numProcs; ++p) {Noff[p] = Noff[p-1] + Nc[p-1];}
36584961fc4SMatthew G. Knepley       if (debug) {
366*e1d83109SMatthew G. Knepley         PetscInt n;
367*e1d83109SMatthew G. Knepley 
368*e1d83109SMatthew G. Knepley         for (p = 0, off = 0; p < numProcs; ++p) {
369*e1d83109SMatthew G. Knepley           for (c = 0; c < Nc[p]; ++c) {
370*e1d83109SMatthew G. Knepley             ierr = PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %d:\n", p, c);
371*e1d83109SMatthew G. Knepley             for (n = 0; n < N[Noff[p]+c]; ++n, ++off) {
372*e1d83109SMatthew G. Knepley               ierr = PetscPrintf(PETSC_COMM_SELF, "  edge (%d, %d) (%d):\n", adj[off].rank, adj[off].index, val[off]);
373*e1d83109SMatthew G. Knepley             }
37484961fc4SMatthew G. Knepley           }
37584961fc4SMatthew G. Knepley         }
37684961fc4SMatthew G. Knepley       }
37784961fc4SMatthew G. Knepley       /* Symmetrize the graph */
37884961fc4SMatthew G. Knepley       ierr = MatCreate(PETSC_COMM_SELF, &G);CHKERRQ(ierr);
379*e1d83109SMatthew G. Knepley       ierr = MatSetSizes(G, Noff[numProcs], Noff[numProcs], Noff[numProcs], Noff[numProcs]);CHKERRQ(ierr);
38084961fc4SMatthew G. Knepley       ierr = MatSetUp(G);CHKERRQ(ierr);
381*e1d83109SMatthew G. Knepley       for (p = 0, off = 0; p < numProcs; ++p) {
382*e1d83109SMatthew G. Knepley         for (c = 0; c < Nc[p]; ++c) {
383*e1d83109SMatthew G. Knepley           const PetscInt r = Noff[p]+c;
384*e1d83109SMatthew G. Knepley           PetscInt       n;
385*e1d83109SMatthew G. Knepley 
386*e1d83109SMatthew G. Knepley           for (n = 0; n < N[r]; ++n, ++off) {
387*e1d83109SMatthew G. Knepley             const PetscInt    q = Noff[adj[off].rank] + adj[off].index;
388*e1d83109SMatthew G. Knepley             const PetscScalar o = val[off] ? 1.0 : 0.0;
38984961fc4SMatthew G. Knepley 
39084961fc4SMatthew G. Knepley             ierr = MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES);CHKERRQ(ierr);
39184961fc4SMatthew G. Knepley             ierr = MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES);CHKERRQ(ierr);
39284961fc4SMatthew G. Knepley           }
39384961fc4SMatthew G. Knepley         }
394*e1d83109SMatthew G. Knepley       }
39584961fc4SMatthew G. Knepley       ierr = MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
39684961fc4SMatthew G. Knepley       ierr = MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
39784961fc4SMatthew G. Knepley 
398*e1d83109SMatthew G. Knepley       ierr = PetscBTCreate(Noff[numProcs], &seenProcs);CHKERRQ(ierr);
399*e1d83109SMatthew G. Knepley       ierr = PetscBTMemzero(Noff[numProcs], seenProcs);CHKERRQ(ierr);
400*e1d83109SMatthew G. Knepley       ierr = PetscBTCreate(Noff[numProcs], &flippedProcs);CHKERRQ(ierr);
401*e1d83109SMatthew G. Knepley       ierr = PetscBTMemzero(Noff[numProcs], flippedProcs);CHKERRQ(ierr);
402*e1d83109SMatthew G. Knepley       ierr = PetscMalloc1(Noff[numProcs], &procFIFO);CHKERRQ(ierr);
40384961fc4SMatthew G. Knepley       pTop = pBottom = 0;
404*e1d83109SMatthew G. Knepley       for (p = 0; p < Noff[numProcs]; ++p) {
40584961fc4SMatthew G. Knepley         if (PetscBTLookup(seenProcs, p)) continue;
40684961fc4SMatthew G. Knepley         /* Initialize FIFO with next proc */
40784961fc4SMatthew G. Knepley         procFIFO[pBottom++] = p;
40884961fc4SMatthew G. Knepley         ierr = PetscBTSet(seenProcs, p);CHKERRQ(ierr);
40984961fc4SMatthew G. Knepley         /* Consider each proc in FIFO */
41084961fc4SMatthew G. Knepley         while (pTop < pBottom) {
41184961fc4SMatthew G. Knepley           const PetscScalar *ornt;
41284961fc4SMatthew G. Knepley           const PetscInt    *neighbors;
413*e1d83109SMatthew G. Knepley           PetscInt           proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n;
41484961fc4SMatthew G. Knepley 
41584961fc4SMatthew G. Knepley           proc     = procFIFO[pTop++];
41684961fc4SMatthew G. Knepley           flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0;
41784961fc4SMatthew G. Knepley           ierr = MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt);CHKERRQ(ierr);
41884961fc4SMatthew G. Knepley           /* Loop over neighboring procs */
41984961fc4SMatthew G. Knepley           for (n = 0; n < numNeighbors; ++n) {
42084961fc4SMatthew G. Knepley             nproc    = neighbors[n];
42184961fc4SMatthew G. Knepley             mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1;
42284961fc4SMatthew G. Knepley             seen     = PetscBTLookup(seenProcs, nproc);
42384961fc4SMatthew G. Knepley             flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0;
42484961fc4SMatthew G. Knepley 
42584961fc4SMatthew G. Knepley             if (mismatch ^ (flippedA ^ flippedB)) {
42684961fc4SMatthew G. Knepley               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);
42784961fc4SMatthew G. Knepley               if (!flippedB) {
42884961fc4SMatthew G. Knepley                 ierr = PetscBTSet(flippedProcs, nproc);CHKERRQ(ierr);
42984961fc4SMatthew G. Knepley               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
43084961fc4SMatthew G. Knepley             } else if (mismatch && flippedA && flippedB) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
43184961fc4SMatthew G. Knepley             if (!seen) {
43284961fc4SMatthew G. Knepley               procFIFO[pBottom++] = nproc;
43384961fc4SMatthew G. Knepley               ierr = PetscBTSet(seenProcs, nproc);CHKERRQ(ierr);
43484961fc4SMatthew G. Knepley             }
43584961fc4SMatthew G. Knepley           }
43684961fc4SMatthew G. Knepley         }
43784961fc4SMatthew G. Knepley       }
43884961fc4SMatthew G. Knepley       ierr = PetscFree(procFIFO);CHKERRQ(ierr);
43984961fc4SMatthew G. Knepley       ierr = MatDestroy(&G);CHKERRQ(ierr);
44084961fc4SMatthew G. Knepley       ierr = PetscFree2(adj, val);CHKERRQ(ierr);
441*e1d83109SMatthew G. Knepley       ierr = PetscBTDestroy(&seenProcs);CHKERRQ(ierr);
44284961fc4SMatthew G. Knepley     }
443*e1d83109SMatthew G. Knepley     /* Scatter flip flags */
444*e1d83109SMatthew G. Knepley     {
445*e1d83109SMatthew G. Knepley       PetscBool *flips = NULL;
446*e1d83109SMatthew G. Knepley 
447*e1d83109SMatthew G. Knepley       if (!rank) {
448*e1d83109SMatthew G. Knepley         ierr = PetscMalloc1(Noff[numProcs], &flips);CHKERRQ(ierr);
449*e1d83109SMatthew G. Knepley         for (p = 0; p < Noff[numProcs]; ++p) {
450*e1d83109SMatthew G. Knepley           flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE;
451*e1d83109SMatthew G. Knepley           if (debug && flips[p]) {ierr = PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p);}
452*e1d83109SMatthew G. Knepley         }
453*e1d83109SMatthew G. Knepley         for (p = 0; p < numProcs; ++p) {
454*e1d83109SMatthew G. Knepley           displs[p+1] = displs[p] + Nc[p];
455*e1d83109SMatthew G. Knepley         }
456*e1d83109SMatthew G. Knepley       }
457*e1d83109SMatthew G. Knepley       ierr = MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, numComponents, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
45884961fc4SMatthew G. Knepley       ierr = PetscFree(flips);CHKERRQ(ierr);
45984961fc4SMatthew G. Knepley     }
460*e1d83109SMatthew G. Knepley     if (!rank) {ierr = PetscBTDestroy(&flippedProcs);CHKERRQ(ierr);}
461*e1d83109SMatthew G. Knepley     ierr = PetscFree(N);CHKERRQ(ierr);
462*e1d83109SMatthew G. Knepley     ierr = PetscFree4(recvcounts, displs, Nc, Noff);CHKERRQ(ierr);
463*e1d83109SMatthew G. Knepley     ierr = PetscFree2(nrankComp, match);CHKERRQ(ierr);
464*e1d83109SMatthew G. Knepley 
465*e1d83109SMatthew G. Knepley     /* Decide whether to flip cells in each component */
466*e1d83109SMatthew G. Knepley     for (c = 0; c < cEnd-cStart; ++c) {if (flipped[cellComp[c]]) {ierr = PetscBTNegate(flippedCells, c);CHKERRQ(ierr);}}
467*e1d83109SMatthew G. Knepley     ierr = PetscFree(flipped);CHKERRQ(ierr);
46884961fc4SMatthew G. Knepley   }
46984961fc4SMatthew G. Knepley   if (flg) {
47084961fc4SMatthew G. Knepley     PetscViewer v;
47184961fc4SMatthew G. Knepley     PetscMPIInt rank;
47284961fc4SMatthew G. Knepley 
47384961fc4SMatthew G. Knepley     ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
4747e09831bSMatthew G. Knepley     ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr);
47584961fc4SMatthew G. Knepley     ierr = PetscViewerASCIISynchronizedAllow(v, PETSC_TRUE);CHKERRQ(ierr);
47684961fc4SMatthew G. Knepley     ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank);CHKERRQ(ierr);
47784961fc4SMatthew G. Knepley     ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr);
47884961fc4SMatthew G. Knepley   }
47984961fc4SMatthew G. Knepley   /* Reverse flipped cells in the mesh */
48084961fc4SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
48184961fc4SMatthew G. Knepley     if (PetscBTLookup(flippedCells, c-cStart)) {ierr = DMPlexReverseCell(dm, c);CHKERRQ(ierr);}
48284961fc4SMatthew G. Knepley   }
48384961fc4SMatthew G. Knepley   ierr = PetscBTDestroy(&seenCells);CHKERRQ(ierr);
48484961fc4SMatthew G. Knepley   ierr = PetscBTDestroy(&flippedCells);CHKERRQ(ierr);
48584961fc4SMatthew G. Knepley   ierr = PetscBTDestroy(&seenFaces);CHKERRQ(ierr);
486*e1d83109SMatthew G. Knepley   ierr = PetscFree2(numNeighbors, neighbors);CHKERRQ(ierr);
487*e1d83109SMatthew G. Knepley   ierr = PetscFree2(rorntComp, lorntComp);CHKERRQ(ierr);
488*e1d83109SMatthew G. Knepley   ierr = PetscFree3(faceFIFO, cellComp, faceComp);CHKERRQ(ierr);
48984961fc4SMatthew G. Knepley   PetscFunctionReturn(0);
49084961fc4SMatthew G. Knepley }
491