xref: /petsc/src/dm/impls/plex/plexorient.c (revision 9852e1235e6de1345fdd93ba94fdd27144d16762)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
284961fc4SMatthew G. Knepley #include <petscsf.h>
384961fc4SMatthew G. Knepley 
484961fc4SMatthew G. Knepley /*@
584961fc4SMatthew G. Knepley   DMPlexReverseCell - Give a mesh cell the opposite orientation
684961fc4SMatthew G. Knepley 
784961fc4SMatthew G. Knepley   Input Parameters:
884961fc4SMatthew G. Knepley + dm   - The DM
984961fc4SMatthew G. Knepley - cell - The cell number
1084961fc4SMatthew G. Knepley 
1184961fc4SMatthew G. Knepley   Note: The modification of the DM is done in-place.
1284961fc4SMatthew G. Knepley 
1384961fc4SMatthew G. Knepley   Level: advanced
1484961fc4SMatthew G. Knepley 
1584961fc4SMatthew G. Knepley .seealso: DMPlexOrient(), DMCreate(), DMPLEX
1684961fc4SMatthew G. Knepley @*/
1784961fc4SMatthew G. Knepley PetscErrorCode DMPlexReverseCell(DM dm, PetscInt cell)
1884961fc4SMatthew G. Knepley {
1984961fc4SMatthew G. Knepley   /* Note that the reverse orientation ro of a face with orientation o is:
2084961fc4SMatthew G. Knepley 
2184961fc4SMatthew G. Knepley        ro = o >= 0 ? -(faceSize - o) : faceSize + o
2284961fc4SMatthew G. Knepley 
2384961fc4SMatthew G. Knepley      where faceSize is the size of the cone for the face.
2484961fc4SMatthew G. Knepley   */
2584961fc4SMatthew G. Knepley   const PetscInt *cone,    *coneO, *support;
2684961fc4SMatthew G. Knepley   PetscInt       *revcone, *revconeO;
2784961fc4SMatthew G. Knepley   PetscInt        maxConeSize, coneSize, supportSize, faceSize, cp, sp;
2884961fc4SMatthew G. Knepley   PetscErrorCode  ierr;
2984961fc4SMatthew G. Knepley 
3084961fc4SMatthew G. Knepley   PetscFunctionBegin;
3184961fc4SMatthew G. Knepley   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
3284961fc4SMatthew G. Knepley   ierr = DMGetWorkArray(dm, maxConeSize, PETSC_INT, &revcone);CHKERRQ(ierr);
3384961fc4SMatthew G. Knepley   ierr = DMGetWorkArray(dm, maxConeSize, PETSC_INT, &revconeO);CHKERRQ(ierr);
3484961fc4SMatthew G. Knepley   /* Reverse cone, and reverse orientations of faces */
3584961fc4SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
3684961fc4SMatthew G. Knepley   ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
3784961fc4SMatthew G. Knepley   ierr = DMPlexGetConeOrientation(dm, cell, &coneO);CHKERRQ(ierr);
3884961fc4SMatthew G. Knepley   for (cp = 0; cp < coneSize; ++cp) {
3984961fc4SMatthew G. Knepley     const PetscInt rcp = coneSize-cp-1;
4084961fc4SMatthew G. Knepley 
4184961fc4SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, cone[rcp], &faceSize);CHKERRQ(ierr);
4284961fc4SMatthew G. Knepley     revcone[cp]  = cone[rcp];
4384961fc4SMatthew G. Knepley     revconeO[cp] = coneO[rcp] >= 0 ? -(faceSize-coneO[rcp]) : faceSize+coneO[rcp];
4484961fc4SMatthew G. Knepley   }
4584961fc4SMatthew G. Knepley   ierr = DMPlexSetCone(dm, cell, revcone);CHKERRQ(ierr);
4684961fc4SMatthew G. Knepley   ierr = DMPlexSetConeOrientation(dm, cell, revconeO);CHKERRQ(ierr);
4784961fc4SMatthew G. Knepley   /* Reverse orientation of this cell in the support hypercells */
4884961fc4SMatthew G. Knepley   faceSize = coneSize;
4984961fc4SMatthew G. Knepley   ierr = DMPlexGetSupportSize(dm, cell, &supportSize);CHKERRQ(ierr);
5084961fc4SMatthew G. Knepley   ierr = DMPlexGetSupport(dm, cell, &support);CHKERRQ(ierr);
5184961fc4SMatthew G. Knepley   for (sp = 0; sp < supportSize; ++sp) {
5284961fc4SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, support[sp], &coneSize);CHKERRQ(ierr);
5384961fc4SMatthew G. Knepley     ierr = DMPlexGetCone(dm, support[sp], &cone);CHKERRQ(ierr);
5484961fc4SMatthew G. Knepley     ierr = DMPlexGetConeOrientation(dm, support[sp], &coneO);CHKERRQ(ierr);
5584961fc4SMatthew G. Knepley     for (cp = 0; cp < coneSize; ++cp) {
5684961fc4SMatthew G. Knepley       if (cone[cp] != cell) continue;
5784961fc4SMatthew G. Knepley       ierr = DMPlexInsertConeOrientation(dm, support[sp], cp, coneO[cp] >= 0 ? -(faceSize-coneO[cp]) : faceSize+coneO[cp]);CHKERRQ(ierr);
5884961fc4SMatthew G. Knepley     }
5984961fc4SMatthew G. Knepley   }
6084961fc4SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, maxConeSize, PETSC_INT, &revcone);CHKERRQ(ierr);
6184961fc4SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, maxConeSize, PETSC_INT, &revconeO);CHKERRQ(ierr);
6284961fc4SMatthew G. Knepley   PetscFunctionReturn(0);
6384961fc4SMatthew G. Knepley }
6484961fc4SMatthew G. Knepley 
657b310ce2SMatthew G. Knepley /*
667b310ce2SMatthew G. Knepley   - Checks face match
677b310ce2SMatthew G. Knepley     - Flips non-matching
687b310ce2SMatthew G. Knepley   - Inserts faces of support cells in FIFO
697b310ce2SMatthew G. Knepley */
707b2de0fdSMatthew 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)
717b310ce2SMatthew G. Knepley {
727b310ce2SMatthew G. Knepley   const PetscInt *support, *coneA, *coneB, *coneOA, *coneOB;
737b310ce2SMatthew G. Knepley   PetscInt        supportSize, coneSizeA, coneSizeB, posA = -1, posB = -1;
747b2de0fdSMatthew G. Knepley   PetscInt        face, dim, seenA, flippedA, seenB, flippedB, mismatch, c;
757b310ce2SMatthew G. Knepley   PetscErrorCode  ierr;
767b310ce2SMatthew G. Knepley 
777b310ce2SMatthew G. Knepley   PetscFunctionBegin;
787b310ce2SMatthew G. Knepley   face = faceFIFO[(*fTop)++];
797b310ce2SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
807b310ce2SMatthew G. Knepley   ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
817b310ce2SMatthew G. Knepley   ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr);
827b310ce2SMatthew G. Knepley   if (supportSize < 2) PetscFunctionReturn(0);
837b310ce2SMatthew G. Knepley   if (supportSize != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %d", supportSize);
847b310ce2SMatthew G. Knepley   seenA    = PetscBTLookup(seenCells,    support[0]-cStart);
857b310ce2SMatthew G. Knepley   flippedA = PetscBTLookup(flippedCells, support[0]-cStart) ? 1 : 0;
867b310ce2SMatthew G. Knepley   seenB    = PetscBTLookup(seenCells,    support[1]-cStart);
877b310ce2SMatthew G. Knepley   flippedB = PetscBTLookup(flippedCells, support[1]-cStart) ? 1 : 0;
887b310ce2SMatthew G. Knepley 
897b310ce2SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, support[0], &coneSizeA);CHKERRQ(ierr);
907b310ce2SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, support[1], &coneSizeB);CHKERRQ(ierr);
917b310ce2SMatthew G. Knepley   ierr = DMPlexGetCone(dm, support[0], &coneA);CHKERRQ(ierr);
927b310ce2SMatthew G. Knepley   ierr = DMPlexGetCone(dm, support[1], &coneB);CHKERRQ(ierr);
937b310ce2SMatthew G. Knepley   ierr = DMPlexGetConeOrientation(dm, support[0], &coneOA);CHKERRQ(ierr);
947b310ce2SMatthew G. Knepley   ierr = DMPlexGetConeOrientation(dm, support[1], &coneOB);CHKERRQ(ierr);
957b310ce2SMatthew G. Knepley   for (c = 0; c < coneSizeA; ++c) {
967b310ce2SMatthew G. Knepley     if (!PetscBTLookup(seenFaces, coneA[c]-fStart)) {
977b310ce2SMatthew G. Knepley       faceFIFO[(*fBottom)++] = coneA[c];
987b310ce2SMatthew G. Knepley       ierr = PetscBTSet(seenFaces, coneA[c]-fStart);CHKERRQ(ierr);
997b310ce2SMatthew G. Knepley     }
1007b310ce2SMatthew G. Knepley     if (coneA[c] == face) posA = c;
1017b310ce2SMatthew 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);
1027b310ce2SMatthew G. Knepley   }
1037b310ce2SMatthew 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]);
1047b310ce2SMatthew G. Knepley   for (c = 0; c < coneSizeB; ++c) {
1057b310ce2SMatthew G. Knepley     if (!PetscBTLookup(seenFaces, coneB[c]-fStart)) {
1067b310ce2SMatthew G. Knepley       faceFIFO[(*fBottom)++] = coneB[c];
1077b310ce2SMatthew G. Knepley       ierr = PetscBTSet(seenFaces, coneB[c]-fStart);CHKERRQ(ierr);
1087b310ce2SMatthew G. Knepley     }
1097b310ce2SMatthew G. Knepley     if (coneB[c] == face) posB = c;
1107b310ce2SMatthew 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);
1117b310ce2SMatthew G. Knepley   }
1127b310ce2SMatthew 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]);
1137b310ce2SMatthew G. Knepley 
1147b310ce2SMatthew G. Knepley   if (dim == 1) {
1157b310ce2SMatthew G. Knepley     mismatch = posA == posB;
1167b310ce2SMatthew G. Knepley   } else {
1177b310ce2SMatthew G. Knepley     mismatch = coneOA[posA] == coneOB[posB];
1187b310ce2SMatthew G. Knepley   }
1197b310ce2SMatthew G. Knepley 
1207b310ce2SMatthew G. Knepley   if (mismatch ^ (flippedA ^ flippedB)) {
1217b310ce2SMatthew 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]);
1227b310ce2SMatthew G. Knepley     if (!seenA && !flippedA) {
1237b310ce2SMatthew G. Knepley       ierr = PetscBTSet(flippedCells, support[0]-cStart);CHKERRQ(ierr);
1247b310ce2SMatthew G. Knepley     } else if (!seenB && !flippedB) {
1257b310ce2SMatthew G. Knepley       ierr = PetscBTSet(flippedCells, support[1]-cStart);CHKERRQ(ierr);
1267b310ce2SMatthew G. Knepley     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
1277b310ce2SMatthew 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");
1287b310ce2SMatthew G. Knepley   ierr = PetscBTSet(seenCells, support[0]-cStart);CHKERRQ(ierr);
1297b310ce2SMatthew G. Knepley   ierr = PetscBTSet(seenCells, support[1]-cStart);CHKERRQ(ierr);
1307b310ce2SMatthew G. Knepley   PetscFunctionReturn(0);
1317b310ce2SMatthew G. Knepley }
1327b310ce2SMatthew G. Knepley 
13384961fc4SMatthew G. Knepley /*@
13484961fc4SMatthew G. Knepley   DMPlexOrient - Give a consistent orientation to the input mesh
13584961fc4SMatthew G. Knepley 
13684961fc4SMatthew G. Knepley   Input Parameters:
13784961fc4SMatthew G. Knepley . dm - The DM
13884961fc4SMatthew G. Knepley 
13984961fc4SMatthew G. Knepley   Note: The orientation data for the DM are change in-place.
14084961fc4SMatthew G. Knepley $ This routine will fail for non-orientable surfaces, such as the Moebius strip.
14184961fc4SMatthew G. Knepley 
14284961fc4SMatthew G. Knepley   Level: advanced
14384961fc4SMatthew G. Knepley 
14484961fc4SMatthew G. Knepley .seealso: DMCreate(), DMPLEX
14584961fc4SMatthew G. Knepley @*/
14684961fc4SMatthew G. Knepley PetscErrorCode DMPlexOrient(DM dm)
14784961fc4SMatthew G. Knepley {
14884961fc4SMatthew G. Knepley   MPI_Comm           comm;
149e1d83109SMatthew G. Knepley   PetscSF            sf;
150e1d83109SMatthew G. Knepley   const PetscInt    *lpoints;
151e1d83109SMatthew G. Knepley   const PetscSFNode *rpoints;
152e1d83109SMatthew G. Knepley   PetscSFNode       *rorntComp = NULL, *lorntComp = NULL;
153e1d83109SMatthew G. Knepley   PetscInt          *numNeighbors, **neighbors;
154e1d83109SMatthew G. Knepley   PetscSFNode       *nrankComp;
155e1d83109SMatthew G. Knepley   PetscBool         *match, *flipped;
15684961fc4SMatthew G. Knepley   PetscBT            seenCells, flippedCells, seenFaces;
157e1d83109SMatthew G. Knepley   PetscInt          *faceFIFO, fTop, fBottom, *cellComp, *faceComp;
1587cadcfe8SMatthew G. Knepley   PetscInt           numLeaves, numRoots, dim, h, cStart, cEnd, c, cell, fStart, fEnd, face, off, totNeighbors = 0;
159a83982f3SMatthew G. Knepley   PetscMPIInt        rank, numComponents, comp = 0;
16084961fc4SMatthew G. Knepley   PetscBool          flg;
16184961fc4SMatthew G. Knepley   PetscErrorCode     ierr;
16284961fc4SMatthew G. Knepley 
16384961fc4SMatthew G. Knepley   PetscFunctionBegin;
16484961fc4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
165a83982f3SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
166c5929fdfSBarry Smith   ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view", &flg);CHKERRQ(ierr);
167e1d83109SMatthew G. Knepley   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
168e1d83109SMatthew G. Knepley   ierr = PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints);CHKERRQ(ierr);
16984961fc4SMatthew G. Knepley   /* Truth Table
17084961fc4SMatthew G. Knepley      mismatch    flips   do action   mismatch   flipA ^ flipB   action
17184961fc4SMatthew G. Knepley          F       0 flips     no         F             F           F
17284961fc4SMatthew G. Knepley          F       1 flip      yes        F             T           T
17384961fc4SMatthew G. Knepley          F       2 flips     no         T             F           T
17484961fc4SMatthew G. Knepley          T       0 flips     yes        T             T           F
17584961fc4SMatthew G. Knepley          T       1 flip      no
17684961fc4SMatthew G. Knepley          T       2 flips     yes
17784961fc4SMatthew G. Knepley   */
17884961fc4SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
17984961fc4SMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(dm, &h);CHKERRQ(ierr);
18084961fc4SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, h,   &cStart, &cEnd);CHKERRQ(ierr);
18184961fc4SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, h+1, &fStart, &fEnd);CHKERRQ(ierr);
18284961fc4SMatthew G. Knepley   ierr = PetscBTCreate(cEnd - cStart, &seenCells);CHKERRQ(ierr);
18384961fc4SMatthew G. Knepley   ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr);
18484961fc4SMatthew G. Knepley   ierr = PetscBTCreate(cEnd - cStart, &flippedCells);CHKERRQ(ierr);
18584961fc4SMatthew G. Knepley   ierr = PetscBTMemzero(cEnd - cStart, flippedCells);CHKERRQ(ierr);
18684961fc4SMatthew G. Knepley   ierr = PetscBTCreate(fEnd - fStart, &seenFaces);CHKERRQ(ierr);
18784961fc4SMatthew G. Knepley   ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr);
188e1d83109SMatthew G. Knepley   ierr = PetscCalloc3(fEnd - fStart, &faceFIFO, cEnd-cStart, &cellComp, fEnd-fStart, &faceComp);CHKERRQ(ierr);
189e1d83109SMatthew G. Knepley   /*
190e1d83109SMatthew G. Knepley    OLD STYLE
191e1d83109SMatthew G. Knepley    - Add an integer array over cells and faces (component) for connected component number
192e1d83109SMatthew G. Knepley    Foreach component
193e1d83109SMatthew G. Knepley      - Mark the initial cell as seen
194e1d83109SMatthew G. Knepley      - Process component as usual
195e1d83109SMatthew G. Knepley      - Set component for all seenCells
196e1d83109SMatthew G. Knepley      - Wipe seenCells and seenFaces (flippedCells can stay)
197e1d83109SMatthew G. Knepley    - Generate parallel adjacency for component using SF and seenFaces
198e1d83109SMatthew G. Knepley    - Collect numComponents adj data from each proc to 0
199e1d83109SMatthew G. Knepley    - Build same serial graph
200e1d83109SMatthew G. Knepley    - Use same solver
201e1d83109SMatthew G. Knepley    - Use Scatterv to to send back flipped flags for each component
202e1d83109SMatthew G. Knepley    - Negate flippedCells by component
203e1d83109SMatthew G. Knepley 
204e1d83109SMatthew G. Knepley    NEW STYLE
205e1d83109SMatthew G. Knepley    - Create the adj on each process
206e1d83109SMatthew G. Knepley    - Bootstrap to complete graph on proc 0
207e1d83109SMatthew G. Knepley   */
208e1d83109SMatthew G. Knepley   /* Loop over components */
209e1d83109SMatthew G. Knepley   for (cell = cStart; cell < cEnd; ++cell) cellComp[cell-cStart] = -1;
210e1d83109SMatthew G. Knepley   do {
211e1d83109SMatthew G. Knepley     /* Look for first unmarked cell */
212e1d83109SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) if (cellComp[cell-cStart] < 0) break;
213e1d83109SMatthew G. Knepley     if (cell >= cEnd) break;
214e1d83109SMatthew G. Knepley     /* Initialize FIFO with first cell in component */
215e1d83109SMatthew G. Knepley     {
21684961fc4SMatthew G. Knepley       const PetscInt *cone;
21784961fc4SMatthew G. Knepley       PetscInt        coneSize;
21884961fc4SMatthew G. Knepley 
219e1d83109SMatthew G. Knepley       fTop = fBottom = 0;
220e1d83109SMatthew G. Knepley       ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
221e1d83109SMatthew G. Knepley       ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
22284961fc4SMatthew G. Knepley       for (c = 0; c < coneSize; ++c) {
22384961fc4SMatthew G. Knepley         faceFIFO[fBottom++] = cone[c];
22484961fc4SMatthew G. Knepley         ierr = PetscBTSet(seenFaces, cone[c]-fStart);CHKERRQ(ierr);
22584961fc4SMatthew G. Knepley       }
22631c8331aSMatthew G. Knepley       ierr = PetscBTSet(seenCells, cell-cStart);CHKERRQ(ierr);
22784961fc4SMatthew G. Knepley     }
22884961fc4SMatthew G. Knepley     /* Consider each face in FIFO */
22984961fc4SMatthew G. Knepley     while (fTop < fBottom) {
2307b2de0fdSMatthew G. Knepley       ierr = DMPlexCheckFace_Internal(dm, faceFIFO, &fTop, &fBottom, cStart, fStart, fEnd, seenCells, flippedCells, seenFaces);CHKERRQ(ierr);
23184961fc4SMatthew G. Knepley     }
232e1d83109SMatthew G. Knepley     /* Set component for cells and faces */
233e1d83109SMatthew G. Knepley     for (cell = 0; cell < cEnd-cStart; ++cell) {
234e1d83109SMatthew G. Knepley       if (PetscBTLookup(seenCells, cell)) cellComp[cell] = comp;
235e1d83109SMatthew G. Knepley     }
236e1d83109SMatthew G. Knepley     for (face = 0; face < fEnd-fStart; ++face) {
237e1d83109SMatthew G. Knepley       if (PetscBTLookup(seenFaces, face)) faceComp[face] = comp;
238e1d83109SMatthew G. Knepley     }
239e1d83109SMatthew G. Knepley     /* Wipe seenCells and seenFaces for next component */
240e1d83109SMatthew G. Knepley     ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr);
241e1d83109SMatthew G. Knepley     ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr);
242e1d83109SMatthew G. Knepley     ++comp;
243e1d83109SMatthew G. Knepley   } while (1);
244e1d83109SMatthew G. Knepley   numComponents = comp;
24584961fc4SMatthew G. Knepley   if (flg) {
24684961fc4SMatthew G. Knepley     PetscViewer v;
24784961fc4SMatthew G. Knepley 
2487e09831bSMatthew G. Knepley     ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr);
2491575c14dSBarry Smith     ierr = PetscViewerASCIIPushSynchronized(v);CHKERRQ(ierr);
25084961fc4SMatthew G. Knepley     ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank);CHKERRQ(ierr);
25184961fc4SMatthew G. Knepley     ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr);
2524d4c343aSBarry Smith     ierr = PetscViewerFlush(v);CHKERRQ(ierr);
2531575c14dSBarry Smith     ierr = PetscViewerASCIIPopSynchronized(v);CHKERRQ(ierr);
25484961fc4SMatthew G. Knepley   }
25584961fc4SMatthew G. Knepley   /* Now all subdomains are oriented, but we need a consistent parallel orientation */
25684961fc4SMatthew G. Knepley   if (numLeaves >= 0) {
257e1d83109SMatthew G. Knepley     /* Store orientations of boundary faces*/
258e1d83109SMatthew G. Knepley     ierr = PetscCalloc2(numRoots,&rorntComp,numRoots,&lorntComp);CHKERRQ(ierr);
25984961fc4SMatthew G. Knepley     for (face = fStart; face < fEnd; ++face) {
260e1d83109SMatthew G. Knepley       const PetscInt *cone, *support, *ornt;
261e1d83109SMatthew G. Knepley       PetscInt        coneSize, supportSize;
262e1d83109SMatthew G. Knepley 
26384961fc4SMatthew G. Knepley       ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
26484961fc4SMatthew G. Knepley       if (supportSize != 1) continue;
26584961fc4SMatthew G. Knepley       ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr);
26684961fc4SMatthew G. Knepley 
26784961fc4SMatthew G. Knepley       ierr = DMPlexGetCone(dm, support[0], &cone);CHKERRQ(ierr);
26884961fc4SMatthew G. Knepley       ierr = DMPlexGetConeSize(dm, support[0], &coneSize);CHKERRQ(ierr);
26984961fc4SMatthew G. Knepley       ierr = DMPlexGetConeOrientation(dm, support[0], &ornt);CHKERRQ(ierr);
27084961fc4SMatthew G. Knepley       for (c = 0; c < coneSize; ++c) if (cone[c] == face) break;
27184961fc4SMatthew G. Knepley       if (dim == 1) {
27284961fc4SMatthew G. Knepley         /* Use cone position instead, shifted to -1 or 1 */
2736e1eb25cSMatthew G. Knepley         if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = 1-c*2;
2746e1eb25cSMatthew G. Knepley         else                                                rorntComp[face].rank = c*2-1;
27584961fc4SMatthew G. Knepley       } else {
276e1d83109SMatthew G. Knepley         if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 :  1;
277e1d83109SMatthew G. Knepley         else                                                rorntComp[face].rank = ornt[c] < 0 ?  1 : -1;
27884961fc4SMatthew G. Knepley       }
279e1d83109SMatthew G. Knepley       rorntComp[face].index = faceComp[face-fStart];
28084961fc4SMatthew G. Knepley     }
281e1d83109SMatthew G. Knepley     /* Communicate boundary edge orientations */
282e1d83109SMatthew G. Knepley     ierr = PetscSFBcastBegin(sf, MPIU_2INT, rorntComp, lorntComp);CHKERRQ(ierr);
283e1d83109SMatthew G. Knepley     ierr = PetscSFBcastEnd(sf, MPIU_2INT, rorntComp, lorntComp);CHKERRQ(ierr);
284e1d83109SMatthew G. Knepley   }
285e1d83109SMatthew G. Knepley   /* Get process adjacency */
286e1d83109SMatthew G. Knepley   ierr = PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors);CHKERRQ(ierr);
287e1d83109SMatthew G. Knepley   for (comp = 0; comp < numComponents; ++comp) {
288e1d83109SMatthew G. Knepley     PetscInt  l, n;
28984961fc4SMatthew G. Knepley 
290e1d83109SMatthew G. Knepley     numNeighbors[comp] = 0;
29139ea070eSMatthew G. Knepley     ierr = PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp]);CHKERRQ(ierr);
292e1d83109SMatthew G. Knepley     /* I know this is p^2 time in general, but for bounded degree its alright */
293e1d83109SMatthew G. Knepley     for (l = 0; l < numLeaves; ++l) {
294e1d83109SMatthew G. Knepley       const PetscInt face = lpoints[l];
295e1d83109SMatthew G. Knepley 
296e1d83109SMatthew G. Knepley       /* Find a representative face (edge) separating pairs of procs */
297e1d83109SMatthew G. Knepley       if ((face >= fStart) && (face < fEnd) && (faceComp[face-fStart] == comp)) {
298269a49d5SMatthew G. Knepley         const PetscInt rrank = rpoints[l].rank;
299e1d83109SMatthew G. Knepley         const PetscInt rcomp = lorntComp[face].index;
300e1d83109SMatthew G. Knepley 
301269a49d5SMatthew G. Knepley         for (n = 0; n < numNeighbors[comp]; ++n) if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break;
302e1d83109SMatthew G. Knepley         if (n >= numNeighbors[comp]) {
303e1d83109SMatthew G. Knepley           PetscInt supportSize;
304e1d83109SMatthew G. Knepley 
305e1d83109SMatthew G. Knepley           ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
306e1d83109SMatthew G. Knepley           if (supportSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %d", supportSize);
307a83982f3SMatthew G. Knepley           if (flg) {ierr = PetscPrintf(PETSC_COMM_SELF, "[%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);CHKERRQ(ierr);}
308e1d83109SMatthew G. Knepley           neighbors[comp][numNeighbors[comp]++] = l;
309e1d83109SMatthew G. Knepley         }
310e1d83109SMatthew G. Knepley       }
311e1d83109SMatthew G. Knepley     }
312e1d83109SMatthew G. Knepley     totNeighbors += numNeighbors[comp];
313e1d83109SMatthew G. Knepley   }
314e1d83109SMatthew G. Knepley   ierr = PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match);CHKERRQ(ierr);
315e1d83109SMatthew G. Knepley   for (comp = 0, off = 0; comp < numComponents; ++comp) {
316e1d83109SMatthew G. Knepley     PetscInt n;
317e1d83109SMatthew G. Knepley 
318e1d83109SMatthew G. Knepley     for (n = 0; n < numNeighbors[comp]; ++n, ++off) {
319e1d83109SMatthew G. Knepley       const PetscInt face = lpoints[neighbors[comp][n]];
320e1d83109SMatthew G. Knepley       const PetscInt o    = rorntComp[face].rank*lorntComp[face].rank;
321e1d83109SMatthew G. Knepley 
322e1d83109SMatthew G. Knepley       if      (o < 0) match[off] = PETSC_TRUE;
323e1d83109SMatthew G. Knepley       else if (o > 0) match[off] = PETSC_FALSE;
324e1d83109SMatthew 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);
325e1d83109SMatthew G. Knepley       nrankComp[off].rank  = rpoints[neighbors[comp][n]].rank;
326e1d83109SMatthew G. Knepley       nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index;
327e1d83109SMatthew G. Knepley     }
328e1d83109SMatthew G. Knepley     ierr = PetscFree(neighbors[comp]);CHKERRQ(ierr);
32984961fc4SMatthew G. Knepley   }
33084961fc4SMatthew G. Knepley   /* Collect the graph on 0 */
331e1d83109SMatthew G. Knepley   if (numLeaves >= 0) {
33284961fc4SMatthew G. Knepley     Mat          G;
33384961fc4SMatthew G. Knepley     PetscBT      seenProcs, flippedProcs;
33484961fc4SMatthew G. Knepley     PetscInt    *procFIFO, pTop, pBottom;
3357cadcfe8SMatthew G. Knepley     PetscInt    *N   = NULL, *Noff;
336e1d83109SMatthew G. Knepley     PetscSFNode *adj = NULL;
33784961fc4SMatthew G. Knepley     PetscBool   *val = NULL;
3387cadcfe8SMatthew G. Knepley     PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc, p, o;
339*9852e123SBarry Smith     PetscMPIInt  size = 0;
34084961fc4SMatthew G. Knepley 
341e1d83109SMatthew G. Knepley     ierr = PetscCalloc1(numComponents, &flipped);CHKERRQ(ierr);
342*9852e123SBarry Smith     if (!rank) {ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);}
343*9852e123SBarry Smith     ierr = PetscCalloc4(size, &recvcounts, size+1, &displs, size, &Nc, size+1, &Noff);CHKERRQ(ierr);
344e1d83109SMatthew G. Knepley     ierr = MPI_Gather(&numComponents, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm);CHKERRQ(ierr);
345*9852e123SBarry Smith     for (p = 0; p < size; ++p) {
346e1d83109SMatthew G. Knepley       displs[p+1] = displs[p] + Nc[p];
347e1d83109SMatthew G. Knepley     }
348*9852e123SBarry Smith     if (!rank) {ierr = PetscMalloc1(displs[size],&N);CHKERRQ(ierr);}
3497cadcfe8SMatthew G. Knepley     ierr = MPI_Gatherv(numNeighbors, numComponents, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm);CHKERRQ(ierr);
350*9852e123SBarry Smith     for (p = 0, o = 0; p < size; ++p) {
351e1d83109SMatthew G. Knepley       recvcounts[p] = 0;
352e1d83109SMatthew G. Knepley       for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o];
35384961fc4SMatthew G. Knepley       displs[p+1] = displs[p] + recvcounts[p];
35484961fc4SMatthew G. Knepley     }
355*9852e123SBarry Smith     if (!rank) {ierr = PetscMalloc2(displs[size], &adj, displs[size], &val);CHKERRQ(ierr);}
356e1d83109SMatthew G. Knepley     ierr = MPI_Gatherv(nrankComp, totNeighbors, MPIU_2INT, adj, recvcounts, displs, MPIU_2INT, 0, comm);CHKERRQ(ierr);
357e1d83109SMatthew G. Knepley     ierr = MPI_Gatherv(match, totNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
358e1d83109SMatthew G. Knepley     ierr = PetscFree2(numNeighbors, neighbors);CHKERRQ(ierr);
359e1d83109SMatthew G. Knepley     if (!rank) {
360*9852e123SBarry Smith       for (p = 1; p <= size; ++p) {Noff[p] = Noff[p-1] + Nc[p-1];}
361a83982f3SMatthew G. Knepley       if (flg) {
362e1d83109SMatthew G. Knepley         PetscInt n;
363e1d83109SMatthew G. Knepley 
364*9852e123SBarry Smith         for (p = 0, off = 0; p < size; ++p) {
365e1d83109SMatthew G. Knepley           for (c = 0; c < Nc[p]; ++c) {
366302440fdSBarry Smith             ierr = PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %d:\n", p, c);CHKERRQ(ierr);
367e1d83109SMatthew G. Knepley             for (n = 0; n < N[Noff[p]+c]; ++n, ++off) {
368302440fdSBarry Smith               ierr = PetscPrintf(PETSC_COMM_SELF, "  edge (%d, %d) (%d):\n", adj[off].rank, adj[off].index, val[off]);CHKERRQ(ierr);
369e1d83109SMatthew G. Knepley             }
37084961fc4SMatthew G. Knepley           }
37184961fc4SMatthew G. Knepley         }
37284961fc4SMatthew G. Knepley       }
37384961fc4SMatthew G. Knepley       /* Symmetrize the graph */
37484961fc4SMatthew G. Knepley       ierr = MatCreate(PETSC_COMM_SELF, &G);CHKERRQ(ierr);
375*9852e123SBarry Smith       ierr = MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size]);CHKERRQ(ierr);
37684961fc4SMatthew G. Knepley       ierr = MatSetUp(G);CHKERRQ(ierr);
377*9852e123SBarry Smith       for (p = 0, off = 0; p < size; ++p) {
378e1d83109SMatthew G. Knepley         for (c = 0; c < Nc[p]; ++c) {
379e1d83109SMatthew G. Knepley           const PetscInt r = Noff[p]+c;
380e1d83109SMatthew G. Knepley           PetscInt       n;
381e1d83109SMatthew G. Knepley 
382e1d83109SMatthew G. Knepley           for (n = 0; n < N[r]; ++n, ++off) {
383e1d83109SMatthew G. Knepley             const PetscInt    q = Noff[adj[off].rank] + adj[off].index;
384e1d83109SMatthew G. Knepley             const PetscScalar o = val[off] ? 1.0 : 0.0;
38584961fc4SMatthew G. Knepley 
38684961fc4SMatthew G. Knepley             ierr = MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES);CHKERRQ(ierr);
38784961fc4SMatthew G. Knepley             ierr = MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES);CHKERRQ(ierr);
38884961fc4SMatthew G. Knepley           }
38984961fc4SMatthew G. Knepley         }
390e1d83109SMatthew G. Knepley       }
39184961fc4SMatthew G. Knepley       ierr = MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
39284961fc4SMatthew G. Knepley       ierr = MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
39384961fc4SMatthew G. Knepley 
394*9852e123SBarry Smith       ierr = PetscBTCreate(Noff[size], &seenProcs);CHKERRQ(ierr);
395*9852e123SBarry Smith       ierr = PetscBTMemzero(Noff[size], seenProcs);CHKERRQ(ierr);
396*9852e123SBarry Smith       ierr = PetscBTCreate(Noff[size], &flippedProcs);CHKERRQ(ierr);
397*9852e123SBarry Smith       ierr = PetscBTMemzero(Noff[size], flippedProcs);CHKERRQ(ierr);
398*9852e123SBarry Smith       ierr = PetscMalloc1(Noff[size], &procFIFO);CHKERRQ(ierr);
39984961fc4SMatthew G. Knepley       pTop = pBottom = 0;
400*9852e123SBarry Smith       for (p = 0; p < Noff[size]; ++p) {
40184961fc4SMatthew G. Knepley         if (PetscBTLookup(seenProcs, p)) continue;
40284961fc4SMatthew G. Knepley         /* Initialize FIFO with next proc */
40384961fc4SMatthew G. Knepley         procFIFO[pBottom++] = p;
40484961fc4SMatthew G. Knepley         ierr = PetscBTSet(seenProcs, p);CHKERRQ(ierr);
40584961fc4SMatthew G. Knepley         /* Consider each proc in FIFO */
40684961fc4SMatthew G. Knepley         while (pTop < pBottom) {
40784961fc4SMatthew G. Knepley           const PetscScalar *ornt;
40884961fc4SMatthew G. Knepley           const PetscInt    *neighbors;
409e1d83109SMatthew G. Knepley           PetscInt           proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n;
41084961fc4SMatthew G. Knepley 
41184961fc4SMatthew G. Knepley           proc     = procFIFO[pTop++];
41284961fc4SMatthew G. Knepley           flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0;
41384961fc4SMatthew G. Knepley           ierr = MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt);CHKERRQ(ierr);
41484961fc4SMatthew G. Knepley           /* Loop over neighboring procs */
41584961fc4SMatthew G. Knepley           for (n = 0; n < numNeighbors; ++n) {
41684961fc4SMatthew G. Knepley             nproc    = neighbors[n];
41784961fc4SMatthew G. Knepley             mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1;
41884961fc4SMatthew G. Knepley             seen     = PetscBTLookup(seenProcs, nproc);
41984961fc4SMatthew G. Knepley             flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0;
42084961fc4SMatthew G. Knepley 
42184961fc4SMatthew G. Knepley             if (mismatch ^ (flippedA ^ flippedB)) {
42284961fc4SMatthew 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);
42384961fc4SMatthew G. Knepley               if (!flippedB) {
42484961fc4SMatthew G. Knepley                 ierr = PetscBTSet(flippedProcs, nproc);CHKERRQ(ierr);
42584961fc4SMatthew G. Knepley               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
42684961fc4SMatthew 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");
42784961fc4SMatthew G. Knepley             if (!seen) {
42884961fc4SMatthew G. Knepley               procFIFO[pBottom++] = nproc;
42984961fc4SMatthew G. Knepley               ierr = PetscBTSet(seenProcs, nproc);CHKERRQ(ierr);
43084961fc4SMatthew G. Knepley             }
43184961fc4SMatthew G. Knepley           }
43284961fc4SMatthew G. Knepley         }
43384961fc4SMatthew G. Knepley       }
43484961fc4SMatthew G. Knepley       ierr = PetscFree(procFIFO);CHKERRQ(ierr);
43584961fc4SMatthew G. Knepley       ierr = MatDestroy(&G);CHKERRQ(ierr);
43684961fc4SMatthew G. Knepley       ierr = PetscFree2(adj, val);CHKERRQ(ierr);
437e1d83109SMatthew G. Knepley       ierr = PetscBTDestroy(&seenProcs);CHKERRQ(ierr);
43884961fc4SMatthew G. Knepley     }
439e1d83109SMatthew G. Knepley     /* Scatter flip flags */
440e1d83109SMatthew G. Knepley     {
441e1d83109SMatthew G. Knepley       PetscBool *flips = NULL;
442e1d83109SMatthew G. Knepley 
443e1d83109SMatthew G. Knepley       if (!rank) {
444*9852e123SBarry Smith         ierr = PetscMalloc1(Noff[size], &flips);CHKERRQ(ierr);
445*9852e123SBarry Smith         for (p = 0; p < Noff[size]; ++p) {
446e1d83109SMatthew G. Knepley           flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE;
447302440fdSBarry Smith           if (flg && flips[p]) {ierr = PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p);CHKERRQ(ierr);}
448e1d83109SMatthew G. Knepley         }
449*9852e123SBarry Smith         for (p = 0; p < size; ++p) {
450e1d83109SMatthew G. Knepley           displs[p+1] = displs[p] + Nc[p];
451e1d83109SMatthew G. Knepley         }
452e1d83109SMatthew G. Knepley       }
453e1d83109SMatthew G. Knepley       ierr = MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, numComponents, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
45484961fc4SMatthew G. Knepley       ierr = PetscFree(flips);CHKERRQ(ierr);
45584961fc4SMatthew G. Knepley     }
456e1d83109SMatthew G. Knepley     if (!rank) {ierr = PetscBTDestroy(&flippedProcs);CHKERRQ(ierr);}
457e1d83109SMatthew G. Knepley     ierr = PetscFree(N);CHKERRQ(ierr);
458e1d83109SMatthew G. Knepley     ierr = PetscFree4(recvcounts, displs, Nc, Noff);CHKERRQ(ierr);
459e1d83109SMatthew G. Knepley     ierr = PetscFree2(nrankComp, match);CHKERRQ(ierr);
460e1d83109SMatthew G. Knepley 
461e1d83109SMatthew G. Knepley     /* Decide whether to flip cells in each component */
462e1d83109SMatthew G. Knepley     for (c = 0; c < cEnd-cStart; ++c) {if (flipped[cellComp[c]]) {ierr = PetscBTNegate(flippedCells, c);CHKERRQ(ierr);}}
463e1d83109SMatthew G. Knepley     ierr = PetscFree(flipped);CHKERRQ(ierr);
46484961fc4SMatthew G. Knepley   }
46584961fc4SMatthew G. Knepley   if (flg) {
46684961fc4SMatthew G. Knepley     PetscViewer v;
46784961fc4SMatthew G. Knepley 
4687e09831bSMatthew G. Knepley     ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr);
4691575c14dSBarry Smith     ierr = PetscViewerASCIIPushSynchronized(v);CHKERRQ(ierr);
47084961fc4SMatthew G. Knepley     ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank);CHKERRQ(ierr);
47184961fc4SMatthew G. Knepley     ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr);
4724d4c343aSBarry Smith     ierr = PetscViewerFlush(v);CHKERRQ(ierr);
4731575c14dSBarry Smith     ierr = PetscViewerASCIIPopSynchronized(v);CHKERRQ(ierr);
47484961fc4SMatthew G. Knepley   }
47584961fc4SMatthew G. Knepley   /* Reverse flipped cells in the mesh */
47684961fc4SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
47784961fc4SMatthew G. Knepley     if (PetscBTLookup(flippedCells, c-cStart)) {ierr = DMPlexReverseCell(dm, c);CHKERRQ(ierr);}
47884961fc4SMatthew G. Knepley   }
47984961fc4SMatthew G. Knepley   ierr = PetscBTDestroy(&seenCells);CHKERRQ(ierr);
48084961fc4SMatthew G. Knepley   ierr = PetscBTDestroy(&flippedCells);CHKERRQ(ierr);
48184961fc4SMatthew G. Knepley   ierr = PetscBTDestroy(&seenFaces);CHKERRQ(ierr);
482e1d83109SMatthew G. Knepley   ierr = PetscFree2(numNeighbors, neighbors);CHKERRQ(ierr);
483e1d83109SMatthew G. Knepley   ierr = PetscFree2(rorntComp, lorntComp);CHKERRQ(ierr);
484e1d83109SMatthew G. Knepley   ierr = PetscFree3(faceFIFO, cellComp, faceComp);CHKERRQ(ierr);
48584961fc4SMatthew G. Knepley   PetscFunctionReturn(0);
48684961fc4SMatthew G. Knepley }
487