xref: /petsc/src/dm/impls/plex/plexorient.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
284961fc4SMatthew G. Knepley #include <petscsf.h>
384961fc4SMatthew G. Knepley 
484961fc4SMatthew G. Knepley /*@
5b5a892a1SMatthew G. Knepley   DMPlexOrientPoint - Act with the given orientation on the cone points of this mesh point, and update its use in the mesh.
627d0e99aSVaclav Hapla 
727d0e99aSVaclav Hapla   Not Collective
827d0e99aSVaclav Hapla 
927d0e99aSVaclav Hapla   Input Parameters:
1027d0e99aSVaclav Hapla + dm - The DM
11b5a892a1SMatthew G. Knepley . p  - The mesh point
12b5a892a1SMatthew G. Knepley - o  - The orientation
1327d0e99aSVaclav Hapla 
1427d0e99aSVaclav Hapla   Level: intermediate
1527d0e99aSVaclav Hapla 
1627d0e99aSVaclav Hapla .seealso: DMPlexOrient(), DMPlexGetCone(), DMPlexGetConeOrientation(), DMPlexInterpolate(), DMPlexGetChart()
1727d0e99aSVaclav Hapla @*/
18b5a892a1SMatthew G. Knepley PetscErrorCode DMPlexOrientPoint(DM dm, PetscInt p, PetscInt o)
1927d0e99aSVaclav Hapla {
20b5a892a1SMatthew G. Knepley   DMPolytopeType  ct;
21b5a892a1SMatthew G. Knepley   const PetscInt *arr, *cone, *ornt, *support;
22b5a892a1SMatthew G. Knepley   PetscInt       *newcone, *newornt;
23b5a892a1SMatthew G. Knepley   PetscInt        coneSize, c, supportSize, s;
2427d0e99aSVaclav Hapla 
2527d0e99aSVaclav Hapla   PetscFunctionBegin;
2627d0e99aSVaclav Hapla   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetCellType(dm, p, &ct));
28b5a892a1SMatthew G. Knepley   arr  = DMPolytopeTypeGetArrangment(ct, o);
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetConeSize(dm, p, &coneSize));
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetCone(dm, p, &cone));
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetConeOrientation(dm, p, &ornt));
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetWorkArray(dm, coneSize, MPIU_INT, &newcone));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetWorkArray(dm, coneSize, MPIU_INT, &newornt));
34b5a892a1SMatthew G. Knepley   for (c = 0; c < coneSize; ++c) {
35b5a892a1SMatthew G. Knepley     DMPolytopeType ft;
36b5a892a1SMatthew G. Knepley     PetscInt       nO;
3727d0e99aSVaclav Hapla 
38*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetCellType(dm, cone[c], &ft));
39b5a892a1SMatthew G. Knepley     nO   = DMPolytopeTypeGetNumArrangments(ft)/2;
40b5a892a1SMatthew G. Knepley     newcone[c] = cone[arr[c*2+0]];
41b5a892a1SMatthew G. Knepley     newornt[c] = DMPolytopeTypeComposeOrientation(ft, arr[c*2+1], ornt[arr[c*2+0]]);
422c71b3e2SJacob Faibussowitsch     PetscCheckFalse(newornt[c] && (newornt[c] >= nO || newornt[c] < -nO),PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid orientation %D not in [%D,%D) for %s %D", newornt[c], -nO, nO, DMPolytopeTypes[ft], cone[c]);
4327d0e99aSVaclav Hapla   }
44*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexSetCone(dm, p, newcone));
45*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexSetConeOrientation(dm, p, newornt));
46*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreWorkArray(dm, coneSize, MPIU_INT, &newcone));
47*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreWorkArray(dm, coneSize, MPIU_INT, &newornt));
48b5a892a1SMatthew G. Knepley   /* Update orientation of this point in the support points */
49*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetSupportSize(dm, p, &supportSize));
50*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetSupport(dm, p, &support));
51b5a892a1SMatthew G. Knepley   for (s = 0; s < supportSize; ++s) {
52*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetConeSize(dm, support[s], &coneSize));
53*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetCone(dm, support[s], &cone));
54*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetConeOrientation(dm, support[s], &ornt));
55b5a892a1SMatthew G. Knepley     for (c = 0; c < coneSize; ++c) {
56b5a892a1SMatthew G. Knepley       PetscInt po;
5727d0e99aSVaclav Hapla 
58b5a892a1SMatthew G. Knepley       if (cone[c] != p) continue;
59b5a892a1SMatthew G. Knepley       /* ornt[c] * 0 = target = po * o so that po = ornt[c] * o^{-1} */
60b5a892a1SMatthew G. Knepley       po   = DMPolytopeTypeComposeOrientationInv(ct, ornt[c], o);
61*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexInsertConeOrientation(dm, support[s], c, po));
6284961fc4SMatthew G. Knepley     }
6384961fc4SMatthew G. Knepley   }
6484961fc4SMatthew G. Knepley   PetscFunctionReturn(0);
6584961fc4SMatthew G. Knepley }
6684961fc4SMatthew G. Knepley 
677b310ce2SMatthew G. Knepley /*
687b310ce2SMatthew G. Knepley   - Checks face match
697b310ce2SMatthew G. Knepley     - Flips non-matching
707b310ce2SMatthew G. Knepley   - Inserts faces of support cells in FIFO
717b310ce2SMatthew G. Knepley */
727b2de0fdSMatthew G. Knepley static PetscErrorCode DMPlexCheckFace_Internal(DM dm, PetscInt *faceFIFO, PetscInt *fTop, PetscInt *fBottom, PetscInt cStart, PetscInt fStart, PetscInt fEnd, PetscBT seenCells, PetscBT flippedCells, PetscBT seenFaces)
737b310ce2SMatthew G. Knepley {
747b310ce2SMatthew G. Knepley   const PetscInt *support, *coneA, *coneB, *coneOA, *coneOB;
757b310ce2SMatthew G. Knepley   PetscInt        supportSize, coneSizeA, coneSizeB, posA = -1, posB = -1;
767b2de0fdSMatthew G. Knepley   PetscInt        face, dim, seenA, flippedA, seenB, flippedB, mismatch, c;
777b310ce2SMatthew G. Knepley 
787b310ce2SMatthew G. Knepley   PetscFunctionBegin;
797b310ce2SMatthew G. Knepley   face = faceFIFO[(*fTop)++];
80*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
81*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetSupportSize(dm, face, &supportSize));
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetSupport(dm, face, &support));
837b310ce2SMatthew G. Knepley   if (supportSize < 2) PetscFunctionReturn(0);
842c71b3e2SJacob Faibussowitsch   PetscCheckFalse(supportSize != 2,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %d", supportSize);
857b310ce2SMatthew G. Knepley   seenA    = PetscBTLookup(seenCells,    support[0]-cStart);
867b310ce2SMatthew G. Knepley   flippedA = PetscBTLookup(flippedCells, support[0]-cStart) ? 1 : 0;
877b310ce2SMatthew G. Knepley   seenB    = PetscBTLookup(seenCells,    support[1]-cStart);
887b310ce2SMatthew G. Knepley   flippedB = PetscBTLookup(flippedCells, support[1]-cStart) ? 1 : 0;
897b310ce2SMatthew G. Knepley 
90*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetConeSize(dm, support[0], &coneSizeA));
91*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetConeSize(dm, support[1], &coneSizeB));
92*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetCone(dm, support[0], &coneA));
93*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetCone(dm, support[1], &coneB));
94*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetConeOrientation(dm, support[0], &coneOA));
95*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetConeOrientation(dm, support[1], &coneOB));
967b310ce2SMatthew G. Knepley   for (c = 0; c < coneSizeA; ++c) {
977b310ce2SMatthew G. Knepley     if (!PetscBTLookup(seenFaces, coneA[c]-fStart)) {
987b310ce2SMatthew G. Knepley       faceFIFO[(*fBottom)++] = coneA[c];
99*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscBTSet(seenFaces, coneA[c]-fStart));
1007b310ce2SMatthew G. Knepley     }
1017b310ce2SMatthew G. Knepley     if (coneA[c] == face) posA = c;
1022c71b3e2SJacob Faibussowitsch     PetscCheckFalse(*fBottom > fEnd-fStart,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %d was pushed exceeding capacity %d > %d", coneA[c], *fBottom, fEnd-fStart);
1037b310ce2SMatthew G. Knepley   }
1042c71b3e2SJacob Faibussowitsch   PetscCheckFalse(posA < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d could not be located in cell %d", face, support[0]);
1057b310ce2SMatthew G. Knepley   for (c = 0; c < coneSizeB; ++c) {
1067b310ce2SMatthew G. Knepley     if (!PetscBTLookup(seenFaces, coneB[c]-fStart)) {
1077b310ce2SMatthew G. Knepley       faceFIFO[(*fBottom)++] = coneB[c];
108*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscBTSet(seenFaces, coneB[c]-fStart));
1097b310ce2SMatthew G. Knepley     }
1107b310ce2SMatthew G. Knepley     if (coneB[c] == face) posB = c;
1112c71b3e2SJacob Faibussowitsch     PetscCheckFalse(*fBottom > fEnd-fStart,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %d was pushed exceeding capacity %d > %d", coneA[c], *fBottom, fEnd-fStart);
1127b310ce2SMatthew G. Knepley   }
1132c71b3e2SJacob Faibussowitsch   PetscCheckFalse(posB < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d could not be located in cell %d", face, support[1]);
1147b310ce2SMatthew G. Knepley 
1157b310ce2SMatthew G. Knepley   if (dim == 1) {
1167b310ce2SMatthew G. Knepley     mismatch = posA == posB;
1177b310ce2SMatthew G. Knepley   } else {
1187b310ce2SMatthew G. Knepley     mismatch = coneOA[posA] == coneOB[posB];
1197b310ce2SMatthew G. Knepley   }
1207b310ce2SMatthew G. Knepley 
1217b310ce2SMatthew G. Knepley   if (mismatch ^ (flippedA ^ flippedB)) {
1222c71b3e2SJacob Faibussowitsch     PetscCheckFalse(seenA && seenB,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]);
1237b310ce2SMatthew G. Knepley     if (!seenA && !flippedA) {
124*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscBTSet(flippedCells, support[0]-cStart));
1257b310ce2SMatthew G. Knepley     } else if (!seenB && !flippedB) {
126*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscBTSet(flippedCells, support[1]-cStart));
1277b310ce2SMatthew G. Knepley     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
1282c71b3e2SJacob Faibussowitsch   } else PetscCheckFalse(mismatch && flippedA && flippedB,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
129*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBTSet(seenCells, support[0]-cStart));
130*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBTSet(seenCells, support[1]-cStart));
1317b310ce2SMatthew G. Knepley   PetscFunctionReturn(0);
1327b310ce2SMatthew G. Knepley }
1337b310ce2SMatthew G. Knepley 
13484961fc4SMatthew G. Knepley /*@
13584961fc4SMatthew G. Knepley   DMPlexOrient - Give a consistent orientation to the input mesh
13684961fc4SMatthew G. Knepley 
13784961fc4SMatthew G. Knepley   Input Parameters:
13884961fc4SMatthew G. Knepley . dm - The DM
13984961fc4SMatthew G. Knepley 
14084961fc4SMatthew G. Knepley   Note: The orientation data for the DM are change in-place.
14184961fc4SMatthew G. Knepley $ This routine will fail for non-orientable surfaces, such as the Moebius strip.
14284961fc4SMatthew G. Knepley 
14384961fc4SMatthew G. Knepley   Level: advanced
14484961fc4SMatthew G. Knepley 
14584961fc4SMatthew G. Knepley .seealso: DMCreate(), DMPLEX
14684961fc4SMatthew G. Knepley @*/
14784961fc4SMatthew G. Knepley PetscErrorCode DMPlexOrient(DM dm)
14884961fc4SMatthew G. Knepley {
14984961fc4SMatthew G. Knepley   MPI_Comm           comm;
150e1d83109SMatthew G. Knepley   PetscSF            sf;
151e1d83109SMatthew G. Knepley   const PetscInt    *lpoints;
152e1d83109SMatthew G. Knepley   const PetscSFNode *rpoints;
153e1d83109SMatthew G. Knepley   PetscSFNode       *rorntComp = NULL, *lorntComp = NULL;
154e1d83109SMatthew G. Knepley   PetscInt          *numNeighbors, **neighbors;
155e1d83109SMatthew G. Knepley   PetscSFNode       *nrankComp;
156e1d83109SMatthew G. Knepley   PetscBool         *match, *flipped;
15784961fc4SMatthew G. Knepley   PetscBT            seenCells, flippedCells, seenFaces;
158e1d83109SMatthew G. Knepley   PetscInt          *faceFIFO, fTop, fBottom, *cellComp, *faceComp;
1597cadcfe8SMatthew G. Knepley   PetscInt           numLeaves, numRoots, dim, h, cStart, cEnd, c, cell, fStart, fEnd, face, off, totNeighbors = 0;
16039526728SToby Isaac   PetscMPIInt        rank, size, numComponents, comp = 0;
16139526728SToby Isaac   PetscBool          flg, flg2;
162a17aa47eSToby Isaac   PetscViewer        viewer = NULL, selfviewer = NULL;
16384961fc4SMatthew G. Knepley 
16484961fc4SMatthew G. Knepley   PetscFunctionBegin;
165*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm));
166*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm, &rank));
167*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm, &size));
168*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view", &flg));
169*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view_synchronized", &flg2));
170*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetPointSF(dm, &sf));
171*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints));
17284961fc4SMatthew G. Knepley   /* Truth Table
17384961fc4SMatthew G. Knepley      mismatch    flips   do action   mismatch   flipA ^ flipB   action
17484961fc4SMatthew G. Knepley          F       0 flips     no         F             F           F
17584961fc4SMatthew G. Knepley          F       1 flip      yes        F             T           T
17684961fc4SMatthew G. Knepley          F       2 flips     no         T             F           T
17784961fc4SMatthew G. Knepley          T       0 flips     yes        T             T           F
17884961fc4SMatthew G. Knepley          T       1 flip      no
17984961fc4SMatthew G. Knepley          T       2 flips     yes
18084961fc4SMatthew G. Knepley   */
181*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
182*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetVTKCellHeight(dm, &h));
183*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, h,   &cStart, &cEnd));
184*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, h+1, &fStart, &fEnd));
185*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBTCreate(cEnd - cStart, &seenCells));
186*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBTMemzero(cEnd - cStart, seenCells));
187*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBTCreate(cEnd - cStart, &flippedCells));
188*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBTMemzero(cEnd - cStart, flippedCells));
189*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBTCreate(fEnd - fStart, &seenFaces));
190*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBTMemzero(fEnd - fStart, seenFaces));
191*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc3(fEnd - fStart, &faceFIFO, cEnd-cStart, &cellComp, fEnd-fStart, &faceComp));
192e1d83109SMatthew G. Knepley   /*
193e1d83109SMatthew G. Knepley    OLD STYLE
194e1d83109SMatthew G. Knepley    - Add an integer array over cells and faces (component) for connected component number
195e1d83109SMatthew G. Knepley    Foreach component
196e1d83109SMatthew G. Knepley      - Mark the initial cell as seen
197e1d83109SMatthew G. Knepley      - Process component as usual
198e1d83109SMatthew G. Knepley      - Set component for all seenCells
199e1d83109SMatthew G. Knepley      - Wipe seenCells and seenFaces (flippedCells can stay)
200e1d83109SMatthew G. Knepley    - Generate parallel adjacency for component using SF and seenFaces
201e1d83109SMatthew G. Knepley    - Collect numComponents adj data from each proc to 0
202e1d83109SMatthew G. Knepley    - Build same serial graph
203e1d83109SMatthew G. Knepley    - Use same solver
204e1d83109SMatthew G. Knepley    - Use Scatterv to to send back flipped flags for each component
205e1d83109SMatthew G. Knepley    - Negate flippedCells by component
206e1d83109SMatthew G. Knepley 
207e1d83109SMatthew G. Knepley    NEW STYLE
208e1d83109SMatthew G. Knepley    - Create the adj on each process
209e1d83109SMatthew G. Knepley    - Bootstrap to complete graph on proc 0
210e1d83109SMatthew G. Knepley   */
211e1d83109SMatthew G. Knepley   /* Loop over components */
212e1d83109SMatthew G. Knepley   for (cell = cStart; cell < cEnd; ++cell) cellComp[cell-cStart] = -1;
213e1d83109SMatthew G. Knepley   do {
214e1d83109SMatthew G. Knepley     /* Look for first unmarked cell */
215e1d83109SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) if (cellComp[cell-cStart] < 0) break;
216e1d83109SMatthew G. Knepley     if (cell >= cEnd) break;
217e1d83109SMatthew G. Knepley     /* Initialize FIFO with first cell in component */
218e1d83109SMatthew G. Knepley     {
21984961fc4SMatthew G. Knepley       const PetscInt *cone;
22084961fc4SMatthew G. Knepley       PetscInt        coneSize;
22184961fc4SMatthew G. Knepley 
222e1d83109SMatthew G. Knepley       fTop = fBottom = 0;
223*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetConeSize(dm, cell, &coneSize));
224*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetCone(dm, cell, &cone));
22584961fc4SMatthew G. Knepley       for (c = 0; c < coneSize; ++c) {
22684961fc4SMatthew G. Knepley         faceFIFO[fBottom++] = cone[c];
227*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscBTSet(seenFaces, cone[c]-fStart));
22884961fc4SMatthew G. Knepley       }
229*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscBTSet(seenCells, cell-cStart));
23084961fc4SMatthew G. Knepley     }
23184961fc4SMatthew G. Knepley     /* Consider each face in FIFO */
23284961fc4SMatthew G. Knepley     while (fTop < fBottom) {
233*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexCheckFace_Internal(dm, faceFIFO, &fTop, &fBottom, cStart, fStart, fEnd, seenCells, flippedCells, seenFaces));
23484961fc4SMatthew G. Knepley     }
235e1d83109SMatthew G. Knepley     /* Set component for cells and faces */
236e1d83109SMatthew G. Knepley     for (cell = 0; cell < cEnd-cStart; ++cell) {
237e1d83109SMatthew G. Knepley       if (PetscBTLookup(seenCells, cell)) cellComp[cell] = comp;
238e1d83109SMatthew G. Knepley     }
239e1d83109SMatthew G. Knepley     for (face = 0; face < fEnd-fStart; ++face) {
240e1d83109SMatthew G. Knepley       if (PetscBTLookup(seenFaces, face)) faceComp[face] = comp;
241e1d83109SMatthew G. Knepley     }
242e1d83109SMatthew G. Knepley     /* Wipe seenCells and seenFaces for next component */
243*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBTMemzero(fEnd - fStart, seenFaces));
244*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBTMemzero(cEnd - cStart, seenCells));
245e1d83109SMatthew G. Knepley     ++comp;
246e1d83109SMatthew G. Knepley   } while (1);
247e1d83109SMatthew G. Knepley   numComponents = comp;
24884961fc4SMatthew G. Knepley   if (flg) {
24984961fc4SMatthew G. Knepley     PetscViewer v;
25084961fc4SMatthew G. Knepley 
251*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIGetStdout(comm, &v));
252*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPushSynchronized(v));
253*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank));
254*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBTView(cEnd-cStart, flippedCells, v));
255*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerFlush(v));
256*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPopSynchronized(v));
25784961fc4SMatthew G. Knepley   }
25884961fc4SMatthew G. Knepley   /* Now all subdomains are oriented, but we need a consistent parallel orientation */
25984961fc4SMatthew G. Knepley   if (numLeaves >= 0) {
260e1d83109SMatthew G. Knepley     /* Store orientations of boundary faces*/
261*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc2(numRoots,&rorntComp,numRoots,&lorntComp));
26284961fc4SMatthew G. Knepley     for (face = fStart; face < fEnd; ++face) {
263e1d83109SMatthew G. Knepley       const PetscInt *cone, *support, *ornt;
264e1d83109SMatthew G. Knepley       PetscInt        coneSize, supportSize;
265e1d83109SMatthew G. Knepley 
266*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetSupportSize(dm, face, &supportSize));
26784961fc4SMatthew G. Knepley       if (supportSize != 1) continue;
268*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetSupport(dm, face, &support));
26984961fc4SMatthew G. Knepley 
270*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetCone(dm, support[0], &cone));
271*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetConeSize(dm, support[0], &coneSize));
272*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetConeOrientation(dm, support[0], &ornt));
27384961fc4SMatthew G. Knepley       for (c = 0; c < coneSize; ++c) if (cone[c] == face) break;
27484961fc4SMatthew G. Knepley       if (dim == 1) {
27584961fc4SMatthew G. Knepley         /* Use cone position instead, shifted to -1 or 1 */
2766e1eb25cSMatthew G. Knepley         if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = 1-c*2;
2776e1eb25cSMatthew G. Knepley         else                                                rorntComp[face].rank = c*2-1;
27884961fc4SMatthew G. Knepley       } else {
279e1d83109SMatthew G. Knepley         if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 :  1;
280e1d83109SMatthew G. Knepley         else                                                rorntComp[face].rank = ornt[c] < 0 ?  1 : -1;
28184961fc4SMatthew G. Knepley       }
282e1d83109SMatthew G. Knepley       rorntComp[face].index = faceComp[face-fStart];
28384961fc4SMatthew G. Knepley     }
284e1d83109SMatthew G. Knepley     /* Communicate boundary edge orientations */
285*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFBcastBegin(sf, MPIU_2INT, rorntComp, lorntComp,MPI_REPLACE));
286*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFBcastEnd(sf, MPIU_2INT, rorntComp, lorntComp,MPI_REPLACE));
287e1d83109SMatthew G. Knepley   }
288e1d83109SMatthew G. Knepley   /* Get process adjacency */
289*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors));
290a17aa47eSToby Isaac   viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm));
291*5f80ce2aSJacob Faibussowitsch   if (flg2) CHKERRQ(PetscViewerASCIIPushSynchronized(viewer));
292*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&selfviewer));
293e1d83109SMatthew G. Knepley   for (comp = 0; comp < numComponents; ++comp) {
294e1d83109SMatthew G. Knepley     PetscInt  l, n;
29584961fc4SMatthew G. Knepley 
296e1d83109SMatthew G. Knepley     numNeighbors[comp] = 0;
297*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp]));
298e1d83109SMatthew G. Knepley     /* I know this is p^2 time in general, but for bounded degree its alright */
299e1d83109SMatthew G. Knepley     for (l = 0; l < numLeaves; ++l) {
300e1d83109SMatthew G. Knepley       const PetscInt face = lpoints[l];
301e1d83109SMatthew G. Knepley 
302e1d83109SMatthew G. Knepley       /* Find a representative face (edge) separating pairs of procs */
303e1d83109SMatthew G. Knepley       if ((face >= fStart) && (face < fEnd) && (faceComp[face-fStart] == comp)) {
304269a49d5SMatthew G. Knepley         const PetscInt rrank = rpoints[l].rank;
305e1d83109SMatthew G. Knepley         const PetscInt rcomp = lorntComp[face].index;
306e1d83109SMatthew G. Knepley 
307269a49d5SMatthew G. Knepley         for (n = 0; n < numNeighbors[comp]; ++n) if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break;
308e1d83109SMatthew G. Knepley         if (n >= numNeighbors[comp]) {
309e1d83109SMatthew G. Knepley           PetscInt supportSize;
310e1d83109SMatthew G. Knepley 
311*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMPlexGetSupportSize(dm, face, &supportSize));
3122c71b3e2SJacob Faibussowitsch           PetscCheckFalse(supportSize != 1,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %d", supportSize);
313*5f80ce2aSJacob Faibussowitsch           if (flg) CHKERRQ(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));
314e1d83109SMatthew G. Knepley           neighbors[comp][numNeighbors[comp]++] = l;
315e1d83109SMatthew G. Knepley         }
316e1d83109SMatthew G. Knepley       }
317e1d83109SMatthew G. Knepley     }
318e1d83109SMatthew G. Knepley     totNeighbors += numNeighbors[comp];
319e1d83109SMatthew G. Knepley   }
320*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&selfviewer));
321*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerFlush(viewer));
322*5f80ce2aSJacob Faibussowitsch   if (flg2) CHKERRQ(PetscViewerASCIIPopSynchronized(viewer));
323*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match));
324e1d83109SMatthew G. Knepley   for (comp = 0, off = 0; comp < numComponents; ++comp) {
325e1d83109SMatthew G. Knepley     PetscInt n;
326e1d83109SMatthew G. Knepley 
327e1d83109SMatthew G. Knepley     for (n = 0; n < numNeighbors[comp]; ++n, ++off) {
328e1d83109SMatthew G. Knepley       const PetscInt face = lpoints[neighbors[comp][n]];
329e1d83109SMatthew G. Knepley       const PetscInt o    = rorntComp[face].rank*lorntComp[face].rank;
330e1d83109SMatthew G. Knepley 
331e1d83109SMatthew G. Knepley       if      (o < 0) match[off] = PETSC_TRUE;
332e1d83109SMatthew G. Knepley       else if (o > 0) match[off] = PETSC_FALSE;
33398921bdaSJacob Faibussowitsch       else SETERRQ(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);
334e1d83109SMatthew G. Knepley       nrankComp[off].rank  = rpoints[neighbors[comp][n]].rank;
335e1d83109SMatthew G. Knepley       nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index;
336e1d83109SMatthew G. Knepley     }
337*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(neighbors[comp]));
33884961fc4SMatthew G. Knepley   }
33984961fc4SMatthew G. Knepley   /* Collect the graph on 0 */
340e1d83109SMatthew G. Knepley   if (numLeaves >= 0) {
34184961fc4SMatthew G. Knepley     Mat          G;
34284961fc4SMatthew G. Knepley     PetscBT      seenProcs, flippedProcs;
34384961fc4SMatthew G. Knepley     PetscInt    *procFIFO, pTop, pBottom;
3447cadcfe8SMatthew G. Knepley     PetscInt    *N   = NULL, *Noff;
345e1d83109SMatthew G. Knepley     PetscSFNode *adj = NULL;
34684961fc4SMatthew G. Knepley     PetscBool   *val = NULL;
3477cadcfe8SMatthew G. Knepley     PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc, p, o;
3489852e123SBarry Smith     PetscMPIInt  size = 0;
34984961fc4SMatthew G. Knepley 
350*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc1(numComponents, &flipped));
351*5f80ce2aSJacob Faibussowitsch     if (rank == 0) CHKERRMPI(MPI_Comm_size(comm, &size));
352*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc4(size, &recvcounts, size+1, &displs, size, &Nc, size+1, &Noff));
353*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Gather(&numComponents, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm));
3549852e123SBarry Smith     for (p = 0; p < size; ++p) {
355e1d83109SMatthew G. Knepley       displs[p+1] = displs[p] + Nc[p];
356e1d83109SMatthew G. Knepley     }
357*5f80ce2aSJacob Faibussowitsch     if (rank == 0) CHKERRQ(PetscMalloc1(displs[size],&N));
358*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Gatherv(numNeighbors, numComponents, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm));
3599852e123SBarry Smith     for (p = 0, o = 0; p < size; ++p) {
360e1d83109SMatthew G. Knepley       recvcounts[p] = 0;
361e1d83109SMatthew G. Knepley       for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o];
36284961fc4SMatthew G. Knepley       displs[p+1] = displs[p] + recvcounts[p];
36384961fc4SMatthew G. Knepley     }
364*5f80ce2aSJacob Faibussowitsch     if (rank == 0) CHKERRQ(PetscMalloc2(displs[size], &adj, displs[size], &val));
365*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Gatherv(nrankComp, totNeighbors, MPIU_2INT, adj, recvcounts, displs, MPIU_2INT, 0, comm));
366*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Gatherv(match, totNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm));
367*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree2(numNeighbors, neighbors));
368dd400576SPatrick Sanan     if (rank == 0) {
3699852e123SBarry Smith       for (p = 1; p <= size; ++p) {Noff[p] = Noff[p-1] + Nc[p-1];}
370a83982f3SMatthew G. Knepley       if (flg) {
371e1d83109SMatthew G. Knepley         PetscInt n;
372e1d83109SMatthew G. Knepley 
3739852e123SBarry Smith         for (p = 0, off = 0; p < size; ++p) {
374e1d83109SMatthew G. Knepley           for (c = 0; c < Nc[p]; ++c) {
375*5f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %d:\n", p, c));
376e1d83109SMatthew G. Knepley             for (n = 0; n < N[Noff[p]+c]; ++n, ++off) {
377*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "  edge (%d, %d) (%d):\n", adj[off].rank, adj[off].index, val[off]));
378e1d83109SMatthew G. Knepley             }
37984961fc4SMatthew G. Knepley           }
38084961fc4SMatthew G. Knepley         }
38184961fc4SMatthew G. Knepley       }
38284961fc4SMatthew G. Knepley       /* Symmetrize the graph */
383*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreate(PETSC_COMM_SELF, &G));
384*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size]));
385*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetUp(G));
3869852e123SBarry Smith       for (p = 0, off = 0; p < size; ++p) {
387e1d83109SMatthew G. Knepley         for (c = 0; c < Nc[p]; ++c) {
388e1d83109SMatthew G. Knepley           const PetscInt r = Noff[p]+c;
389e1d83109SMatthew G. Knepley           PetscInt       n;
390e1d83109SMatthew G. Knepley 
391e1d83109SMatthew G. Knepley           for (n = 0; n < N[r]; ++n, ++off) {
392e1d83109SMatthew G. Knepley             const PetscInt    q = Noff[adj[off].rank] + adj[off].index;
393e1d83109SMatthew G. Knepley             const PetscScalar o = val[off] ? 1.0 : 0.0;
39484961fc4SMatthew G. Knepley 
395*5f80ce2aSJacob Faibussowitsch             CHKERRQ(MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES));
396*5f80ce2aSJacob Faibussowitsch             CHKERRQ(MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES));
39784961fc4SMatthew G. Knepley           }
39884961fc4SMatthew G. Knepley         }
399e1d83109SMatthew G. Knepley       }
400*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY));
401*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY));
40284961fc4SMatthew G. Knepley 
403*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscBTCreate(Noff[size], &seenProcs));
404*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscBTMemzero(Noff[size], seenProcs));
405*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscBTCreate(Noff[size], &flippedProcs));
406*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscBTMemzero(Noff[size], flippedProcs));
407*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(Noff[size], &procFIFO));
40884961fc4SMatthew G. Knepley       pTop = pBottom = 0;
4099852e123SBarry Smith       for (p = 0; p < Noff[size]; ++p) {
41084961fc4SMatthew G. Knepley         if (PetscBTLookup(seenProcs, p)) continue;
41184961fc4SMatthew G. Knepley         /* Initialize FIFO with next proc */
41284961fc4SMatthew G. Knepley         procFIFO[pBottom++] = p;
413*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscBTSet(seenProcs, p));
41484961fc4SMatthew G. Knepley         /* Consider each proc in FIFO */
41584961fc4SMatthew G. Knepley         while (pTop < pBottom) {
41684961fc4SMatthew G. Knepley           const PetscScalar *ornt;
41784961fc4SMatthew G. Knepley           const PetscInt    *neighbors;
418e1d83109SMatthew G. Knepley           PetscInt           proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n;
41984961fc4SMatthew G. Knepley 
42084961fc4SMatthew G. Knepley           proc     = procFIFO[pTop++];
42184961fc4SMatthew G. Knepley           flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0;
422*5f80ce2aSJacob Faibussowitsch           CHKERRQ(MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt));
42384961fc4SMatthew G. Knepley           /* Loop over neighboring procs */
42484961fc4SMatthew G. Knepley           for (n = 0; n < numNeighbors; ++n) {
42584961fc4SMatthew G. Knepley             nproc    = neighbors[n];
42684961fc4SMatthew G. Knepley             mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1;
42784961fc4SMatthew G. Knepley             seen     = PetscBTLookup(seenProcs, nproc);
42884961fc4SMatthew G. Knepley             flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0;
42984961fc4SMatthew G. Knepley 
43084961fc4SMatthew G. Knepley             if (mismatch ^ (flippedA ^ flippedB)) {
4312c71b3e2SJacob Faibussowitsch               PetscCheckFalse(seen,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Previously seen procs %d and %d do not match: Fault mesh is non-orientable", proc, nproc);
43284961fc4SMatthew G. Knepley               if (!flippedB) {
433*5f80ce2aSJacob Faibussowitsch                 CHKERRQ(PetscBTSet(flippedProcs, nproc));
43484961fc4SMatthew G. Knepley               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
4352c71b3e2SJacob Faibussowitsch             } else PetscCheckFalse(mismatch && flippedA && flippedB,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
43684961fc4SMatthew G. Knepley             if (!seen) {
43784961fc4SMatthew G. Knepley               procFIFO[pBottom++] = nproc;
438*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscBTSet(seenProcs, nproc));
43984961fc4SMatthew G. Knepley             }
44084961fc4SMatthew G. Knepley           }
44184961fc4SMatthew G. Knepley         }
44284961fc4SMatthew G. Knepley       }
443*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(procFIFO));
444*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&G));
445*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree2(adj, val));
446*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscBTDestroy(&seenProcs));
44784961fc4SMatthew G. Knepley     }
448e1d83109SMatthew G. Knepley     /* Scatter flip flags */
449e1d83109SMatthew G. Knepley     {
450e1d83109SMatthew G. Knepley       PetscBool *flips = NULL;
451e1d83109SMatthew G. Knepley 
452dd400576SPatrick Sanan       if (rank == 0) {
453*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscMalloc1(Noff[size], &flips));
4549852e123SBarry Smith         for (p = 0; p < Noff[size]; ++p) {
455e1d83109SMatthew G. Knepley           flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE;
456*5f80ce2aSJacob Faibussowitsch           if (flg && flips[p]) CHKERRQ(PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p));
457e1d83109SMatthew G. Knepley         }
4589852e123SBarry Smith         for (p = 0; p < size; ++p) {
459e1d83109SMatthew G. Knepley           displs[p+1] = displs[p] + Nc[p];
460e1d83109SMatthew G. Knepley         }
461e1d83109SMatthew G. Knepley       }
462*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, numComponents, MPIU_BOOL, 0, comm));
463*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(flips));
46484961fc4SMatthew G. Knepley     }
465*5f80ce2aSJacob Faibussowitsch     if (rank == 0) CHKERRQ(PetscBTDestroy(&flippedProcs));
466*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(N));
467*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree4(recvcounts, displs, Nc, Noff));
468*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree2(nrankComp, match));
469e1d83109SMatthew G. Knepley 
470e1d83109SMatthew G. Knepley     /* Decide whether to flip cells in each component */
471*5f80ce2aSJacob Faibussowitsch     for (c = 0; c < cEnd-cStart; ++c) {if (flipped[cellComp[c]]) CHKERRQ(PetscBTNegate(flippedCells, c));}
472*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(flipped));
47384961fc4SMatthew G. Knepley   }
47484961fc4SMatthew G. Knepley   if (flg) {
47584961fc4SMatthew G. Knepley     PetscViewer v;
47684961fc4SMatthew G. Knepley 
477*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIGetStdout(comm, &v));
478*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPushSynchronized(v));
479*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank));
480*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBTView(cEnd-cStart, flippedCells, v));
481*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerFlush(v));
482*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPopSynchronized(v));
48384961fc4SMatthew G. Knepley   }
48484961fc4SMatthew G. Knepley   /* Reverse flipped cells in the mesh */
48584961fc4SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
4860d366550SMatthew G. Knepley     if (PetscBTLookup(flippedCells, c-cStart)) {
487*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexOrientPoint(dm, c, -1));
4880d366550SMatthew G. Knepley     }
48984961fc4SMatthew G. Knepley   }
490*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBTDestroy(&seenCells));
491*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBTDestroy(&flippedCells));
492*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBTDestroy(&seenFaces));
493*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(numNeighbors, neighbors));
494*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(rorntComp, lorntComp));
495*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree3(faceFIFO, cellComp, faceComp));
49684961fc4SMatthew G. Knepley   PetscFunctionReturn(0);
49784961fc4SMatthew G. Knepley }
498