xref: /petsc/src/dm/impls/plex/plexorient.c (revision 1dca8a0504492127e77eac64bc165d7372dd6d63)
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);
279566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, p, &ct));
28b5a892a1SMatthew G. Knepley   arr  = DMPolytopeTypeGetArrangment(ct, o);
299566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
309566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCone(dm, p, &cone));
319566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeOrientation(dm, p, &ornt));
329566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, coneSize, MPIU_INT, &newcone));
339566063dSJacob Faibussowitsch   PetscCall(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 
389566063dSJacob Faibussowitsch     PetscCall(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]]);
4263a3b9bcSJacob Faibussowitsch     PetscCheck(!newornt[c] || !(newornt[c] >= nO || newornt[c] < -nO),PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid orientation %" PetscInt_FMT " not in [%" PetscInt_FMT ",%" PetscInt_FMT ") for %s %" PetscInt_FMT, newornt[c], -nO, nO, DMPolytopeTypes[ft], cone[c]);
4327d0e99aSVaclav Hapla   }
449566063dSJacob Faibussowitsch   PetscCall(DMPlexSetCone(dm, p, newcone));
459566063dSJacob Faibussowitsch   PetscCall(DMPlexSetConeOrientation(dm, p, newornt));
469566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm, coneSize, MPIU_INT, &newcone));
479566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm, coneSize, MPIU_INT, &newornt));
48b5a892a1SMatthew G. Knepley   /* Update orientation of this point in the support points */
499566063dSJacob Faibussowitsch   PetscCall(DMPlexGetSupportSize(dm, p, &supportSize));
509566063dSJacob Faibussowitsch   PetscCall(DMPlexGetSupport(dm, p, &support));
51b5a892a1SMatthew G. Knepley   for (s = 0; s < supportSize; ++s) {
529566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeSize(dm, support[s], &coneSize));
539566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCone(dm, support[s], &cone));
549566063dSJacob Faibussowitsch     PetscCall(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);
619566063dSJacob Faibussowitsch       PetscCall(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)++];
809566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
819566063dSJacob Faibussowitsch   PetscCall(DMPlexGetSupportSize(dm, face, &supportSize));
829566063dSJacob Faibussowitsch   PetscCall(DMPlexGetSupport(dm, face, &support));
837b310ce2SMatthew G. Knepley   if (supportSize < 2) PetscFunctionReturn(0);
8463a3b9bcSJacob Faibussowitsch   PetscCheck(supportSize == 2,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %" PetscInt_FMT, 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 
909566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeSize(dm, support[0], &coneSizeA));
919566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeSize(dm, support[1], &coneSizeB));
929566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCone(dm, support[0], &coneA));
939566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCone(dm, support[1], &coneB));
949566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeOrientation(dm, support[0], &coneOA));
959566063dSJacob Faibussowitsch   PetscCall(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];
999566063dSJacob Faibussowitsch       PetscCall(PetscBTSet(seenFaces, coneA[c]-fStart));
1007b310ce2SMatthew G. Knepley     }
1017b310ce2SMatthew G. Knepley     if (coneA[c] == face) posA = c;
10263a3b9bcSJacob Faibussowitsch     PetscCheck(*fBottom <= fEnd-fStart,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %" PetscInt_FMT " was pushed exceeding capacity %" PetscInt_FMT " > %" PetscInt_FMT, coneA[c], *fBottom, fEnd-fStart);
1037b310ce2SMatthew G. Knepley   }
10463a3b9bcSJacob Faibussowitsch   PetscCheck(posA >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " could not be located in cell %" PetscInt_FMT, 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];
1089566063dSJacob Faibussowitsch       PetscCall(PetscBTSet(seenFaces, coneB[c]-fStart));
1097b310ce2SMatthew G. Knepley     }
1107b310ce2SMatthew G. Knepley     if (coneB[c] == face) posB = c;
11163a3b9bcSJacob Faibussowitsch     PetscCheck(*fBottom <= fEnd-fStart,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %" PetscInt_FMT " was pushed exceeding capacity %" PetscInt_FMT " > %" PetscInt_FMT, coneA[c], *fBottom, fEnd-fStart);
1127b310ce2SMatthew G. Knepley   }
11363a3b9bcSJacob Faibussowitsch   PetscCheck(posB >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " could not be located in cell %" PetscInt_FMT, 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)) {
12263a3b9bcSJacob Faibussowitsch     PetscCheck(!seenA || !seenB,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Previously seen cells %" PetscInt_FMT " and %" PetscInt_FMT " do not match: Fault mesh is non-orientable", support[0], support[1]);
1237b310ce2SMatthew G. Knepley     if (!seenA && !flippedA) {
1249566063dSJacob Faibussowitsch       PetscCall(PetscBTSet(flippedCells, support[0]-cStart));
1257b310ce2SMatthew G. Knepley     } else if (!seenB && !flippedB) {
1269566063dSJacob Faibussowitsch       PetscCall(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");
128*1dca8a05SBarry Smith   } else PetscCheck(!mismatch || !flippedA || !flippedB,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
1299566063dSJacob Faibussowitsch   PetscCall(PetscBTSet(seenCells, support[0]-cStart));
1309566063dSJacob Faibussowitsch   PetscCall(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;
1543d6051bcSMatthew G. Knepley   PetscInt          *numNeighbors, **neighbors, *locSupport = NULL;
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;
1659566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
1669566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1679566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
1689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view", &flg));
1699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view_synchronized", &flg2));
1709566063dSJacob Faibussowitsch   PetscCall(DMGetPointSF(dm, &sf));
1719566063dSJacob Faibussowitsch   PetscCall(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   */
1819566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1829566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(dm, &h));
1839566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, h,   &cStart, &cEnd));
1849566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, h+1, &fStart, &fEnd));
1859566063dSJacob Faibussowitsch   PetscCall(PetscBTCreate(cEnd - cStart, &seenCells));
1869566063dSJacob Faibussowitsch   PetscCall(PetscBTMemzero(cEnd - cStart, seenCells));
1879566063dSJacob Faibussowitsch   PetscCall(PetscBTCreate(cEnd - cStart, &flippedCells));
1889566063dSJacob Faibussowitsch   PetscCall(PetscBTMemzero(cEnd - cStart, flippedCells));
1899566063dSJacob Faibussowitsch   PetscCall(PetscBTCreate(fEnd - fStart, &seenFaces));
1909566063dSJacob Faibussowitsch   PetscCall(PetscBTMemzero(fEnd - fStart, seenFaces));
1919566063dSJacob Faibussowitsch   PetscCall(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;
2239566063dSJacob Faibussowitsch       PetscCall(DMPlexGetConeSize(dm, cell, &coneSize));
2249566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(dm, cell, &cone));
22584961fc4SMatthew G. Knepley       for (c = 0; c < coneSize; ++c) {
22684961fc4SMatthew G. Knepley         faceFIFO[fBottom++] = cone[c];
2279566063dSJacob Faibussowitsch         PetscCall(PetscBTSet(seenFaces, cone[c]-fStart));
22884961fc4SMatthew G. Knepley       }
2299566063dSJacob Faibussowitsch       PetscCall(PetscBTSet(seenCells, cell-cStart));
23084961fc4SMatthew G. Knepley     }
23184961fc4SMatthew G. Knepley     /* Consider each face in FIFO */
23284961fc4SMatthew G. Knepley     while (fTop < fBottom) {
2339566063dSJacob Faibussowitsch       PetscCall(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 */
2439566063dSJacob Faibussowitsch     PetscCall(PetscBTMemzero(fEnd - fStart, seenFaces));
2449566063dSJacob Faibussowitsch     PetscCall(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 
2519566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIGetStdout(comm, &v));
2529566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(v));
2539566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank));
2549566063dSJacob Faibussowitsch     PetscCall(PetscBTView(cEnd-cStart, flippedCells, v));
2559566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(v));
2569566063dSJacob Faibussowitsch     PetscCall(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) {
2603d6051bcSMatthew G. Knepley     PetscInt maxSupportSize, neighbor;
2613d6051bcSMatthew G. Knepley 
262e1d83109SMatthew G. Knepley     /* Store orientations of boundary faces*/
2633d6051bcSMatthew G. Knepley     PetscCall(DMPlexGetMaxSizes(dm, NULL, &maxSupportSize));
2643d6051bcSMatthew G. Knepley     PetscCall(PetscCalloc3(numRoots, &rorntComp, numRoots, &lorntComp, maxSupportSize, &locSupport));
26584961fc4SMatthew G. Knepley     for (face = fStart; face < fEnd; ++face) {
266e1d83109SMatthew G. Knepley       const PetscInt *cone, *support, *ornt;
2673d6051bcSMatthew G. Knepley       PetscInt        coneSize, supportSize, Ns = 0, s, l;
268e1d83109SMatthew G. Knepley 
2699566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupportSize(dm, face, &supportSize));
2703d6051bcSMatthew G. Knepley       /* Ignore overlapping cells */
2719566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupport(dm, face, &support));
2723d6051bcSMatthew G. Knepley       for (s = 0; s < supportSize; ++s) {
2733d6051bcSMatthew G. Knepley         PetscCall(PetscFindInt(support[s], numLeaves, lpoints, &l));
2743d6051bcSMatthew G. Knepley         if (l >= 0) continue;
2753d6051bcSMatthew G. Knepley         locSupport[Ns++] = support[s];
2763d6051bcSMatthew G. Knepley       }
2773d6051bcSMatthew G. Knepley       if (Ns != 1) continue;
2783d6051bcSMatthew G. Knepley       neighbor = locSupport[0];
2793d6051bcSMatthew G. Knepley       PetscCall(DMPlexGetCone(dm, neighbor, &cone));
2803d6051bcSMatthew G. Knepley       PetscCall(DMPlexGetConeSize(dm, neighbor, &coneSize));
2813d6051bcSMatthew G. Knepley       PetscCall(DMPlexGetConeOrientation(dm, neighbor, &ornt));
28284961fc4SMatthew G. Knepley       for (c = 0; c < coneSize; ++c) if (cone[c] == face) break;
28384961fc4SMatthew G. Knepley       if (dim == 1) {
28484961fc4SMatthew G. Knepley         /* Use cone position instead, shifted to -1 or 1 */
2853d6051bcSMatthew G. Knepley         if (PetscBTLookup(flippedCells, neighbor-cStart)) rorntComp[face].rank = 1-c*2;
2866e1eb25cSMatthew G. Knepley         else                                              rorntComp[face].rank = c*2-1;
28784961fc4SMatthew G. Knepley       } else {
2883d6051bcSMatthew G. Knepley         if (PetscBTLookup(flippedCells, neighbor-cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 :  1;
289e1d83109SMatthew G. Knepley         else                                              rorntComp[face].rank = ornt[c] < 0 ?  1 : -1;
29084961fc4SMatthew G. Knepley       }
291e1d83109SMatthew G. Knepley       rorntComp[face].index = faceComp[face-fStart];
29284961fc4SMatthew G. Knepley     }
293e1d83109SMatthew G. Knepley     /* Communicate boundary edge orientations */
2949566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf, MPIU_2INT, rorntComp, lorntComp,MPI_REPLACE));
2959566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf, MPIU_2INT, rorntComp, lorntComp,MPI_REPLACE));
296e1d83109SMatthew G. Knepley   }
297e1d83109SMatthew G. Knepley   /* Get process adjacency */
2989566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors));
299a17aa47eSToby Isaac   viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm));
3009566063dSJacob Faibussowitsch   if (flg2) PetscCall(PetscViewerASCIIPushSynchronized(viewer));
3019566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&selfviewer));
302e1d83109SMatthew G. Knepley   for (comp = 0; comp < numComponents; ++comp) {
303e1d83109SMatthew G. Knepley     PetscInt  l, n;
30484961fc4SMatthew G. Knepley 
305e1d83109SMatthew G. Knepley     numNeighbors[comp] = 0;
3069566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp]));
307e1d83109SMatthew G. Knepley     /* I know this is p^2 time in general, but for bounded degree its alright */
308e1d83109SMatthew G. Knepley     for (l = 0; l < numLeaves; ++l) {
309e1d83109SMatthew G. Knepley       const PetscInt face = lpoints[l];
310e1d83109SMatthew G. Knepley 
311e1d83109SMatthew G. Knepley       /* Find a representative face (edge) separating pairs of procs */
3123d6051bcSMatthew G. Knepley       if ((face >= fStart) && (face < fEnd) && (faceComp[face-fStart] == comp) && rorntComp[face].rank) {
313269a49d5SMatthew G. Knepley         const PetscInt rrank = rpoints[l].rank;
314e1d83109SMatthew G. Knepley         const PetscInt rcomp = lorntComp[face].index;
315e1d83109SMatthew G. Knepley 
316269a49d5SMatthew G. Knepley         for (n = 0; n < numNeighbors[comp]; ++n) if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break;
317e1d83109SMatthew G. Knepley         if (n >= numNeighbors[comp]) {
318e1d83109SMatthew G. Knepley           PetscInt supportSize;
319e1d83109SMatthew G. Knepley 
3209566063dSJacob Faibussowitsch           PetscCall(DMPlexGetSupportSize(dm, face, &supportSize));
32163a3b9bcSJacob Faibussowitsch           PetscCheck(supportSize == 1,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %" PetscInt_FMT, supportSize);
32263a3b9bcSJacob Faibussowitsch           if (flg) PetscCall(PetscViewerASCIIPrintf(selfviewer, "[%d]: component %d, Found representative leaf %" PetscInt_FMT " (face %" PetscInt_FMT ") connecting to face %" PetscInt_FMT " on (%" PetscInt_FMT ", %" PetscInt_FMT ") with orientation %" PetscInt_FMT "\n", rank, comp, l, face, rpoints[l].index, rrank, rcomp, lorntComp[face].rank));
323e1d83109SMatthew G. Knepley           neighbors[comp][numNeighbors[comp]++] = l;
324e1d83109SMatthew G. Knepley         }
325e1d83109SMatthew G. Knepley       }
326e1d83109SMatthew G. Knepley     }
327e1d83109SMatthew G. Knepley     totNeighbors += numNeighbors[comp];
328e1d83109SMatthew G. Knepley   }
3299566063dSJacob Faibussowitsch   PetscCall(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&selfviewer));
3309566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
3319566063dSJacob Faibussowitsch   if (flg2) PetscCall(PetscViewerASCIIPopSynchronized(viewer));
3329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match));
333e1d83109SMatthew G. Knepley   for (comp = 0, off = 0; comp < numComponents; ++comp) {
334e1d83109SMatthew G. Knepley     PetscInt n;
335e1d83109SMatthew G. Knepley 
336e1d83109SMatthew G. Knepley     for (n = 0; n < numNeighbors[comp]; ++n, ++off) {
337e1d83109SMatthew G. Knepley       const PetscInt face = lpoints[neighbors[comp][n]];
338e1d83109SMatthew G. Knepley       const PetscInt o    = rorntComp[face].rank*lorntComp[face].rank;
339e1d83109SMatthew G. Knepley 
340e1d83109SMatthew G. Knepley       if      (o < 0) match[off] = PETSC_TRUE;
341e1d83109SMatthew G. Knepley       else if (o > 0) match[off] = PETSC_FALSE;
34263a3b9bcSJacob Faibussowitsch       else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid face %" PetscInt_FMT " (%" PetscInt_FMT ", %" PetscInt_FMT ") neighbor: %" PetscInt_FMT " comp: %d", face, rorntComp[face].rank, lorntComp[face].rank, neighbors[comp][n], comp);
343e1d83109SMatthew G. Knepley       nrankComp[off].rank  = rpoints[neighbors[comp][n]].rank;
344e1d83109SMatthew G. Knepley       nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index;
345e1d83109SMatthew G. Knepley     }
3469566063dSJacob Faibussowitsch     PetscCall(PetscFree(neighbors[comp]));
34784961fc4SMatthew G. Knepley   }
34884961fc4SMatthew G. Knepley   /* Collect the graph on 0 */
349e1d83109SMatthew G. Knepley   if (numLeaves >= 0) {
35084961fc4SMatthew G. Knepley     Mat          G;
35184961fc4SMatthew G. Knepley     PetscBT      seenProcs, flippedProcs;
35284961fc4SMatthew G. Knepley     PetscInt    *procFIFO, pTop, pBottom;
3537cadcfe8SMatthew G. Knepley     PetscInt    *N   = NULL, *Noff;
354e1d83109SMatthew G. Knepley     PetscSFNode *adj = NULL;
35584961fc4SMatthew G. Knepley     PetscBool   *val = NULL;
3567cadcfe8SMatthew G. Knepley     PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc, p, o;
3579852e123SBarry Smith     PetscMPIInt  size = 0;
35884961fc4SMatthew G. Knepley 
3599566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(numComponents, &flipped));
3609566063dSJacob Faibussowitsch     if (rank == 0) PetscCallMPI(MPI_Comm_size(comm, &size));
3619566063dSJacob Faibussowitsch     PetscCall(PetscCalloc4(size, &recvcounts, size+1, &displs, size, &Nc, size+1, &Noff));
3629566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Gather(&numComponents, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm));
3639852e123SBarry Smith     for (p = 0; p < size; ++p) {
364e1d83109SMatthew G. Knepley       displs[p+1] = displs[p] + Nc[p];
365e1d83109SMatthew G. Knepley     }
3669566063dSJacob Faibussowitsch     if (rank == 0) PetscCall(PetscMalloc1(displs[size],&N));
3679566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Gatherv(numNeighbors, numComponents, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm));
3689852e123SBarry Smith     for (p = 0, o = 0; p < size; ++p) {
369e1d83109SMatthew G. Knepley       recvcounts[p] = 0;
370e1d83109SMatthew G. Knepley       for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o];
37184961fc4SMatthew G. Knepley       displs[p+1] = displs[p] + recvcounts[p];
37284961fc4SMatthew G. Knepley     }
3739566063dSJacob Faibussowitsch     if (rank == 0) PetscCall(PetscMalloc2(displs[size], &adj, displs[size], &val));
3749566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Gatherv(nrankComp, totNeighbors, MPIU_2INT, adj, recvcounts, displs, MPIU_2INT, 0, comm));
3759566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Gatherv(match, totNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm));
3769566063dSJacob Faibussowitsch     PetscCall(PetscFree2(numNeighbors, neighbors));
377dd400576SPatrick Sanan     if (rank == 0) {
3789852e123SBarry Smith       for (p = 1; p <= size; ++p) {Noff[p] = Noff[p-1] + Nc[p-1];}
379a83982f3SMatthew G. Knepley       if (flg) {
380e1d83109SMatthew G. Knepley         PetscInt n;
381e1d83109SMatthew G. Knepley 
3829852e123SBarry Smith         for (p = 0, off = 0; p < size; ++p) {
383e1d83109SMatthew G. Knepley           for (c = 0; c < Nc[p]; ++c) {
38463a3b9bcSJacob Faibussowitsch             PetscCall(PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %" PetscInt_FMT ":\n", p, c));
385e1d83109SMatthew G. Knepley             for (n = 0; n < N[Noff[p]+c]; ++n, ++off) {
38663a3b9bcSJacob Faibussowitsch               PetscCall(PetscPrintf(PETSC_COMM_SELF, "  edge (%" PetscInt_FMT ", %" PetscInt_FMT ") (%s):\n", adj[off].rank, adj[off].index, PetscBools[val[off]]));
387e1d83109SMatthew G. Knepley             }
38884961fc4SMatthew G. Knepley           }
38984961fc4SMatthew G. Knepley         }
39084961fc4SMatthew G. Knepley       }
39184961fc4SMatthew G. Knepley       /* Symmetrize the graph */
3929566063dSJacob Faibussowitsch       PetscCall(MatCreate(PETSC_COMM_SELF, &G));
3939566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size]));
3949566063dSJacob Faibussowitsch       PetscCall(MatSetUp(G));
3959852e123SBarry Smith       for (p = 0, off = 0; p < size; ++p) {
396e1d83109SMatthew G. Knepley         for (c = 0; c < Nc[p]; ++c) {
397e1d83109SMatthew G. Knepley           const PetscInt r = Noff[p]+c;
398e1d83109SMatthew G. Knepley           PetscInt       n;
399e1d83109SMatthew G. Knepley 
400e1d83109SMatthew G. Knepley           for (n = 0; n < N[r]; ++n, ++off) {
401e1d83109SMatthew G. Knepley             const PetscInt    q = Noff[adj[off].rank] + adj[off].index;
402e1d83109SMatthew G. Knepley             const PetscScalar o = val[off] ? 1.0 : 0.0;
40384961fc4SMatthew G. Knepley 
4049566063dSJacob Faibussowitsch             PetscCall(MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES));
4059566063dSJacob Faibussowitsch             PetscCall(MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES));
40684961fc4SMatthew G. Knepley           }
40784961fc4SMatthew G. Knepley         }
408e1d83109SMatthew G. Knepley       }
4099566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY));
4109566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY));
41184961fc4SMatthew G. Knepley 
4129566063dSJacob Faibussowitsch       PetscCall(PetscBTCreate(Noff[size], &seenProcs));
4139566063dSJacob Faibussowitsch       PetscCall(PetscBTMemzero(Noff[size], seenProcs));
4149566063dSJacob Faibussowitsch       PetscCall(PetscBTCreate(Noff[size], &flippedProcs));
4159566063dSJacob Faibussowitsch       PetscCall(PetscBTMemzero(Noff[size], flippedProcs));
4169566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(Noff[size], &procFIFO));
41784961fc4SMatthew G. Knepley       pTop = pBottom = 0;
4189852e123SBarry Smith       for (p = 0; p < Noff[size]; ++p) {
41984961fc4SMatthew G. Knepley         if (PetscBTLookup(seenProcs, p)) continue;
42084961fc4SMatthew G. Knepley         /* Initialize FIFO with next proc */
42184961fc4SMatthew G. Knepley         procFIFO[pBottom++] = p;
4229566063dSJacob Faibussowitsch         PetscCall(PetscBTSet(seenProcs, p));
42384961fc4SMatthew G. Knepley         /* Consider each proc in FIFO */
42484961fc4SMatthew G. Knepley         while (pTop < pBottom) {
42584961fc4SMatthew G. Knepley           const PetscScalar *ornt;
42684961fc4SMatthew G. Knepley           const PetscInt    *neighbors;
427e1d83109SMatthew G. Knepley           PetscInt           proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n;
42884961fc4SMatthew G. Knepley 
42984961fc4SMatthew G. Knepley           proc     = procFIFO[pTop++];
43084961fc4SMatthew G. Knepley           flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0;
4319566063dSJacob Faibussowitsch           PetscCall(MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt));
43284961fc4SMatthew G. Knepley           /* Loop over neighboring procs */
43384961fc4SMatthew G. Knepley           for (n = 0; n < numNeighbors; ++n) {
43484961fc4SMatthew G. Knepley             nproc    = neighbors[n];
43584961fc4SMatthew G. Knepley             mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1;
43684961fc4SMatthew G. Knepley             seen     = PetscBTLookup(seenProcs, nproc);
43784961fc4SMatthew G. Knepley             flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0;
43884961fc4SMatthew G. Knepley 
43984961fc4SMatthew G. Knepley             if (mismatch ^ (flippedA ^ flippedB)) {
44063a3b9bcSJacob Faibussowitsch               PetscCheck(!seen,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Previously seen procs %" PetscInt_FMT " and %" PetscInt_FMT " do not match: Fault mesh is non-orientable", proc, nproc);
44184961fc4SMatthew G. Knepley               if (!flippedB) {
4429566063dSJacob Faibussowitsch                 PetscCall(PetscBTSet(flippedProcs, nproc));
44384961fc4SMatthew G. Knepley               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
444*1dca8a05SBarry Smith             } else PetscCheck(!mismatch || !flippedA || !flippedB,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
44584961fc4SMatthew G. Knepley             if (!seen) {
44684961fc4SMatthew G. Knepley               procFIFO[pBottom++] = nproc;
4479566063dSJacob Faibussowitsch               PetscCall(PetscBTSet(seenProcs, nproc));
44884961fc4SMatthew G. Knepley             }
44984961fc4SMatthew G. Knepley           }
45084961fc4SMatthew G. Knepley         }
45184961fc4SMatthew G. Knepley       }
4529566063dSJacob Faibussowitsch       PetscCall(PetscFree(procFIFO));
4539566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&G));
4549566063dSJacob Faibussowitsch       PetscCall(PetscFree2(adj, val));
4559566063dSJacob Faibussowitsch       PetscCall(PetscBTDestroy(&seenProcs));
45684961fc4SMatthew G. Knepley     }
457e1d83109SMatthew G. Knepley     /* Scatter flip flags */
458e1d83109SMatthew G. Knepley     {
459e1d83109SMatthew G. Knepley       PetscBool *flips = NULL;
460e1d83109SMatthew G. Knepley 
461dd400576SPatrick Sanan       if (rank == 0) {
4629566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(Noff[size], &flips));
4639852e123SBarry Smith         for (p = 0; p < Noff[size]; ++p) {
464e1d83109SMatthew G. Knepley           flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE;
4659566063dSJacob Faibussowitsch           if (flg && flips[p]) PetscCall(PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p));
466e1d83109SMatthew G. Knepley         }
4679852e123SBarry Smith         for (p = 0; p < size; ++p) {
468e1d83109SMatthew G. Knepley           displs[p+1] = displs[p] + Nc[p];
469e1d83109SMatthew G. Knepley         }
470e1d83109SMatthew G. Knepley       }
4719566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, numComponents, MPIU_BOOL, 0, comm));
4729566063dSJacob Faibussowitsch       PetscCall(PetscFree(flips));
47384961fc4SMatthew G. Knepley     }
4749566063dSJacob Faibussowitsch     if (rank == 0) PetscCall(PetscBTDestroy(&flippedProcs));
4759566063dSJacob Faibussowitsch     PetscCall(PetscFree(N));
4769566063dSJacob Faibussowitsch     PetscCall(PetscFree4(recvcounts, displs, Nc, Noff));
4779566063dSJacob Faibussowitsch     PetscCall(PetscFree2(nrankComp, match));
478e1d83109SMatthew G. Knepley 
479e1d83109SMatthew G. Knepley     /* Decide whether to flip cells in each component */
4809566063dSJacob Faibussowitsch     for (c = 0; c < cEnd-cStart; ++c) {if (flipped[cellComp[c]]) PetscCall(PetscBTNegate(flippedCells, c));}
4819566063dSJacob Faibussowitsch     PetscCall(PetscFree(flipped));
48284961fc4SMatthew G. Knepley   }
48384961fc4SMatthew G. Knepley   if (flg) {
48484961fc4SMatthew G. Knepley     PetscViewer v;
48584961fc4SMatthew G. Knepley 
4869566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIGetStdout(comm, &v));
4879566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(v));
4889566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank));
4899566063dSJacob Faibussowitsch     PetscCall(PetscBTView(cEnd-cStart, flippedCells, v));
4909566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(v));
4919566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(v));
49284961fc4SMatthew G. Knepley   }
49384961fc4SMatthew G. Knepley   /* Reverse flipped cells in the mesh */
49484961fc4SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
4950d366550SMatthew G. Knepley     if (PetscBTLookup(flippedCells, c-cStart)) {
4969566063dSJacob Faibussowitsch       PetscCall(DMPlexOrientPoint(dm, c, -1));
4970d366550SMatthew G. Knepley     }
49884961fc4SMatthew G. Knepley   }
4999566063dSJacob Faibussowitsch   PetscCall(PetscBTDestroy(&seenCells));
5009566063dSJacob Faibussowitsch   PetscCall(PetscBTDestroy(&flippedCells));
5019566063dSJacob Faibussowitsch   PetscCall(PetscBTDestroy(&seenFaces));
5029566063dSJacob Faibussowitsch   PetscCall(PetscFree2(numNeighbors, neighbors));
5033d6051bcSMatthew G. Knepley   PetscCall(PetscFree3(rorntComp, lorntComp, locSupport));
5049566063dSJacob Faibussowitsch   PetscCall(PetscFree3(faceFIFO, cellComp, faceComp));
50584961fc4SMatthew G. Knepley   PetscFunctionReturn(0);
50684961fc4SMatthew G. Knepley }
507