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 16db781477SPatrick Sanan .seealso: `DMPlexOrient()`, `DMPlexGetCone()`, `DMPlexGetConeOrientation()`, `DMPlexInterpolate()`, `DMPlexGetChart()` 1727d0e99aSVaclav Hapla @*/ 189371c9d4SSatish Balay PetscErrorCode DMPlexOrientPoint(DM dm, PetscInt p, PetscInt o) { 19b5a892a1SMatthew G. Knepley DMPolytopeType ct; 20b5a892a1SMatthew G. Knepley const PetscInt *arr, *cone, *ornt, *support; 21b5a892a1SMatthew G. Knepley PetscInt *newcone, *newornt; 22b5a892a1SMatthew G. Knepley PetscInt coneSize, c, supportSize, s; 2327d0e99aSVaclav Hapla 2427d0e99aSVaclav Hapla PetscFunctionBegin; 2527d0e99aSVaclav Hapla PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 269566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, p, &ct)); 27b5a892a1SMatthew G. Knepley arr = DMPolytopeTypeGetArrangment(ct, o); 289566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, p, &coneSize)); 299566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, p, &cone)); 309566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeOrientation(dm, p, &ornt)); 319566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, coneSize, MPIU_INT, &newcone)); 329566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, coneSize, MPIU_INT, &newornt)); 33b5a892a1SMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 34b5a892a1SMatthew G. Knepley DMPolytopeType ft; 35b5a892a1SMatthew G. Knepley PetscInt nO; 3627d0e99aSVaclav Hapla 379566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cone[c], &ft)); 38b5a892a1SMatthew G. Knepley nO = DMPolytopeTypeGetNumArrangments(ft) / 2; 39b5a892a1SMatthew G. Knepley newcone[c] = cone[arr[c * 2 + 0]]; 40b5a892a1SMatthew G. Knepley newornt[c] = DMPolytopeTypeComposeOrientation(ft, arr[c * 2 + 1], ornt[arr[c * 2 + 0]]); 4163a3b9bcSJacob 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]); 4227d0e99aSVaclav Hapla } 439566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(dm, p, newcone)); 449566063dSJacob Faibussowitsch PetscCall(DMPlexSetConeOrientation(dm, p, newornt)); 459566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, coneSize, MPIU_INT, &newcone)); 469566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, coneSize, MPIU_INT, &newornt)); 47b5a892a1SMatthew G. Knepley /* Update orientation of this point in the support points */ 489566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, p, &supportSize)); 499566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, p, &support)); 50b5a892a1SMatthew G. Knepley for (s = 0; s < supportSize; ++s) { 519566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, support[s], &coneSize)); 529566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, support[s], &cone)); 539566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeOrientation(dm, support[s], &ornt)); 54b5a892a1SMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 55b5a892a1SMatthew G. Knepley PetscInt po; 5627d0e99aSVaclav Hapla 57b5a892a1SMatthew G. Knepley if (cone[c] != p) continue; 58b5a892a1SMatthew G. Knepley /* ornt[c] * 0 = target = po * o so that po = ornt[c] * o^{-1} */ 59b5a892a1SMatthew G. Knepley po = DMPolytopeTypeComposeOrientationInv(ct, ornt[c], o); 609566063dSJacob Faibussowitsch PetscCall(DMPlexInsertConeOrientation(dm, support[s], c, po)); 6184961fc4SMatthew G. Knepley } 6284961fc4SMatthew G. Knepley } 6384961fc4SMatthew G. Knepley PetscFunctionReturn(0); 6484961fc4SMatthew G. Knepley } 6584961fc4SMatthew G. Knepley 667b310ce2SMatthew G. Knepley /* 677b310ce2SMatthew G. Knepley - Checks face match 687b310ce2SMatthew G. Knepley - Flips non-matching 697b310ce2SMatthew G. Knepley - Inserts faces of support cells in FIFO 707b310ce2SMatthew G. Knepley */ 719371c9d4SSatish Balay static PetscErrorCode DMPlexCheckFace_Internal(DM dm, PetscInt *faceFIFO, PetscInt *fTop, PetscInt *fBottom, PetscInt cStart, PetscInt fStart, PetscInt fEnd, PetscBT seenCells, PetscBT flippedCells, PetscBT seenFaces) { 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 767b310ce2SMatthew G. Knepley PetscFunctionBegin; 777b310ce2SMatthew G. Knepley face = faceFIFO[(*fTop)++]; 789566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 799566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, face, &supportSize)); 809566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, face, &support)); 817b310ce2SMatthew G. Knepley if (supportSize < 2) PetscFunctionReturn(0); 8263a3b9bcSJacob Faibussowitsch PetscCheck(supportSize == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %" PetscInt_FMT, supportSize); 837b310ce2SMatthew G. Knepley seenA = PetscBTLookup(seenCells, support[0] - cStart); 847b310ce2SMatthew G. Knepley flippedA = PetscBTLookup(flippedCells, support[0] - cStart) ? 1 : 0; 857b310ce2SMatthew G. Knepley seenB = PetscBTLookup(seenCells, support[1] - cStart); 867b310ce2SMatthew G. Knepley flippedB = PetscBTLookup(flippedCells, support[1] - cStart) ? 1 : 0; 877b310ce2SMatthew G. Knepley 889566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, support[0], &coneSizeA)); 899566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, support[1], &coneSizeB)); 909566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, support[0], &coneA)); 919566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, support[1], &coneB)); 929566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeOrientation(dm, support[0], &coneOA)); 939566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeOrientation(dm, support[1], &coneOB)); 947b310ce2SMatthew G. Knepley for (c = 0; c < coneSizeA; ++c) { 957b310ce2SMatthew G. Knepley if (!PetscBTLookup(seenFaces, coneA[c] - fStart)) { 967b310ce2SMatthew G. Knepley faceFIFO[(*fBottom)++] = coneA[c]; 979566063dSJacob Faibussowitsch PetscCall(PetscBTSet(seenFaces, coneA[c] - fStart)); 987b310ce2SMatthew G. Knepley } 997b310ce2SMatthew G. Knepley if (coneA[c] == face) posA = c; 10063a3b9bcSJacob 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); 1017b310ce2SMatthew G. Knepley } 10263a3b9bcSJacob 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]); 1037b310ce2SMatthew G. Knepley for (c = 0; c < coneSizeB; ++c) { 1047b310ce2SMatthew G. Knepley if (!PetscBTLookup(seenFaces, coneB[c] - fStart)) { 1057b310ce2SMatthew G. Knepley faceFIFO[(*fBottom)++] = coneB[c]; 1069566063dSJacob Faibussowitsch PetscCall(PetscBTSet(seenFaces, coneB[c] - fStart)); 1077b310ce2SMatthew G. Knepley } 1087b310ce2SMatthew G. Knepley if (coneB[c] == face) posB = c; 10963a3b9bcSJacob 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); 1107b310ce2SMatthew G. Knepley } 11163a3b9bcSJacob 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]); 1127b310ce2SMatthew G. Knepley 1137b310ce2SMatthew G. Knepley if (dim == 1) { 1147b310ce2SMatthew G. Knepley mismatch = posA == posB; 1157b310ce2SMatthew G. Knepley } else { 1167b310ce2SMatthew G. Knepley mismatch = coneOA[posA] == coneOB[posB]; 1177b310ce2SMatthew G. Knepley } 1187b310ce2SMatthew G. Knepley 1197b310ce2SMatthew G. Knepley if (mismatch ^ (flippedA ^ flippedB)) { 12063a3b9bcSJacob 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]); 1217b310ce2SMatthew G. Knepley if (!seenA && !flippedA) { 1229566063dSJacob Faibussowitsch PetscCall(PetscBTSet(flippedCells, support[0] - cStart)); 1237b310ce2SMatthew G. Knepley } else if (!seenB && !flippedB) { 1249566063dSJacob Faibussowitsch PetscCall(PetscBTSet(flippedCells, support[1] - cStart)); 1257b310ce2SMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable"); 1261dca8a05SBarry Smith } else PetscCheck(!mismatch || !flippedA || !flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable"); 1279566063dSJacob Faibussowitsch PetscCall(PetscBTSet(seenCells, support[0] - cStart)); 1289566063dSJacob Faibussowitsch PetscCall(PetscBTSet(seenCells, support[1] - cStart)); 1297b310ce2SMatthew G. Knepley PetscFunctionReturn(0); 1307b310ce2SMatthew G. Knepley } 1317b310ce2SMatthew G. Knepley 13284961fc4SMatthew G. Knepley /*@ 13384961fc4SMatthew G. Knepley DMPlexOrient - Give a consistent orientation to the input mesh 13484961fc4SMatthew G. Knepley 13584961fc4SMatthew G. Knepley Input Parameters: 13684961fc4SMatthew G. Knepley . dm - The DM 13784961fc4SMatthew G. Knepley 13884961fc4SMatthew G. Knepley Note: The orientation data for the DM are change in-place. 13984961fc4SMatthew G. Knepley $ This routine will fail for non-orientable surfaces, such as the Moebius strip. 14084961fc4SMatthew G. Knepley 14184961fc4SMatthew G. Knepley Level: advanced 14284961fc4SMatthew G. Knepley 143db781477SPatrick Sanan .seealso: `DMCreate()`, `DMPLEX` 14484961fc4SMatthew G. Knepley @*/ 1459371c9d4SSatish Balay PetscErrorCode DMPlexOrient(DM dm) { 14684961fc4SMatthew G. Knepley MPI_Comm comm; 147e1d83109SMatthew G. Knepley PetscSF sf; 148e1d83109SMatthew G. Knepley const PetscInt *lpoints; 149e1d83109SMatthew G. Knepley const PetscSFNode *rpoints; 150e1d83109SMatthew G. Knepley PetscSFNode *rorntComp = NULL, *lorntComp = NULL; 1513d6051bcSMatthew G. Knepley PetscInt *numNeighbors, **neighbors, *locSupport = NULL; 152e1d83109SMatthew G. Knepley PetscSFNode *nrankComp; 153e1d83109SMatthew G. Knepley PetscBool *match, *flipped; 15484961fc4SMatthew G. Knepley PetscBT seenCells, flippedCells, seenFaces; 155e1d83109SMatthew G. Knepley PetscInt *faceFIFO, fTop, fBottom, *cellComp, *faceComp; 1567cadcfe8SMatthew G. Knepley PetscInt numLeaves, numRoots, dim, h, cStart, cEnd, c, cell, fStart, fEnd, face, off, totNeighbors = 0; 15739526728SToby Isaac PetscMPIInt rank, size, numComponents, comp = 0; 15839526728SToby Isaac PetscBool flg, flg2; 159a17aa47eSToby Isaac PetscViewer viewer = NULL, selfviewer = NULL; 16084961fc4SMatthew G. Knepley 16184961fc4SMatthew G. Knepley PetscFunctionBegin; 1629566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 1639566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1649566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 1659566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-orientation_view", &flg)); 1669566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-orientation_view_synchronized", &flg2)); 1679566063dSJacob Faibussowitsch PetscCall(DMGetPointSF(dm, &sf)); 1689566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints)); 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 */ 1789566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1799566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(dm, &h)); 1809566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, h, &cStart, &cEnd)); 1819566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, h + 1, &fStart, &fEnd)); 1829566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(cEnd - cStart, &seenCells)); 1839566063dSJacob Faibussowitsch PetscCall(PetscBTMemzero(cEnd - cStart, seenCells)); 1849566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(cEnd - cStart, &flippedCells)); 1859566063dSJacob Faibussowitsch PetscCall(PetscBTMemzero(cEnd - cStart, flippedCells)); 1869566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(fEnd - fStart, &seenFaces)); 1879566063dSJacob Faibussowitsch PetscCall(PetscBTMemzero(fEnd - fStart, seenFaces)); 1889566063dSJacob Faibussowitsch PetscCall(PetscCalloc3(fEnd - fStart, &faceFIFO, cEnd - cStart, &cellComp, fEnd - fStart, &faceComp)); 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 */ 2129371c9d4SSatish Balay for (cell = cStart; cell < cEnd; ++cell) 2139371c9d4SSatish Balay if (cellComp[cell - cStart] < 0) break; 214e1d83109SMatthew G. Knepley if (cell >= cEnd) break; 215e1d83109SMatthew G. Knepley /* Initialize FIFO with first cell in component */ 216e1d83109SMatthew G. Knepley { 21784961fc4SMatthew G. Knepley const PetscInt *cone; 21884961fc4SMatthew G. Knepley PetscInt coneSize; 21984961fc4SMatthew G. Knepley 220e1d83109SMatthew G. Knepley fTop = fBottom = 0; 2219566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 2229566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, cell, &cone)); 22384961fc4SMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 22484961fc4SMatthew G. Knepley faceFIFO[fBottom++] = cone[c]; 2259566063dSJacob Faibussowitsch PetscCall(PetscBTSet(seenFaces, cone[c] - fStart)); 22684961fc4SMatthew G. Knepley } 2279566063dSJacob Faibussowitsch PetscCall(PetscBTSet(seenCells, cell - cStart)); 22884961fc4SMatthew G. Knepley } 22984961fc4SMatthew G. Knepley /* Consider each face in FIFO */ 230*48a46eb9SPierre Jolivet while (fTop < fBottom) PetscCall(DMPlexCheckFace_Internal(dm, faceFIFO, &fTop, &fBottom, cStart, fStart, fEnd, seenCells, flippedCells, seenFaces)); 231e1d83109SMatthew G. Knepley /* Set component for cells and faces */ 232e1d83109SMatthew G. Knepley for (cell = 0; cell < cEnd - cStart; ++cell) { 233e1d83109SMatthew G. Knepley if (PetscBTLookup(seenCells, cell)) cellComp[cell] = comp; 234e1d83109SMatthew G. Knepley } 235e1d83109SMatthew G. Knepley for (face = 0; face < fEnd - fStart; ++face) { 236e1d83109SMatthew G. Knepley if (PetscBTLookup(seenFaces, face)) faceComp[face] = comp; 237e1d83109SMatthew G. Knepley } 238e1d83109SMatthew G. Knepley /* Wipe seenCells and seenFaces for next component */ 2399566063dSJacob Faibussowitsch PetscCall(PetscBTMemzero(fEnd - fStart, seenFaces)); 2409566063dSJacob Faibussowitsch PetscCall(PetscBTMemzero(cEnd - cStart, seenCells)); 241e1d83109SMatthew G. Knepley ++comp; 242e1d83109SMatthew G. Knepley } while (1); 243e1d83109SMatthew G. Knepley numComponents = comp; 24484961fc4SMatthew G. Knepley if (flg) { 24584961fc4SMatthew G. Knepley PetscViewer v; 24684961fc4SMatthew G. Knepley 2479566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetStdout(comm, &v)); 2489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(v)); 2499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank)); 2509566063dSJacob Faibussowitsch PetscCall(PetscBTView(cEnd - cStart, flippedCells, v)); 2519566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(v)); 2529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(v)); 25384961fc4SMatthew G. Knepley } 25484961fc4SMatthew G. Knepley /* Now all subdomains are oriented, but we need a consistent parallel orientation */ 25584961fc4SMatthew G. Knepley if (numLeaves >= 0) { 2563d6051bcSMatthew G. Knepley PetscInt maxSupportSize, neighbor; 2573d6051bcSMatthew G. Knepley 258e1d83109SMatthew G. Knepley /* Store orientations of boundary faces*/ 2593d6051bcSMatthew G. Knepley PetscCall(DMPlexGetMaxSizes(dm, NULL, &maxSupportSize)); 2603d6051bcSMatthew G. Knepley PetscCall(PetscCalloc3(numRoots, &rorntComp, numRoots, &lorntComp, maxSupportSize, &locSupport)); 26184961fc4SMatthew G. Knepley for (face = fStart; face < fEnd; ++face) { 262e1d83109SMatthew G. Knepley const PetscInt *cone, *support, *ornt; 2633d6051bcSMatthew G. Knepley PetscInt coneSize, supportSize, Ns = 0, s, l; 264e1d83109SMatthew G. Knepley 2659566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, face, &supportSize)); 2663d6051bcSMatthew G. Knepley /* Ignore overlapping cells */ 2679566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, face, &support)); 2683d6051bcSMatthew G. Knepley for (s = 0; s < supportSize; ++s) { 2693d6051bcSMatthew G. Knepley PetscCall(PetscFindInt(support[s], numLeaves, lpoints, &l)); 2703d6051bcSMatthew G. Knepley if (l >= 0) continue; 2713d6051bcSMatthew G. Knepley locSupport[Ns++] = support[s]; 2723d6051bcSMatthew G. Knepley } 2733d6051bcSMatthew G. Knepley if (Ns != 1) continue; 2743d6051bcSMatthew G. Knepley neighbor = locSupport[0]; 2753d6051bcSMatthew G. Knepley PetscCall(DMPlexGetCone(dm, neighbor, &cone)); 2763d6051bcSMatthew G. Knepley PetscCall(DMPlexGetConeSize(dm, neighbor, &coneSize)); 2773d6051bcSMatthew G. Knepley PetscCall(DMPlexGetConeOrientation(dm, neighbor, &ornt)); 2789371c9d4SSatish Balay for (c = 0; c < coneSize; ++c) 2799371c9d4SSatish Balay if (cone[c] == face) break; 28084961fc4SMatthew G. Knepley if (dim == 1) { 28184961fc4SMatthew G. Knepley /* Use cone position instead, shifted to -1 or 1 */ 2823d6051bcSMatthew G. Knepley if (PetscBTLookup(flippedCells, neighbor - cStart)) rorntComp[face].rank = 1 - c * 2; 2836e1eb25cSMatthew G. Knepley else rorntComp[face].rank = c * 2 - 1; 28484961fc4SMatthew G. Knepley } else { 2853d6051bcSMatthew G. Knepley if (PetscBTLookup(flippedCells, neighbor - cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 : 1; 286e1d83109SMatthew G. Knepley else rorntComp[face].rank = ornt[c] < 0 ? 1 : -1; 28784961fc4SMatthew G. Knepley } 288e1d83109SMatthew G. Knepley rorntComp[face].index = faceComp[face - fStart]; 28984961fc4SMatthew G. Knepley } 290e1d83109SMatthew G. Knepley /* Communicate boundary edge orientations */ 2919566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf, MPIU_2INT, rorntComp, lorntComp, MPI_REPLACE)); 2929566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf, MPIU_2INT, rorntComp, lorntComp, MPI_REPLACE)); 293e1d83109SMatthew G. Knepley } 294e1d83109SMatthew G. Knepley /* Get process adjacency */ 2959566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors)); 296a17aa47eSToby Isaac viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm)); 2979566063dSJacob Faibussowitsch if (flg2) PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 2989566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &selfviewer)); 299e1d83109SMatthew G. Knepley for (comp = 0; comp < numComponents; ++comp) { 300e1d83109SMatthew G. Knepley PetscInt l, n; 30184961fc4SMatthew G. Knepley 302e1d83109SMatthew G. Knepley numNeighbors[comp] = 0; 3039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp])); 304e1d83109SMatthew G. Knepley /* I know this is p^2 time in general, but for bounded degree its alright */ 305e1d83109SMatthew G. Knepley for (l = 0; l < numLeaves; ++l) { 306e1d83109SMatthew G. Knepley const PetscInt face = lpoints[l]; 307e1d83109SMatthew G. Knepley 308e1d83109SMatthew G. Knepley /* Find a representative face (edge) separating pairs of procs */ 3093d6051bcSMatthew G. Knepley if ((face >= fStart) && (face < fEnd) && (faceComp[face - fStart] == comp) && rorntComp[face].rank) { 310269a49d5SMatthew G. Knepley const PetscInt rrank = rpoints[l].rank; 311e1d83109SMatthew G. Knepley const PetscInt rcomp = lorntComp[face].index; 312e1d83109SMatthew G. Knepley 3139371c9d4SSatish Balay for (n = 0; n < numNeighbors[comp]; ++n) 3149371c9d4SSatish Balay if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break; 315e1d83109SMatthew G. Knepley if (n >= numNeighbors[comp]) { 316e1d83109SMatthew G. Knepley PetscInt supportSize; 317e1d83109SMatthew G. Knepley 3189566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, face, &supportSize)); 31963a3b9bcSJacob Faibussowitsch PetscCheck(supportSize == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %" PetscInt_FMT, supportSize); 3209371c9d4SSatish Balay if (flg) 3219371c9d4SSatish Balay 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, 3229371c9d4SSatish Balay 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)); 3639371c9d4SSatish Balay for (p = 0; p < size; ++p) { displs[p + 1] = displs[p] + Nc[p]; } 3649566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscMalloc1(displs[size], &N)); 3659566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gatherv(numNeighbors, numComponents, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm)); 3669852e123SBarry Smith for (p = 0, o = 0; p < size; ++p) { 367e1d83109SMatthew G. Knepley recvcounts[p] = 0; 368e1d83109SMatthew G. Knepley for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o]; 36984961fc4SMatthew G. Knepley displs[p + 1] = displs[p] + recvcounts[p]; 37084961fc4SMatthew G. Knepley } 3719566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscMalloc2(displs[size], &adj, displs[size], &val)); 3729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gatherv(nrankComp, totNeighbors, MPIU_2INT, adj, recvcounts, displs, MPIU_2INT, 0, comm)); 3739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gatherv(match, totNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm)); 3749566063dSJacob Faibussowitsch PetscCall(PetscFree2(numNeighbors, neighbors)); 375dd400576SPatrick Sanan if (rank == 0) { 3769852e123SBarry Smith for (p = 1; p <= size; ++p) { Noff[p] = Noff[p - 1] + Nc[p - 1]; } 377a83982f3SMatthew G. Knepley if (flg) { 378e1d83109SMatthew G. Knepley PetscInt n; 379e1d83109SMatthew G. Knepley 3809852e123SBarry Smith for (p = 0, off = 0; p < size; ++p) { 381e1d83109SMatthew G. Knepley for (c = 0; c < Nc[p]; ++c) { 38263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %" PetscInt_FMT ":\n", p, c)); 383*48a46eb9SPierre Jolivet for (n = 0; n < N[Noff[p] + c]; ++n, ++off) PetscCall(PetscPrintf(PETSC_COMM_SELF, " edge (%" PetscInt_FMT ", %" PetscInt_FMT ") (%s):\n", adj[off].rank, adj[off].index, PetscBools[val[off]])); 38484961fc4SMatthew G. Knepley } 38584961fc4SMatthew G. Knepley } 38684961fc4SMatthew G. Knepley } 38784961fc4SMatthew G. Knepley /* Symmetrize the graph */ 3889566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &G)); 3899566063dSJacob Faibussowitsch PetscCall(MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size])); 3909566063dSJacob Faibussowitsch PetscCall(MatSetUp(G)); 3919852e123SBarry Smith for (p = 0, off = 0; p < size; ++p) { 392e1d83109SMatthew G. Knepley for (c = 0; c < Nc[p]; ++c) { 393e1d83109SMatthew G. Knepley const PetscInt r = Noff[p] + c; 394e1d83109SMatthew G. Knepley PetscInt n; 395e1d83109SMatthew G. Knepley 396e1d83109SMatthew G. Knepley for (n = 0; n < N[r]; ++n, ++off) { 397e1d83109SMatthew G. Knepley const PetscInt q = Noff[adj[off].rank] + adj[off].index; 398e1d83109SMatthew G. Knepley const PetscScalar o = val[off] ? 1.0 : 0.0; 39984961fc4SMatthew G. Knepley 4009566063dSJacob Faibussowitsch PetscCall(MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES)); 4019566063dSJacob Faibussowitsch PetscCall(MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES)); 40284961fc4SMatthew G. Knepley } 40384961fc4SMatthew G. Knepley } 404e1d83109SMatthew G. Knepley } 4059566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY)); 4069566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY)); 40784961fc4SMatthew G. Knepley 4089566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(Noff[size], &seenProcs)); 4099566063dSJacob Faibussowitsch PetscCall(PetscBTMemzero(Noff[size], seenProcs)); 4109566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(Noff[size], &flippedProcs)); 4119566063dSJacob Faibussowitsch PetscCall(PetscBTMemzero(Noff[size], flippedProcs)); 4129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Noff[size], &procFIFO)); 41384961fc4SMatthew G. Knepley pTop = pBottom = 0; 4149852e123SBarry Smith for (p = 0; p < Noff[size]; ++p) { 41584961fc4SMatthew G. Knepley if (PetscBTLookup(seenProcs, p)) continue; 41684961fc4SMatthew G. Knepley /* Initialize FIFO with next proc */ 41784961fc4SMatthew G. Knepley procFIFO[pBottom++] = p; 4189566063dSJacob Faibussowitsch PetscCall(PetscBTSet(seenProcs, p)); 41984961fc4SMatthew G. Knepley /* Consider each proc in FIFO */ 42084961fc4SMatthew G. Knepley while (pTop < pBottom) { 42184961fc4SMatthew G. Knepley const PetscScalar *ornt; 42284961fc4SMatthew G. Knepley const PetscInt *neighbors; 423e1d83109SMatthew G. Knepley PetscInt proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n; 42484961fc4SMatthew G. Knepley 42584961fc4SMatthew G. Knepley proc = procFIFO[pTop++]; 42684961fc4SMatthew G. Knepley flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0; 4279566063dSJacob Faibussowitsch PetscCall(MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt)); 42884961fc4SMatthew G. Knepley /* Loop over neighboring procs */ 42984961fc4SMatthew G. Knepley for (n = 0; n < numNeighbors; ++n) { 43084961fc4SMatthew G. Knepley nproc = neighbors[n]; 43184961fc4SMatthew G. Knepley mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1; 43284961fc4SMatthew G. Knepley seen = PetscBTLookup(seenProcs, nproc); 43384961fc4SMatthew G. Knepley flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0; 43484961fc4SMatthew G. Knepley 43584961fc4SMatthew G. Knepley if (mismatch ^ (flippedA ^ flippedB)) { 43663a3b9bcSJacob 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); 43784961fc4SMatthew G. Knepley if (!flippedB) { 4389566063dSJacob Faibussowitsch PetscCall(PetscBTSet(flippedProcs, nproc)); 43984961fc4SMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable"); 4401dca8a05SBarry Smith } else PetscCheck(!mismatch || !flippedA || !flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable"); 44184961fc4SMatthew G. Knepley if (!seen) { 44284961fc4SMatthew G. Knepley procFIFO[pBottom++] = nproc; 4439566063dSJacob Faibussowitsch PetscCall(PetscBTSet(seenProcs, nproc)); 44484961fc4SMatthew G. Knepley } 44584961fc4SMatthew G. Knepley } 44684961fc4SMatthew G. Knepley } 44784961fc4SMatthew G. Knepley } 4489566063dSJacob Faibussowitsch PetscCall(PetscFree(procFIFO)); 4499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&G)); 4509566063dSJacob Faibussowitsch PetscCall(PetscFree2(adj, val)); 4519566063dSJacob Faibussowitsch PetscCall(PetscBTDestroy(&seenProcs)); 45284961fc4SMatthew G. Knepley } 453e1d83109SMatthew G. Knepley /* Scatter flip flags */ 454e1d83109SMatthew G. Knepley { 455e1d83109SMatthew G. Knepley PetscBool *flips = NULL; 456e1d83109SMatthew G. Knepley 457dd400576SPatrick Sanan if (rank == 0) { 4589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Noff[size], &flips)); 4599852e123SBarry Smith for (p = 0; p < Noff[size]; ++p) { 460e1d83109SMatthew G. Knepley flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE; 4619566063dSJacob Faibussowitsch if (flg && flips[p]) PetscCall(PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p)); 462e1d83109SMatthew G. Knepley } 4639371c9d4SSatish Balay for (p = 0; p < size; ++p) { displs[p + 1] = displs[p] + Nc[p]; } 464e1d83109SMatthew G. Knepley } 4659566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, numComponents, MPIU_BOOL, 0, comm)); 4669566063dSJacob Faibussowitsch PetscCall(PetscFree(flips)); 46784961fc4SMatthew G. Knepley } 4689566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscBTDestroy(&flippedProcs)); 4699566063dSJacob Faibussowitsch PetscCall(PetscFree(N)); 4709566063dSJacob Faibussowitsch PetscCall(PetscFree4(recvcounts, displs, Nc, Noff)); 4719566063dSJacob Faibussowitsch PetscCall(PetscFree2(nrankComp, match)); 472e1d83109SMatthew G. Knepley 473e1d83109SMatthew G. Knepley /* Decide whether to flip cells in each component */ 4749371c9d4SSatish Balay for (c = 0; c < cEnd - cStart; ++c) { 4759371c9d4SSatish Balay if (flipped[cellComp[c]]) PetscCall(PetscBTNegate(flippedCells, c)); 4769371c9d4SSatish Balay } 4779566063dSJacob Faibussowitsch PetscCall(PetscFree(flipped)); 47884961fc4SMatthew G. Knepley } 47984961fc4SMatthew G. Knepley if (flg) { 48084961fc4SMatthew G. Knepley PetscViewer v; 48184961fc4SMatthew G. Knepley 4829566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetStdout(comm, &v)); 4839566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(v)); 4849566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank)); 4859566063dSJacob Faibussowitsch PetscCall(PetscBTView(cEnd - cStart, flippedCells, v)); 4869566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(v)); 4879566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(v)); 48884961fc4SMatthew G. Knepley } 48984961fc4SMatthew G. Knepley /* Reverse flipped cells in the mesh */ 49084961fc4SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 491*48a46eb9SPierre Jolivet if (PetscBTLookup(flippedCells, c - cStart)) PetscCall(DMPlexOrientPoint(dm, c, -1)); 49284961fc4SMatthew G. Knepley } 4939566063dSJacob Faibussowitsch PetscCall(PetscBTDestroy(&seenCells)); 4949566063dSJacob Faibussowitsch PetscCall(PetscBTDestroy(&flippedCells)); 4959566063dSJacob Faibussowitsch PetscCall(PetscBTDestroy(&seenFaces)); 4969566063dSJacob Faibussowitsch PetscCall(PetscFree2(numNeighbors, neighbors)); 4973d6051bcSMatthew G. Knepley PetscCall(PetscFree3(rorntComp, lorntComp, locSupport)); 4989566063dSJacob Faibussowitsch PetscCall(PetscFree3(faceFIFO, cellComp, faceComp)); 49984961fc4SMatthew G. Knepley PetscFunctionReturn(0); 50084961fc4SMatthew G. Knepley } 501