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 PetscErrorCode ierr; 2527d0e99aSVaclav Hapla 2627d0e99aSVaclav Hapla PetscFunctionBegin; 2727d0e99aSVaclav Hapla PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 28b5a892a1SMatthew G. Knepley ierr = DMPlexGetCellType(dm, p, &ct);CHKERRQ(ierr); 29b5a892a1SMatthew G. Knepley arr = DMPolytopeTypeGetArrangment(ct, o); 3027d0e99aSVaclav Hapla ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 3127d0e99aSVaclav Hapla ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 32b5a892a1SMatthew G. Knepley ierr = DMPlexGetConeOrientation(dm, p, &ornt);CHKERRQ(ierr); 33b5a892a1SMatthew G. Knepley ierr = DMGetWorkArray(dm, coneSize, MPIU_INT, &newcone);CHKERRQ(ierr); 34b5a892a1SMatthew G. Knepley ierr = DMGetWorkArray(dm, coneSize, MPIU_INT, &newornt);CHKERRQ(ierr); 35b5a892a1SMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 36b5a892a1SMatthew G. Knepley DMPolytopeType ft; 37b5a892a1SMatthew G. Knepley PetscInt nO; 3827d0e99aSVaclav Hapla 39b5a892a1SMatthew G. Knepley ierr = DMPlexGetCellType(dm, cone[c], &ft);CHKERRQ(ierr); 40b5a892a1SMatthew G. Knepley nO = DMPolytopeTypeGetNumArrangments(ft)/2; 41b5a892a1SMatthew G. Knepley newcone[c] = cone[arr[c*2+0]]; 42b5a892a1SMatthew G. Knepley newornt[c] = DMPolytopeTypeComposeOrientation(ft, arr[c*2+1], ornt[arr[c*2+0]]); 43*2c71b3e2SJacob 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]); 4427d0e99aSVaclav Hapla } 4527d0e99aSVaclav Hapla ierr = DMPlexSetCone(dm, p, newcone);CHKERRQ(ierr); 46b5a892a1SMatthew G. Knepley ierr = DMPlexSetConeOrientation(dm, p, newornt);CHKERRQ(ierr); 47b5a892a1SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, coneSize, MPIU_INT, &newcone);CHKERRQ(ierr); 48b5a892a1SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, coneSize, MPIU_INT, &newornt);CHKERRQ(ierr); 49b5a892a1SMatthew G. Knepley /* Update orientation of this point in the support points */ 50b5a892a1SMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, p, &supportSize);CHKERRQ(ierr); 51b5a892a1SMatthew G. Knepley ierr = DMPlexGetSupport(dm, p, &support);CHKERRQ(ierr); 52b5a892a1SMatthew G. Knepley for (s = 0; s < supportSize; ++s) { 53b5a892a1SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, support[s], &coneSize);CHKERRQ(ierr); 54b5a892a1SMatthew G. Knepley ierr = DMPlexGetCone(dm, support[s], &cone);CHKERRQ(ierr); 55b5a892a1SMatthew G. Knepley ierr = DMPlexGetConeOrientation(dm, support[s], &ornt);CHKERRQ(ierr); 56b5a892a1SMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 57b5a892a1SMatthew G. Knepley PetscInt po; 5827d0e99aSVaclav Hapla 59b5a892a1SMatthew G. Knepley if (cone[c] != p) continue; 60b5a892a1SMatthew G. Knepley /* ornt[c] * 0 = target = po * o so that po = ornt[c] * o^{-1} */ 61b5a892a1SMatthew G. Knepley po = DMPolytopeTypeComposeOrientationInv(ct, ornt[c], o); 62b5a892a1SMatthew G. Knepley ierr = DMPlexInsertConeOrientation(dm, support[s], c, po);CHKERRQ(ierr); 6384961fc4SMatthew G. Knepley } 6484961fc4SMatthew G. Knepley } 6584961fc4SMatthew G. Knepley PetscFunctionReturn(0); 6684961fc4SMatthew G. Knepley } 6784961fc4SMatthew G. Knepley 687b310ce2SMatthew G. Knepley /* 697b310ce2SMatthew G. Knepley - Checks face match 707b310ce2SMatthew G. Knepley - Flips non-matching 717b310ce2SMatthew G. Knepley - Inserts faces of support cells in FIFO 727b310ce2SMatthew G. Knepley */ 737b2de0fdSMatthew 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) 747b310ce2SMatthew G. Knepley { 757b310ce2SMatthew G. Knepley const PetscInt *support, *coneA, *coneB, *coneOA, *coneOB; 767b310ce2SMatthew G. Knepley PetscInt supportSize, coneSizeA, coneSizeB, posA = -1, posB = -1; 777b2de0fdSMatthew G. Knepley PetscInt face, dim, seenA, flippedA, seenB, flippedB, mismatch, c; 787b310ce2SMatthew G. Knepley PetscErrorCode ierr; 797b310ce2SMatthew G. Knepley 807b310ce2SMatthew G. Knepley PetscFunctionBegin; 817b310ce2SMatthew G. Knepley face = faceFIFO[(*fTop)++]; 827b310ce2SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 837b310ce2SMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); 847b310ce2SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr); 857b310ce2SMatthew G. Knepley if (supportSize < 2) PetscFunctionReturn(0); 86*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(supportSize != 2,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %d", supportSize); 877b310ce2SMatthew G. Knepley seenA = PetscBTLookup(seenCells, support[0]-cStart); 887b310ce2SMatthew G. Knepley flippedA = PetscBTLookup(flippedCells, support[0]-cStart) ? 1 : 0; 897b310ce2SMatthew G. Knepley seenB = PetscBTLookup(seenCells, support[1]-cStart); 907b310ce2SMatthew G. Knepley flippedB = PetscBTLookup(flippedCells, support[1]-cStart) ? 1 : 0; 917b310ce2SMatthew G. Knepley 927b310ce2SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, support[0], &coneSizeA);CHKERRQ(ierr); 937b310ce2SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, support[1], &coneSizeB);CHKERRQ(ierr); 947b310ce2SMatthew G. Knepley ierr = DMPlexGetCone(dm, support[0], &coneA);CHKERRQ(ierr); 957b310ce2SMatthew G. Knepley ierr = DMPlexGetCone(dm, support[1], &coneB);CHKERRQ(ierr); 967b310ce2SMatthew G. Knepley ierr = DMPlexGetConeOrientation(dm, support[0], &coneOA);CHKERRQ(ierr); 977b310ce2SMatthew G. Knepley ierr = DMPlexGetConeOrientation(dm, support[1], &coneOB);CHKERRQ(ierr); 987b310ce2SMatthew G. Knepley for (c = 0; c < coneSizeA; ++c) { 997b310ce2SMatthew G. Knepley if (!PetscBTLookup(seenFaces, coneA[c]-fStart)) { 1007b310ce2SMatthew G. Knepley faceFIFO[(*fBottom)++] = coneA[c]; 1017b310ce2SMatthew G. Knepley ierr = PetscBTSet(seenFaces, coneA[c]-fStart);CHKERRQ(ierr); 1027b310ce2SMatthew G. Knepley } 1037b310ce2SMatthew G. Knepley if (coneA[c] == face) posA = c; 104*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(*fBottom > fEnd-fStart,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %d was pushed exceeding capacity %d > %d", coneA[c], *fBottom, fEnd-fStart); 1057b310ce2SMatthew G. Knepley } 106*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(posA < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d could not be located in cell %d", face, support[0]); 1077b310ce2SMatthew G. Knepley for (c = 0; c < coneSizeB; ++c) { 1087b310ce2SMatthew G. Knepley if (!PetscBTLookup(seenFaces, coneB[c]-fStart)) { 1097b310ce2SMatthew G. Knepley faceFIFO[(*fBottom)++] = coneB[c]; 1107b310ce2SMatthew G. Knepley ierr = PetscBTSet(seenFaces, coneB[c]-fStart);CHKERRQ(ierr); 1117b310ce2SMatthew G. Knepley } 1127b310ce2SMatthew G. Knepley if (coneB[c] == face) posB = c; 113*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(*fBottom > fEnd-fStart,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %d was pushed exceeding capacity %d > %d", coneA[c], *fBottom, fEnd-fStart); 1147b310ce2SMatthew G. Knepley } 115*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(posB < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d could not be located in cell %d", face, support[1]); 1167b310ce2SMatthew G. Knepley 1177b310ce2SMatthew G. Knepley if (dim == 1) { 1187b310ce2SMatthew G. Knepley mismatch = posA == posB; 1197b310ce2SMatthew G. Knepley } else { 1207b310ce2SMatthew G. Knepley mismatch = coneOA[posA] == coneOB[posB]; 1217b310ce2SMatthew G. Knepley } 1227b310ce2SMatthew G. Knepley 1237b310ce2SMatthew G. Knepley if (mismatch ^ (flippedA ^ flippedB)) { 124*2c71b3e2SJacob 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]); 1257b310ce2SMatthew G. Knepley if (!seenA && !flippedA) { 1267b310ce2SMatthew G. Knepley ierr = PetscBTSet(flippedCells, support[0]-cStart);CHKERRQ(ierr); 1277b310ce2SMatthew G. Knepley } else if (!seenB && !flippedB) { 1287b310ce2SMatthew G. Knepley ierr = PetscBTSet(flippedCells, support[1]-cStart);CHKERRQ(ierr); 1297b310ce2SMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable"); 130*2c71b3e2SJacob Faibussowitsch } else PetscCheckFalse(mismatch && flippedA && flippedB,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable"); 1317b310ce2SMatthew G. Knepley ierr = PetscBTSet(seenCells, support[0]-cStart);CHKERRQ(ierr); 1327b310ce2SMatthew G. Knepley ierr = PetscBTSet(seenCells, support[1]-cStart);CHKERRQ(ierr); 1337b310ce2SMatthew G. Knepley PetscFunctionReturn(0); 1347b310ce2SMatthew G. Knepley } 1357b310ce2SMatthew G. Knepley 13684961fc4SMatthew G. Knepley /*@ 13784961fc4SMatthew G. Knepley DMPlexOrient - Give a consistent orientation to the input mesh 13884961fc4SMatthew G. Knepley 13984961fc4SMatthew G. Knepley Input Parameters: 14084961fc4SMatthew G. Knepley . dm - The DM 14184961fc4SMatthew G. Knepley 14284961fc4SMatthew G. Knepley Note: The orientation data for the DM are change in-place. 14384961fc4SMatthew G. Knepley $ This routine will fail for non-orientable surfaces, such as the Moebius strip. 14484961fc4SMatthew G. Knepley 14584961fc4SMatthew G. Knepley Level: advanced 14684961fc4SMatthew G. Knepley 14784961fc4SMatthew G. Knepley .seealso: DMCreate(), DMPLEX 14884961fc4SMatthew G. Knepley @*/ 14984961fc4SMatthew G. Knepley PetscErrorCode DMPlexOrient(DM dm) 15084961fc4SMatthew G. Knepley { 15184961fc4SMatthew G. Knepley MPI_Comm comm; 152e1d83109SMatthew G. Knepley PetscSF sf; 153e1d83109SMatthew G. Knepley const PetscInt *lpoints; 154e1d83109SMatthew G. Knepley const PetscSFNode *rpoints; 155e1d83109SMatthew G. Knepley PetscSFNode *rorntComp = NULL, *lorntComp = NULL; 156e1d83109SMatthew G. Knepley PetscInt *numNeighbors, **neighbors; 157e1d83109SMatthew G. Knepley PetscSFNode *nrankComp; 158e1d83109SMatthew G. Knepley PetscBool *match, *flipped; 15984961fc4SMatthew G. Knepley PetscBT seenCells, flippedCells, seenFaces; 160e1d83109SMatthew G. Knepley PetscInt *faceFIFO, fTop, fBottom, *cellComp, *faceComp; 1617cadcfe8SMatthew G. Knepley PetscInt numLeaves, numRoots, dim, h, cStart, cEnd, c, cell, fStart, fEnd, face, off, totNeighbors = 0; 16239526728SToby Isaac PetscMPIInt rank, size, numComponents, comp = 0; 16339526728SToby Isaac PetscBool flg, flg2; 164a17aa47eSToby Isaac PetscViewer viewer = NULL, selfviewer = NULL; 16584961fc4SMatthew G. Knepley PetscErrorCode ierr; 16684961fc4SMatthew G. Knepley 16784961fc4SMatthew G. Knepley PetscFunctionBegin; 16884961fc4SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 169ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 170ffc4695bSBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr); 171c5929fdfSBarry Smith ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view", &flg);CHKERRQ(ierr); 17239526728SToby Isaac ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view_synchronized", &flg2);CHKERRQ(ierr); 173e1d83109SMatthew G. Knepley ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 174e1d83109SMatthew G. Knepley ierr = PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints);CHKERRQ(ierr); 17584961fc4SMatthew G. Knepley /* Truth Table 17684961fc4SMatthew G. Knepley mismatch flips do action mismatch flipA ^ flipB action 17784961fc4SMatthew G. Knepley F 0 flips no F F F 17884961fc4SMatthew G. Knepley F 1 flip yes F T T 17984961fc4SMatthew G. Knepley F 2 flips no T F T 18084961fc4SMatthew G. Knepley T 0 flips yes T T F 18184961fc4SMatthew G. Knepley T 1 flip no 18284961fc4SMatthew G. Knepley T 2 flips yes 18384961fc4SMatthew G. Knepley */ 18484961fc4SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 18584961fc4SMatthew G. Knepley ierr = DMPlexGetVTKCellHeight(dm, &h);CHKERRQ(ierr); 18684961fc4SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, h, &cStart, &cEnd);CHKERRQ(ierr); 18784961fc4SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, h+1, &fStart, &fEnd);CHKERRQ(ierr); 18884961fc4SMatthew G. Knepley ierr = PetscBTCreate(cEnd - cStart, &seenCells);CHKERRQ(ierr); 18984961fc4SMatthew G. Knepley ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr); 19084961fc4SMatthew G. Knepley ierr = PetscBTCreate(cEnd - cStart, &flippedCells);CHKERRQ(ierr); 19184961fc4SMatthew G. Knepley ierr = PetscBTMemzero(cEnd - cStart, flippedCells);CHKERRQ(ierr); 19284961fc4SMatthew G. Knepley ierr = PetscBTCreate(fEnd - fStart, &seenFaces);CHKERRQ(ierr); 19384961fc4SMatthew G. Knepley ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr); 194e1d83109SMatthew G. Knepley ierr = PetscCalloc3(fEnd - fStart, &faceFIFO, cEnd-cStart, &cellComp, fEnd-fStart, &faceComp);CHKERRQ(ierr); 195e1d83109SMatthew G. Knepley /* 196e1d83109SMatthew G. Knepley OLD STYLE 197e1d83109SMatthew G. Knepley - Add an integer array over cells and faces (component) for connected component number 198e1d83109SMatthew G. Knepley Foreach component 199e1d83109SMatthew G. Knepley - Mark the initial cell as seen 200e1d83109SMatthew G. Knepley - Process component as usual 201e1d83109SMatthew G. Knepley - Set component for all seenCells 202e1d83109SMatthew G. Knepley - Wipe seenCells and seenFaces (flippedCells can stay) 203e1d83109SMatthew G. Knepley - Generate parallel adjacency for component using SF and seenFaces 204e1d83109SMatthew G. Knepley - Collect numComponents adj data from each proc to 0 205e1d83109SMatthew G. Knepley - Build same serial graph 206e1d83109SMatthew G. Knepley - Use same solver 207e1d83109SMatthew G. Knepley - Use Scatterv to to send back flipped flags for each component 208e1d83109SMatthew G. Knepley - Negate flippedCells by component 209e1d83109SMatthew G. Knepley 210e1d83109SMatthew G. Knepley NEW STYLE 211e1d83109SMatthew G. Knepley - Create the adj on each process 212e1d83109SMatthew G. Knepley - Bootstrap to complete graph on proc 0 213e1d83109SMatthew G. Knepley */ 214e1d83109SMatthew G. Knepley /* Loop over components */ 215e1d83109SMatthew G. Knepley for (cell = cStart; cell < cEnd; ++cell) cellComp[cell-cStart] = -1; 216e1d83109SMatthew G. Knepley do { 217e1d83109SMatthew G. Knepley /* Look for first unmarked cell */ 218e1d83109SMatthew G. Knepley for (cell = cStart; cell < cEnd; ++cell) if (cellComp[cell-cStart] < 0) break; 219e1d83109SMatthew G. Knepley if (cell >= cEnd) break; 220e1d83109SMatthew G. Knepley /* Initialize FIFO with first cell in component */ 221e1d83109SMatthew G. Knepley { 22284961fc4SMatthew G. Knepley const PetscInt *cone; 22384961fc4SMatthew G. Knepley PetscInt coneSize; 22484961fc4SMatthew G. Knepley 225e1d83109SMatthew G. Knepley fTop = fBottom = 0; 226e1d83109SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 227e1d83109SMatthew G. Knepley ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr); 22884961fc4SMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 22984961fc4SMatthew G. Knepley faceFIFO[fBottom++] = cone[c]; 23084961fc4SMatthew G. Knepley ierr = PetscBTSet(seenFaces, cone[c]-fStart);CHKERRQ(ierr); 23184961fc4SMatthew G. Knepley } 23231c8331aSMatthew G. Knepley ierr = PetscBTSet(seenCells, cell-cStart);CHKERRQ(ierr); 23384961fc4SMatthew G. Knepley } 23484961fc4SMatthew G. Knepley /* Consider each face in FIFO */ 23584961fc4SMatthew G. Knepley while (fTop < fBottom) { 2367b2de0fdSMatthew G. Knepley ierr = DMPlexCheckFace_Internal(dm, faceFIFO, &fTop, &fBottom, cStart, fStart, fEnd, seenCells, flippedCells, seenFaces);CHKERRQ(ierr); 23784961fc4SMatthew G. Knepley } 238e1d83109SMatthew G. Knepley /* Set component for cells and faces */ 239e1d83109SMatthew G. Knepley for (cell = 0; cell < cEnd-cStart; ++cell) { 240e1d83109SMatthew G. Knepley if (PetscBTLookup(seenCells, cell)) cellComp[cell] = comp; 241e1d83109SMatthew G. Knepley } 242e1d83109SMatthew G. Knepley for (face = 0; face < fEnd-fStart; ++face) { 243e1d83109SMatthew G. Knepley if (PetscBTLookup(seenFaces, face)) faceComp[face] = comp; 244e1d83109SMatthew G. Knepley } 245e1d83109SMatthew G. Knepley /* Wipe seenCells and seenFaces for next component */ 246e1d83109SMatthew G. Knepley ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr); 247e1d83109SMatthew G. Knepley ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr); 248e1d83109SMatthew G. Knepley ++comp; 249e1d83109SMatthew G. Knepley } while (1); 250e1d83109SMatthew G. Knepley numComponents = comp; 25184961fc4SMatthew G. Knepley if (flg) { 25284961fc4SMatthew G. Knepley PetscViewer v; 25384961fc4SMatthew G. Knepley 2547e09831bSMatthew G. Knepley ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr); 2551575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(v);CHKERRQ(ierr); 25684961fc4SMatthew G. Knepley ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank);CHKERRQ(ierr); 25784961fc4SMatthew G. Knepley ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr); 2584d4c343aSBarry Smith ierr = PetscViewerFlush(v);CHKERRQ(ierr); 2591575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(v);CHKERRQ(ierr); 26084961fc4SMatthew G. Knepley } 26184961fc4SMatthew G. Knepley /* Now all subdomains are oriented, but we need a consistent parallel orientation */ 26284961fc4SMatthew G. Knepley if (numLeaves >= 0) { 263e1d83109SMatthew G. Knepley /* Store orientations of boundary faces*/ 264e1d83109SMatthew G. Knepley ierr = PetscCalloc2(numRoots,&rorntComp,numRoots,&lorntComp);CHKERRQ(ierr); 26584961fc4SMatthew G. Knepley for (face = fStart; face < fEnd; ++face) { 266e1d83109SMatthew G. Knepley const PetscInt *cone, *support, *ornt; 267e1d83109SMatthew G. Knepley PetscInt coneSize, supportSize; 268e1d83109SMatthew G. Knepley 26984961fc4SMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); 27084961fc4SMatthew G. Knepley if (supportSize != 1) continue; 27184961fc4SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr); 27284961fc4SMatthew G. Knepley 27384961fc4SMatthew G. Knepley ierr = DMPlexGetCone(dm, support[0], &cone);CHKERRQ(ierr); 27484961fc4SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, support[0], &coneSize);CHKERRQ(ierr); 27584961fc4SMatthew G. Knepley ierr = DMPlexGetConeOrientation(dm, support[0], &ornt);CHKERRQ(ierr); 27684961fc4SMatthew G. Knepley for (c = 0; c < coneSize; ++c) if (cone[c] == face) break; 27784961fc4SMatthew G. Knepley if (dim == 1) { 27884961fc4SMatthew G. Knepley /* Use cone position instead, shifted to -1 or 1 */ 2796e1eb25cSMatthew G. Knepley if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = 1-c*2; 2806e1eb25cSMatthew G. Knepley else rorntComp[face].rank = c*2-1; 28184961fc4SMatthew G. Knepley } else { 282e1d83109SMatthew G. Knepley if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 : 1; 283e1d83109SMatthew G. Knepley else rorntComp[face].rank = ornt[c] < 0 ? 1 : -1; 28484961fc4SMatthew G. Knepley } 285e1d83109SMatthew G. Knepley rorntComp[face].index = faceComp[face-fStart]; 28684961fc4SMatthew G. Knepley } 287e1d83109SMatthew G. Knepley /* Communicate boundary edge orientations */ 288ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(sf, MPIU_2INT, rorntComp, lorntComp,MPI_REPLACE);CHKERRQ(ierr); 289ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(sf, MPIU_2INT, rorntComp, lorntComp,MPI_REPLACE);CHKERRQ(ierr); 290e1d83109SMatthew G. Knepley } 291e1d83109SMatthew G. Knepley /* Get process adjacency */ 292e1d83109SMatthew G. Knepley ierr = PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors);CHKERRQ(ierr); 293a17aa47eSToby Isaac viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm)); 294a17aa47eSToby Isaac if (flg2) {ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);} 295a17aa47eSToby Isaac ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&selfviewer);CHKERRQ(ierr); 296e1d83109SMatthew G. Knepley for (comp = 0; comp < numComponents; ++comp) { 297e1d83109SMatthew G. Knepley PetscInt l, n; 29884961fc4SMatthew G. Knepley 299e1d83109SMatthew G. Knepley numNeighbors[comp] = 0; 30039ea070eSMatthew G. Knepley ierr = PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp]);CHKERRQ(ierr); 301e1d83109SMatthew G. Knepley /* I know this is p^2 time in general, but for bounded degree its alright */ 302e1d83109SMatthew G. Knepley for (l = 0; l < numLeaves; ++l) { 303e1d83109SMatthew G. Knepley const PetscInt face = lpoints[l]; 304e1d83109SMatthew G. Knepley 305e1d83109SMatthew G. Knepley /* Find a representative face (edge) separating pairs of procs */ 306e1d83109SMatthew G. Knepley if ((face >= fStart) && (face < fEnd) && (faceComp[face-fStart] == comp)) { 307269a49d5SMatthew G. Knepley const PetscInt rrank = rpoints[l].rank; 308e1d83109SMatthew G. Knepley const PetscInt rcomp = lorntComp[face].index; 309e1d83109SMatthew G. Knepley 310269a49d5SMatthew G. Knepley for (n = 0; n < numNeighbors[comp]; ++n) if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break; 311e1d83109SMatthew G. Knepley if (n >= numNeighbors[comp]) { 312e1d83109SMatthew G. Knepley PetscInt supportSize; 313e1d83109SMatthew G. Knepley 314e1d83109SMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); 315*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(supportSize != 1,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %d", supportSize); 31635041eceSToby Isaac if (flg) {ierr = 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);CHKERRQ(ierr);} 317e1d83109SMatthew G. Knepley neighbors[comp][numNeighbors[comp]++] = l; 318e1d83109SMatthew G. Knepley } 319e1d83109SMatthew G. Knepley } 320e1d83109SMatthew G. Knepley } 321e1d83109SMatthew G. Knepley totNeighbors += numNeighbors[comp]; 322e1d83109SMatthew G. Knepley } 323a17aa47eSToby Isaac ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&selfviewer);CHKERRQ(ierr); 324a17aa47eSToby Isaac ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 325a17aa47eSToby Isaac if (flg2) {ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);} 326e1d83109SMatthew G. Knepley ierr = PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match);CHKERRQ(ierr); 327e1d83109SMatthew G. Knepley for (comp = 0, off = 0; comp < numComponents; ++comp) { 328e1d83109SMatthew G. Knepley PetscInt n; 329e1d83109SMatthew G. Knepley 330e1d83109SMatthew G. Knepley for (n = 0; n < numNeighbors[comp]; ++n, ++off) { 331e1d83109SMatthew G. Knepley const PetscInt face = lpoints[neighbors[comp][n]]; 332e1d83109SMatthew G. Knepley const PetscInt o = rorntComp[face].rank*lorntComp[face].rank; 333e1d83109SMatthew G. Knepley 334e1d83109SMatthew G. Knepley if (o < 0) match[off] = PETSC_TRUE; 335e1d83109SMatthew G. Knepley else if (o > 0) match[off] = PETSC_FALSE; 33698921bdaSJacob 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); 337e1d83109SMatthew G. Knepley nrankComp[off].rank = rpoints[neighbors[comp][n]].rank; 338e1d83109SMatthew G. Knepley nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index; 339e1d83109SMatthew G. Knepley } 340e1d83109SMatthew G. Knepley ierr = PetscFree(neighbors[comp]);CHKERRQ(ierr); 34184961fc4SMatthew G. Knepley } 34284961fc4SMatthew G. Knepley /* Collect the graph on 0 */ 343e1d83109SMatthew G. Knepley if (numLeaves >= 0) { 34484961fc4SMatthew G. Knepley Mat G; 34584961fc4SMatthew G. Knepley PetscBT seenProcs, flippedProcs; 34684961fc4SMatthew G. Knepley PetscInt *procFIFO, pTop, pBottom; 3477cadcfe8SMatthew G. Knepley PetscInt *N = NULL, *Noff; 348e1d83109SMatthew G. Knepley PetscSFNode *adj = NULL; 34984961fc4SMatthew G. Knepley PetscBool *val = NULL; 3507cadcfe8SMatthew G. Knepley PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc, p, o; 3519852e123SBarry Smith PetscMPIInt size = 0; 35284961fc4SMatthew G. Knepley 353e1d83109SMatthew G. Knepley ierr = PetscCalloc1(numComponents, &flipped);CHKERRQ(ierr); 354dd400576SPatrick Sanan if (rank == 0) {ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);} 3559852e123SBarry Smith ierr = PetscCalloc4(size, &recvcounts, size+1, &displs, size, &Nc, size+1, &Noff);CHKERRQ(ierr); 356ffc4695bSBarry Smith ierr = MPI_Gather(&numComponents, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm);CHKERRMPI(ierr); 3579852e123SBarry Smith for (p = 0; p < size; ++p) { 358e1d83109SMatthew G. Knepley displs[p+1] = displs[p] + Nc[p]; 359e1d83109SMatthew G. Knepley } 360dd400576SPatrick Sanan if (rank == 0) {ierr = PetscMalloc1(displs[size],&N);CHKERRQ(ierr);} 361ffc4695bSBarry Smith ierr = MPI_Gatherv(numNeighbors, numComponents, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm);CHKERRMPI(ierr); 3629852e123SBarry Smith for (p = 0, o = 0; p < size; ++p) { 363e1d83109SMatthew G. Knepley recvcounts[p] = 0; 364e1d83109SMatthew G. Knepley for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o]; 36584961fc4SMatthew G. Knepley displs[p+1] = displs[p] + recvcounts[p]; 36684961fc4SMatthew G. Knepley } 367dd400576SPatrick Sanan if (rank == 0) {ierr = PetscMalloc2(displs[size], &adj, displs[size], &val);CHKERRQ(ierr);} 368ffc4695bSBarry Smith ierr = MPI_Gatherv(nrankComp, totNeighbors, MPIU_2INT, adj, recvcounts, displs, MPIU_2INT, 0, comm);CHKERRMPI(ierr); 369ffc4695bSBarry Smith ierr = MPI_Gatherv(match, totNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm);CHKERRMPI(ierr); 370e1d83109SMatthew G. Knepley ierr = PetscFree2(numNeighbors, neighbors);CHKERRQ(ierr); 371dd400576SPatrick Sanan if (rank == 0) { 3729852e123SBarry Smith for (p = 1; p <= size; ++p) {Noff[p] = Noff[p-1] + Nc[p-1];} 373a83982f3SMatthew G. Knepley if (flg) { 374e1d83109SMatthew G. Knepley PetscInt n; 375e1d83109SMatthew G. Knepley 3769852e123SBarry Smith for (p = 0, off = 0; p < size; ++p) { 377e1d83109SMatthew G. Knepley for (c = 0; c < Nc[p]; ++c) { 378302440fdSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %d:\n", p, c);CHKERRQ(ierr); 379e1d83109SMatthew G. Knepley for (n = 0; n < N[Noff[p]+c]; ++n, ++off) { 380302440fdSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF, " edge (%d, %d) (%d):\n", adj[off].rank, adj[off].index, val[off]);CHKERRQ(ierr); 381e1d83109SMatthew G. Knepley } 38284961fc4SMatthew G. Knepley } 38384961fc4SMatthew G. Knepley } 38484961fc4SMatthew G. Knepley } 38584961fc4SMatthew G. Knepley /* Symmetrize the graph */ 38684961fc4SMatthew G. Knepley ierr = MatCreate(PETSC_COMM_SELF, &G);CHKERRQ(ierr); 3879852e123SBarry Smith ierr = MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size]);CHKERRQ(ierr); 38884961fc4SMatthew G. Knepley ierr = MatSetUp(G);CHKERRQ(ierr); 3899852e123SBarry Smith for (p = 0, off = 0; p < size; ++p) { 390e1d83109SMatthew G. Knepley for (c = 0; c < Nc[p]; ++c) { 391e1d83109SMatthew G. Knepley const PetscInt r = Noff[p]+c; 392e1d83109SMatthew G. Knepley PetscInt n; 393e1d83109SMatthew G. Knepley 394e1d83109SMatthew G. Knepley for (n = 0; n < N[r]; ++n, ++off) { 395e1d83109SMatthew G. Knepley const PetscInt q = Noff[adj[off].rank] + adj[off].index; 396e1d83109SMatthew G. Knepley const PetscScalar o = val[off] ? 1.0 : 0.0; 39784961fc4SMatthew G. Knepley 39884961fc4SMatthew G. Knepley ierr = MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES);CHKERRQ(ierr); 39984961fc4SMatthew G. Knepley ierr = MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES);CHKERRQ(ierr); 40084961fc4SMatthew G. Knepley } 40184961fc4SMatthew G. Knepley } 402e1d83109SMatthew G. Knepley } 40384961fc4SMatthew G. Knepley ierr = MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 40484961fc4SMatthew G. Knepley ierr = MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 40584961fc4SMatthew G. Knepley 4069852e123SBarry Smith ierr = PetscBTCreate(Noff[size], &seenProcs);CHKERRQ(ierr); 4079852e123SBarry Smith ierr = PetscBTMemzero(Noff[size], seenProcs);CHKERRQ(ierr); 4089852e123SBarry Smith ierr = PetscBTCreate(Noff[size], &flippedProcs);CHKERRQ(ierr); 4099852e123SBarry Smith ierr = PetscBTMemzero(Noff[size], flippedProcs);CHKERRQ(ierr); 4109852e123SBarry Smith ierr = PetscMalloc1(Noff[size], &procFIFO);CHKERRQ(ierr); 41184961fc4SMatthew G. Knepley pTop = pBottom = 0; 4129852e123SBarry Smith for (p = 0; p < Noff[size]; ++p) { 41384961fc4SMatthew G. Knepley if (PetscBTLookup(seenProcs, p)) continue; 41484961fc4SMatthew G. Knepley /* Initialize FIFO with next proc */ 41584961fc4SMatthew G. Knepley procFIFO[pBottom++] = p; 41684961fc4SMatthew G. Knepley ierr = PetscBTSet(seenProcs, p);CHKERRQ(ierr); 41784961fc4SMatthew G. Knepley /* Consider each proc in FIFO */ 41884961fc4SMatthew G. Knepley while (pTop < pBottom) { 41984961fc4SMatthew G. Knepley const PetscScalar *ornt; 42084961fc4SMatthew G. Knepley const PetscInt *neighbors; 421e1d83109SMatthew G. Knepley PetscInt proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n; 42284961fc4SMatthew G. Knepley 42384961fc4SMatthew G. Knepley proc = procFIFO[pTop++]; 42484961fc4SMatthew G. Knepley flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0; 42584961fc4SMatthew G. Knepley ierr = MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt);CHKERRQ(ierr); 42684961fc4SMatthew G. Knepley /* Loop over neighboring procs */ 42784961fc4SMatthew G. Knepley for (n = 0; n < numNeighbors; ++n) { 42884961fc4SMatthew G. Knepley nproc = neighbors[n]; 42984961fc4SMatthew G. Knepley mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1; 43084961fc4SMatthew G. Knepley seen = PetscBTLookup(seenProcs, nproc); 43184961fc4SMatthew G. Knepley flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0; 43284961fc4SMatthew G. Knepley 43384961fc4SMatthew G. Knepley if (mismatch ^ (flippedA ^ flippedB)) { 434*2c71b3e2SJacob 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); 43584961fc4SMatthew G. Knepley if (!flippedB) { 43684961fc4SMatthew G. Knepley ierr = PetscBTSet(flippedProcs, nproc);CHKERRQ(ierr); 43784961fc4SMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable"); 438*2c71b3e2SJacob Faibussowitsch } else PetscCheckFalse(mismatch && flippedA && flippedB,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable"); 43984961fc4SMatthew G. Knepley if (!seen) { 44084961fc4SMatthew G. Knepley procFIFO[pBottom++] = nproc; 44184961fc4SMatthew G. Knepley ierr = PetscBTSet(seenProcs, nproc);CHKERRQ(ierr); 44284961fc4SMatthew G. Knepley } 44384961fc4SMatthew G. Knepley } 44484961fc4SMatthew G. Knepley } 44584961fc4SMatthew G. Knepley } 44684961fc4SMatthew G. Knepley ierr = PetscFree(procFIFO);CHKERRQ(ierr); 44784961fc4SMatthew G. Knepley ierr = MatDestroy(&G);CHKERRQ(ierr); 44884961fc4SMatthew G. Knepley ierr = PetscFree2(adj, val);CHKERRQ(ierr); 449e1d83109SMatthew G. Knepley ierr = PetscBTDestroy(&seenProcs);CHKERRQ(ierr); 45084961fc4SMatthew G. Knepley } 451e1d83109SMatthew G. Knepley /* Scatter flip flags */ 452e1d83109SMatthew G. Knepley { 453e1d83109SMatthew G. Knepley PetscBool *flips = NULL; 454e1d83109SMatthew G. Knepley 455dd400576SPatrick Sanan if (rank == 0) { 4569852e123SBarry Smith ierr = PetscMalloc1(Noff[size], &flips);CHKERRQ(ierr); 4579852e123SBarry Smith for (p = 0; p < Noff[size]; ++p) { 458e1d83109SMatthew G. Knepley flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE; 459302440fdSBarry Smith if (flg && flips[p]) {ierr = PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p);CHKERRQ(ierr);} 460e1d83109SMatthew G. Knepley } 4619852e123SBarry Smith for (p = 0; p < size; ++p) { 462e1d83109SMatthew G. Knepley displs[p+1] = displs[p] + Nc[p]; 463e1d83109SMatthew G. Knepley } 464e1d83109SMatthew G. Knepley } 465ffc4695bSBarry Smith ierr = MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, numComponents, MPIU_BOOL, 0, comm);CHKERRMPI(ierr); 46684961fc4SMatthew G. Knepley ierr = PetscFree(flips);CHKERRQ(ierr); 46784961fc4SMatthew G. Knepley } 468dd400576SPatrick Sanan if (rank == 0) {ierr = PetscBTDestroy(&flippedProcs);CHKERRQ(ierr);} 469e1d83109SMatthew G. Knepley ierr = PetscFree(N);CHKERRQ(ierr); 470e1d83109SMatthew G. Knepley ierr = PetscFree4(recvcounts, displs, Nc, Noff);CHKERRQ(ierr); 471e1d83109SMatthew G. Knepley ierr = PetscFree2(nrankComp, match);CHKERRQ(ierr); 472e1d83109SMatthew G. Knepley 473e1d83109SMatthew G. Knepley /* Decide whether to flip cells in each component */ 474e1d83109SMatthew G. Knepley for (c = 0; c < cEnd-cStart; ++c) {if (flipped[cellComp[c]]) {ierr = PetscBTNegate(flippedCells, c);CHKERRQ(ierr);}} 475e1d83109SMatthew G. Knepley ierr = PetscFree(flipped);CHKERRQ(ierr); 47684961fc4SMatthew G. Knepley } 47784961fc4SMatthew G. Knepley if (flg) { 47884961fc4SMatthew G. Knepley PetscViewer v; 47984961fc4SMatthew G. Knepley 4807e09831bSMatthew G. Knepley ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr); 4811575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(v);CHKERRQ(ierr); 48284961fc4SMatthew G. Knepley ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank);CHKERRQ(ierr); 48384961fc4SMatthew G. Knepley ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr); 4844d4c343aSBarry Smith ierr = PetscViewerFlush(v);CHKERRQ(ierr); 4851575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(v);CHKERRQ(ierr); 48684961fc4SMatthew G. Knepley } 48784961fc4SMatthew G. Knepley /* Reverse flipped cells in the mesh */ 48884961fc4SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 4890d366550SMatthew G. Knepley if (PetscBTLookup(flippedCells, c-cStart)) { 490b5a892a1SMatthew G. Knepley ierr = DMPlexOrientPoint(dm, c, -1);CHKERRQ(ierr); 4910d366550SMatthew G. Knepley } 49284961fc4SMatthew G. Knepley } 49384961fc4SMatthew G. Knepley ierr = PetscBTDestroy(&seenCells);CHKERRQ(ierr); 49484961fc4SMatthew G. Knepley ierr = PetscBTDestroy(&flippedCells);CHKERRQ(ierr); 49584961fc4SMatthew G. Knepley ierr = PetscBTDestroy(&seenFaces);CHKERRQ(ierr); 496e1d83109SMatthew G. Knepley ierr = PetscFree2(numNeighbors, neighbors);CHKERRQ(ierr); 497e1d83109SMatthew G. Knepley ierr = PetscFree2(rorntComp, lorntComp);CHKERRQ(ierr); 498e1d83109SMatthew G. Knepley ierr = PetscFree3(faceFIFO, cellComp, faceComp);CHKERRQ(ierr); 49984961fc4SMatthew G. Knepley PetscFunctionReturn(0); 50084961fc4SMatthew G. Knepley } 501