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: 10*a1cb98faSBarry Smith + dm - The `DM` 11b5a892a1SMatthew G. Knepley . p - The mesh point 12b5a892a1SMatthew G. Knepley - o - The orientation 1327d0e99aSVaclav Hapla 1427d0e99aSVaclav Hapla Level: intermediate 1527d0e99aSVaclav Hapla 16*a1cb98faSBarry Smith .seealso: [](chapter_unstructured), `DM`, `DMPLEX`, `DMPlexOrient()`, `DMPlexGetCone()`, `DMPlexGetConeOrientation()`, `DMPlexInterpolate()`, `DMPlexGetChart()` 1727d0e99aSVaclav Hapla @*/ 18d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexOrientPoint(DM dm, PetscInt p, PetscInt o) 19d71ae5a4SJacob Faibussowitsch { 20b5a892a1SMatthew G. Knepley DMPolytopeType ct; 21b5a892a1SMatthew G. Knepley const PetscInt *arr, *cone, *ornt, *support; 22b5a892a1SMatthew G. Knepley PetscInt *newcone, *newornt; 23b5a892a1SMatthew G. Knepley PetscInt coneSize, c, supportSize, s; 2427d0e99aSVaclav Hapla 2527d0e99aSVaclav Hapla PetscFunctionBegin; 2627d0e99aSVaclav Hapla PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 279566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, p, &ct)); 28b5a892a1SMatthew G. Knepley arr = DMPolytopeTypeGetArrangment(ct, o); 299566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, p, &coneSize)); 309566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, p, &cone)); 319566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeOrientation(dm, p, &ornt)); 329566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, coneSize, MPIU_INT, &newcone)); 339566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, coneSize, MPIU_INT, &newornt)); 34b5a892a1SMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 35b5a892a1SMatthew G. Knepley DMPolytopeType ft; 36b5a892a1SMatthew G. Knepley PetscInt nO; 3727d0e99aSVaclav Hapla 389566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cone[c], &ft)); 39b5a892a1SMatthew G. Knepley nO = DMPolytopeTypeGetNumArrangments(ft) / 2; 40b5a892a1SMatthew G. Knepley newcone[c] = cone[arr[c * 2 + 0]]; 41b5a892a1SMatthew G. Knepley newornt[c] = DMPolytopeTypeComposeOrientation(ft, arr[c * 2 + 1], ornt[arr[c * 2 + 0]]); 4263a3b9bcSJacob Faibussowitsch PetscCheck(!newornt[c] || !(newornt[c] >= nO || newornt[c] < -nO), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid orientation %" PetscInt_FMT " not in [%" PetscInt_FMT ",%" PetscInt_FMT ") for %s %" PetscInt_FMT, newornt[c], -nO, nO, DMPolytopeTypes[ft], cone[c]); 4327d0e99aSVaclav Hapla } 449566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(dm, p, newcone)); 459566063dSJacob Faibussowitsch PetscCall(DMPlexSetConeOrientation(dm, p, newornt)); 469566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, coneSize, MPIU_INT, &newcone)); 479566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, coneSize, MPIU_INT, &newornt)); 48b5a892a1SMatthew G. Knepley /* Update orientation of this point in the support points */ 499566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, p, &supportSize)); 509566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, p, &support)); 51b5a892a1SMatthew G. Knepley for (s = 0; s < supportSize; ++s) { 529566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, support[s], &coneSize)); 539566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, support[s], &cone)); 549566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeOrientation(dm, support[s], &ornt)); 55b5a892a1SMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 56b5a892a1SMatthew G. Knepley PetscInt po; 5727d0e99aSVaclav Hapla 58b5a892a1SMatthew G. Knepley if (cone[c] != p) continue; 59b5a892a1SMatthew G. Knepley /* ornt[c] * 0 = target = po * o so that po = ornt[c] * o^{-1} */ 60b5a892a1SMatthew G. Knepley po = DMPolytopeTypeComposeOrientationInv(ct, ornt[c], o); 619566063dSJacob Faibussowitsch PetscCall(DMPlexInsertConeOrientation(dm, support[s], c, po)); 6284961fc4SMatthew G. Knepley } 6384961fc4SMatthew G. Knepley } 6484961fc4SMatthew G. Knepley PetscFunctionReturn(0); 6584961fc4SMatthew G. Knepley } 6684961fc4SMatthew G. Knepley 677b310ce2SMatthew G. Knepley /* 687b310ce2SMatthew G. Knepley - Checks face match 697b310ce2SMatthew G. Knepley - Flips non-matching 707b310ce2SMatthew G. Knepley - Inserts faces of support cells in FIFO 717b310ce2SMatthew G. Knepley */ 72d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCheckFace_Internal(DM dm, PetscInt *faceFIFO, PetscInt *fTop, PetscInt *fBottom, PetscInt cStart, PetscInt fStart, PetscInt fEnd, PetscBT seenCells, PetscBT flippedCells, PetscBT seenFaces) 73d71ae5a4SJacob Faibussowitsch { 747b310ce2SMatthew G. Knepley const PetscInt *support, *coneA, *coneB, *coneOA, *coneOB; 757b310ce2SMatthew G. Knepley PetscInt supportSize, coneSizeA, coneSizeB, posA = -1, posB = -1; 767b2de0fdSMatthew G. Knepley PetscInt face, dim, seenA, flippedA, seenB, flippedB, mismatch, c; 777b310ce2SMatthew G. Knepley 787b310ce2SMatthew G. Knepley PetscFunctionBegin; 797b310ce2SMatthew G. Knepley face = faceFIFO[(*fTop)++]; 809566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 819566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, face, &supportSize)); 829566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, face, &support)); 837b310ce2SMatthew G. Knepley if (supportSize < 2) PetscFunctionReturn(0); 8463a3b9bcSJacob Faibussowitsch PetscCheck(supportSize == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %" PetscInt_FMT, supportSize); 857b310ce2SMatthew G. Knepley seenA = PetscBTLookup(seenCells, support[0] - cStart); 867b310ce2SMatthew G. Knepley flippedA = PetscBTLookup(flippedCells, support[0] - cStart) ? 1 : 0; 877b310ce2SMatthew G. Knepley seenB = PetscBTLookup(seenCells, support[1] - cStart); 887b310ce2SMatthew G. Knepley flippedB = PetscBTLookup(flippedCells, support[1] - cStart) ? 1 : 0; 897b310ce2SMatthew G. Knepley 909566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, support[0], &coneSizeA)); 919566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, support[1], &coneSizeB)); 929566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, support[0], &coneA)); 939566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, support[1], &coneB)); 949566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeOrientation(dm, support[0], &coneOA)); 959566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeOrientation(dm, support[1], &coneOB)); 967b310ce2SMatthew G. Knepley for (c = 0; c < coneSizeA; ++c) { 977b310ce2SMatthew G. Knepley if (!PetscBTLookup(seenFaces, coneA[c] - fStart)) { 987b310ce2SMatthew G. Knepley faceFIFO[(*fBottom)++] = coneA[c]; 999566063dSJacob Faibussowitsch PetscCall(PetscBTSet(seenFaces, coneA[c] - fStart)); 1007b310ce2SMatthew G. Knepley } 1017b310ce2SMatthew G. Knepley if (coneA[c] == face) posA = c; 10263a3b9bcSJacob Faibussowitsch PetscCheck(*fBottom <= fEnd - fStart, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %" PetscInt_FMT " was pushed exceeding capacity %" PetscInt_FMT " > %" PetscInt_FMT, coneA[c], *fBottom, fEnd - fStart); 1037b310ce2SMatthew G. Knepley } 10463a3b9bcSJacob Faibussowitsch PetscCheck(posA >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " could not be located in cell %" PetscInt_FMT, face, support[0]); 1057b310ce2SMatthew G. Knepley for (c = 0; c < coneSizeB; ++c) { 1067b310ce2SMatthew G. Knepley if (!PetscBTLookup(seenFaces, coneB[c] - fStart)) { 1077b310ce2SMatthew G. Knepley faceFIFO[(*fBottom)++] = coneB[c]; 1089566063dSJacob Faibussowitsch PetscCall(PetscBTSet(seenFaces, coneB[c] - fStart)); 1097b310ce2SMatthew G. Knepley } 1107b310ce2SMatthew G. Knepley if (coneB[c] == face) posB = c; 11163a3b9bcSJacob Faibussowitsch PetscCheck(*fBottom <= fEnd - fStart, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %" PetscInt_FMT " was pushed exceeding capacity %" PetscInt_FMT " > %" PetscInt_FMT, coneA[c], *fBottom, fEnd - fStart); 1127b310ce2SMatthew G. Knepley } 11363a3b9bcSJacob Faibussowitsch PetscCheck(posB >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " could not be located in cell %" PetscInt_FMT, face, support[1]); 1147b310ce2SMatthew G. Knepley 1157b310ce2SMatthew G. Knepley if (dim == 1) { 1167b310ce2SMatthew G. Knepley mismatch = posA == posB; 1177b310ce2SMatthew G. Knepley } else { 1187b310ce2SMatthew G. Knepley mismatch = coneOA[posA] == coneOB[posB]; 1197b310ce2SMatthew G. Knepley } 1207b310ce2SMatthew G. Knepley 1217b310ce2SMatthew G. Knepley if (mismatch ^ (flippedA ^ flippedB)) { 12263a3b9bcSJacob Faibussowitsch PetscCheck(!seenA || !seenB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Previously seen cells %" PetscInt_FMT " and %" PetscInt_FMT " do not match: Fault mesh is non-orientable", support[0], support[1]); 1237b310ce2SMatthew G. Knepley if (!seenA && !flippedA) { 1249566063dSJacob Faibussowitsch PetscCall(PetscBTSet(flippedCells, support[0] - cStart)); 1257b310ce2SMatthew G. Knepley } else if (!seenB && !flippedB) { 1269566063dSJacob Faibussowitsch PetscCall(PetscBTSet(flippedCells, support[1] - cStart)); 1277b310ce2SMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable"); 1281dca8a05SBarry Smith } else PetscCheck(!mismatch || !flippedA || !flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable"); 1299566063dSJacob Faibussowitsch PetscCall(PetscBTSet(seenCells, support[0] - cStart)); 1309566063dSJacob Faibussowitsch PetscCall(PetscBTSet(seenCells, support[1] - cStart)); 1317b310ce2SMatthew G. Knepley PetscFunctionReturn(0); 1327b310ce2SMatthew G. Knepley } 1337b310ce2SMatthew G. Knepley 13484961fc4SMatthew G. Knepley /*@ 13584961fc4SMatthew G. Knepley DMPlexOrient - Give a consistent orientation to the input mesh 13684961fc4SMatthew G. Knepley 13784961fc4SMatthew G. Knepley Input Parameters: 138*a1cb98faSBarry Smith . dm - The `DM` 13984961fc4SMatthew G. Knepley 140*a1cb98faSBarry Smith Note: 141*a1cb98faSBarry Smith The orientation data for the `DM` are change in-place. 142*a1cb98faSBarry Smith 143*a1cb98faSBarry Smith This routine will fail for non-orientable surfaces, such as the Moebius strip. 14484961fc4SMatthew G. Knepley 14584961fc4SMatthew G. Knepley Level: advanced 14684961fc4SMatthew G. Knepley 147*a1cb98faSBarry Smith .seealso: [](chapter_unstructured), `DM`, `DMPLEX`, `DMCreate()`, `DMPLEX` 14884961fc4SMatthew G. Knepley @*/ 149d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexOrient(DM dm) 150d71ae5a4SJacob Faibussowitsch { 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; 1563d6051bcSMatthew G. Knepley PetscInt *numNeighbors, **neighbors, *locSupport = NULL; 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 16684961fc4SMatthew G. Knepley PetscFunctionBegin; 1679566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 1689566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1699566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 1709566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-orientation_view", &flg)); 1719566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-orientation_view_synchronized", &flg2)); 1729566063dSJacob Faibussowitsch PetscCall(DMGetPointSF(dm, &sf)); 1739566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints)); 17484961fc4SMatthew G. Knepley /* Truth Table 17584961fc4SMatthew G. Knepley mismatch flips do action mismatch flipA ^ flipB action 17684961fc4SMatthew G. Knepley F 0 flips no F F F 17784961fc4SMatthew G. Knepley F 1 flip yes F T T 17884961fc4SMatthew G. Knepley F 2 flips no T F T 17984961fc4SMatthew G. Knepley T 0 flips yes T T F 18084961fc4SMatthew G. Knepley T 1 flip no 18184961fc4SMatthew G. Knepley T 2 flips yes 18284961fc4SMatthew G. Knepley */ 1839566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1849566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(dm, &h)); 1859566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, h, &cStart, &cEnd)); 1869566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, h + 1, &fStart, &fEnd)); 1879566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(cEnd - cStart, &seenCells)); 1889566063dSJacob Faibussowitsch PetscCall(PetscBTMemzero(cEnd - cStart, seenCells)); 1899566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(cEnd - cStart, &flippedCells)); 1909566063dSJacob Faibussowitsch PetscCall(PetscBTMemzero(cEnd - cStart, flippedCells)); 1919566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(fEnd - fStart, &seenFaces)); 1929566063dSJacob Faibussowitsch PetscCall(PetscBTMemzero(fEnd - fStart, seenFaces)); 1939566063dSJacob Faibussowitsch PetscCall(PetscCalloc3(fEnd - fStart, &faceFIFO, cEnd - cStart, &cellComp, fEnd - fStart, &faceComp)); 194e1d83109SMatthew G. Knepley /* 195e1d83109SMatthew G. Knepley OLD STYLE 196e1d83109SMatthew G. Knepley - Add an integer array over cells and faces (component) for connected component number 197e1d83109SMatthew G. Knepley Foreach component 198e1d83109SMatthew G. Knepley - Mark the initial cell as seen 199e1d83109SMatthew G. Knepley - Process component as usual 200e1d83109SMatthew G. Knepley - Set component for all seenCells 201e1d83109SMatthew G. Knepley - Wipe seenCells and seenFaces (flippedCells can stay) 202e1d83109SMatthew G. Knepley - Generate parallel adjacency for component using SF and seenFaces 203e1d83109SMatthew G. Knepley - Collect numComponents adj data from each proc to 0 204e1d83109SMatthew G. Knepley - Build same serial graph 205e1d83109SMatthew G. Knepley - Use same solver 206e1d83109SMatthew G. Knepley - Use Scatterv to to send back flipped flags for each component 207e1d83109SMatthew G. Knepley - Negate flippedCells by component 208e1d83109SMatthew G. Knepley 209e1d83109SMatthew G. Knepley NEW STYLE 210e1d83109SMatthew G. Knepley - Create the adj on each process 211e1d83109SMatthew G. Knepley - Bootstrap to complete graph on proc 0 212e1d83109SMatthew G. Knepley */ 213e1d83109SMatthew G. Knepley /* Loop over components */ 214e1d83109SMatthew G. Knepley for (cell = cStart; cell < cEnd; ++cell) cellComp[cell - cStart] = -1; 215e1d83109SMatthew G. Knepley do { 216e1d83109SMatthew G. Knepley /* Look for first unmarked cell */ 2179371c9d4SSatish Balay for (cell = cStart; cell < cEnd; ++cell) 2189371c9d4SSatish Balay 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; 2269566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 2279566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, cell, &cone)); 22884961fc4SMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 22984961fc4SMatthew G. Knepley faceFIFO[fBottom++] = cone[c]; 2309566063dSJacob Faibussowitsch PetscCall(PetscBTSet(seenFaces, cone[c] - fStart)); 23184961fc4SMatthew G. Knepley } 2329566063dSJacob Faibussowitsch PetscCall(PetscBTSet(seenCells, cell - cStart)); 23384961fc4SMatthew G. Knepley } 23484961fc4SMatthew G. Knepley /* Consider each face in FIFO */ 23548a46eb9SPierre Jolivet while (fTop < fBottom) PetscCall(DMPlexCheckFace_Internal(dm, faceFIFO, &fTop, &fBottom, cStart, fStart, fEnd, seenCells, flippedCells, seenFaces)); 236e1d83109SMatthew G. Knepley /* Set component for cells and faces */ 237e1d83109SMatthew G. Knepley for (cell = 0; cell < cEnd - cStart; ++cell) { 238e1d83109SMatthew G. Knepley if (PetscBTLookup(seenCells, cell)) cellComp[cell] = comp; 239e1d83109SMatthew G. Knepley } 240e1d83109SMatthew G. Knepley for (face = 0; face < fEnd - fStart; ++face) { 241e1d83109SMatthew G. Knepley if (PetscBTLookup(seenFaces, face)) faceComp[face] = comp; 242e1d83109SMatthew G. Knepley } 243e1d83109SMatthew G. Knepley /* Wipe seenCells and seenFaces for next component */ 2449566063dSJacob Faibussowitsch PetscCall(PetscBTMemzero(fEnd - fStart, seenFaces)); 2459566063dSJacob Faibussowitsch PetscCall(PetscBTMemzero(cEnd - cStart, seenCells)); 246e1d83109SMatthew G. Knepley ++comp; 247e1d83109SMatthew G. Knepley } while (1); 248e1d83109SMatthew G. Knepley numComponents = comp; 24984961fc4SMatthew G. Knepley if (flg) { 25084961fc4SMatthew G. Knepley PetscViewer v; 25184961fc4SMatthew G. Knepley 2529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetStdout(comm, &v)); 2539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(v)); 2549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank)); 2559566063dSJacob Faibussowitsch PetscCall(PetscBTView(cEnd - cStart, flippedCells, v)); 2569566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(v)); 2579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(v)); 25884961fc4SMatthew G. Knepley } 25984961fc4SMatthew G. Knepley /* Now all subdomains are oriented, but we need a consistent parallel orientation */ 26084961fc4SMatthew G. Knepley if (numLeaves >= 0) { 2613d6051bcSMatthew G. Knepley PetscInt maxSupportSize, neighbor; 2623d6051bcSMatthew G. Knepley 263e1d83109SMatthew G. Knepley /* Store orientations of boundary faces*/ 2643d6051bcSMatthew G. Knepley PetscCall(DMPlexGetMaxSizes(dm, NULL, &maxSupportSize)); 2653d6051bcSMatthew G. Knepley PetscCall(PetscCalloc3(numRoots, &rorntComp, numRoots, &lorntComp, maxSupportSize, &locSupport)); 26684961fc4SMatthew G. Knepley for (face = fStart; face < fEnd; ++face) { 267e1d83109SMatthew G. Knepley const PetscInt *cone, *support, *ornt; 2683d6051bcSMatthew G. Knepley PetscInt coneSize, supportSize, Ns = 0, s, l; 269e1d83109SMatthew G. Knepley 2709566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, face, &supportSize)); 2713d6051bcSMatthew G. Knepley /* Ignore overlapping cells */ 2729566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, face, &support)); 2733d6051bcSMatthew G. Knepley for (s = 0; s < supportSize; ++s) { 2743d6051bcSMatthew G. Knepley PetscCall(PetscFindInt(support[s], numLeaves, lpoints, &l)); 2753d6051bcSMatthew G. Knepley if (l >= 0) continue; 2763d6051bcSMatthew G. Knepley locSupport[Ns++] = support[s]; 2773d6051bcSMatthew G. Knepley } 2783d6051bcSMatthew G. Knepley if (Ns != 1) continue; 2793d6051bcSMatthew G. Knepley neighbor = locSupport[0]; 2803d6051bcSMatthew G. Knepley PetscCall(DMPlexGetCone(dm, neighbor, &cone)); 2813d6051bcSMatthew G. Knepley PetscCall(DMPlexGetConeSize(dm, neighbor, &coneSize)); 2823d6051bcSMatthew G. Knepley PetscCall(DMPlexGetConeOrientation(dm, neighbor, &ornt)); 2839371c9d4SSatish Balay for (c = 0; c < coneSize; ++c) 2849371c9d4SSatish Balay if (cone[c] == face) break; 28584961fc4SMatthew G. Knepley if (dim == 1) { 28684961fc4SMatthew G. Knepley /* Use cone position instead, shifted to -1 or 1 */ 2873d6051bcSMatthew G. Knepley if (PetscBTLookup(flippedCells, neighbor - cStart)) rorntComp[face].rank = 1 - c * 2; 2886e1eb25cSMatthew G. Knepley else rorntComp[face].rank = c * 2 - 1; 28984961fc4SMatthew G. Knepley } else { 2903d6051bcSMatthew G. Knepley if (PetscBTLookup(flippedCells, neighbor - cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 : 1; 291e1d83109SMatthew G. Knepley else rorntComp[face].rank = ornt[c] < 0 ? 1 : -1; 29284961fc4SMatthew G. Knepley } 293e1d83109SMatthew G. Knepley rorntComp[face].index = faceComp[face - fStart]; 29484961fc4SMatthew G. Knepley } 295e1d83109SMatthew G. Knepley /* Communicate boundary edge orientations */ 2969566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf, MPIU_2INT, rorntComp, lorntComp, MPI_REPLACE)); 2979566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf, MPIU_2INT, rorntComp, lorntComp, MPI_REPLACE)); 298e1d83109SMatthew G. Knepley } 299e1d83109SMatthew G. Knepley /* Get process adjacency */ 3009566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors)); 301a17aa47eSToby Isaac viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm)); 3029566063dSJacob Faibussowitsch if (flg2) PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 3039566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &selfviewer)); 304e1d83109SMatthew G. Knepley for (comp = 0; comp < numComponents; ++comp) { 305e1d83109SMatthew G. Knepley PetscInt l, n; 30684961fc4SMatthew G. Knepley 307e1d83109SMatthew G. Knepley numNeighbors[comp] = 0; 3089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp])); 309e1d83109SMatthew G. Knepley /* I know this is p^2 time in general, but for bounded degree its alright */ 310e1d83109SMatthew G. Knepley for (l = 0; l < numLeaves; ++l) { 311e1d83109SMatthew G. Knepley const PetscInt face = lpoints[l]; 312e1d83109SMatthew G. Knepley 313e1d83109SMatthew G. Knepley /* Find a representative face (edge) separating pairs of procs */ 3143d6051bcSMatthew G. Knepley if ((face >= fStart) && (face < fEnd) && (faceComp[face - fStart] == comp) && rorntComp[face].rank) { 315269a49d5SMatthew G. Knepley const PetscInt rrank = rpoints[l].rank; 316e1d83109SMatthew G. Knepley const PetscInt rcomp = lorntComp[face].index; 317e1d83109SMatthew G. Knepley 3189371c9d4SSatish Balay for (n = 0; n < numNeighbors[comp]; ++n) 3199371c9d4SSatish Balay if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break; 320e1d83109SMatthew G. Knepley if (n >= numNeighbors[comp]) { 321e1d83109SMatthew G. Knepley PetscInt supportSize; 322e1d83109SMatthew G. Knepley 3239566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, face, &supportSize)); 32463a3b9bcSJacob Faibussowitsch PetscCheck(supportSize == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %" PetscInt_FMT, supportSize); 3259371c9d4SSatish Balay if (flg) 3269371c9d4SSatish 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, 3279371c9d4SSatish Balay rpoints[l].index, rrank, rcomp, lorntComp[face].rank)); 328e1d83109SMatthew G. Knepley neighbors[comp][numNeighbors[comp]++] = l; 329e1d83109SMatthew G. Knepley } 330e1d83109SMatthew G. Knepley } 331e1d83109SMatthew G. Knepley } 332e1d83109SMatthew G. Knepley totNeighbors += numNeighbors[comp]; 333e1d83109SMatthew G. Knepley } 3349566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &selfviewer)); 3359566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 3369566063dSJacob Faibussowitsch if (flg2) PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 3379566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match)); 338e1d83109SMatthew G. Knepley for (comp = 0, off = 0; comp < numComponents; ++comp) { 339e1d83109SMatthew G. Knepley PetscInt n; 340e1d83109SMatthew G. Knepley 341e1d83109SMatthew G. Knepley for (n = 0; n < numNeighbors[comp]; ++n, ++off) { 342e1d83109SMatthew G. Knepley const PetscInt face = lpoints[neighbors[comp][n]]; 343e1d83109SMatthew G. Knepley const PetscInt o = rorntComp[face].rank * lorntComp[face].rank; 344e1d83109SMatthew G. Knepley 345e1d83109SMatthew G. Knepley if (o < 0) match[off] = PETSC_TRUE; 346e1d83109SMatthew G. Knepley else if (o > 0) match[off] = PETSC_FALSE; 34763a3b9bcSJacob 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); 348e1d83109SMatthew G. Knepley nrankComp[off].rank = rpoints[neighbors[comp][n]].rank; 349e1d83109SMatthew G. Knepley nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index; 350e1d83109SMatthew G. Knepley } 3519566063dSJacob Faibussowitsch PetscCall(PetscFree(neighbors[comp])); 35284961fc4SMatthew G. Knepley } 35384961fc4SMatthew G. Knepley /* Collect the graph on 0 */ 354e1d83109SMatthew G. Knepley if (numLeaves >= 0) { 35584961fc4SMatthew G. Knepley Mat G; 35684961fc4SMatthew G. Knepley PetscBT seenProcs, flippedProcs; 35784961fc4SMatthew G. Knepley PetscInt *procFIFO, pTop, pBottom; 3587cadcfe8SMatthew G. Knepley PetscInt *N = NULL, *Noff; 359e1d83109SMatthew G. Knepley PetscSFNode *adj = NULL; 36084961fc4SMatthew G. Knepley PetscBool *val = NULL; 3617cadcfe8SMatthew G. Knepley PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc, p, o; 3629852e123SBarry Smith PetscMPIInt size = 0; 36384961fc4SMatthew G. Knepley 3649566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(numComponents, &flipped)); 3659566063dSJacob Faibussowitsch if (rank == 0) PetscCallMPI(MPI_Comm_size(comm, &size)); 3669566063dSJacob Faibussowitsch PetscCall(PetscCalloc4(size, &recvcounts, size + 1, &displs, size, &Nc, size + 1, &Noff)); 3679566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gather(&numComponents, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm)); 368ad540459SPierre Jolivet for (p = 0; p < size; ++p) displs[p + 1] = displs[p] + Nc[p]; 3699566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscMalloc1(displs[size], &N)); 3709566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gatherv(numNeighbors, numComponents, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm)); 3719852e123SBarry Smith for (p = 0, o = 0; p < size; ++p) { 372e1d83109SMatthew G. Knepley recvcounts[p] = 0; 373e1d83109SMatthew G. Knepley for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o]; 37484961fc4SMatthew G. Knepley displs[p + 1] = displs[p] + recvcounts[p]; 37584961fc4SMatthew G. Knepley } 3769566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscMalloc2(displs[size], &adj, displs[size], &val)); 3779566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gatherv(nrankComp, totNeighbors, MPIU_2INT, adj, recvcounts, displs, MPIU_2INT, 0, comm)); 3789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gatherv(match, totNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm)); 3799566063dSJacob Faibussowitsch PetscCall(PetscFree2(numNeighbors, neighbors)); 380dd400576SPatrick Sanan if (rank == 0) { 381ad540459SPierre Jolivet for (p = 1; p <= size; ++p) Noff[p] = Noff[p - 1] + Nc[p - 1]; 382a83982f3SMatthew G. Knepley if (flg) { 383e1d83109SMatthew G. Knepley PetscInt n; 384e1d83109SMatthew G. Knepley 3859852e123SBarry Smith for (p = 0, off = 0; p < size; ++p) { 386e1d83109SMatthew G. Knepley for (c = 0; c < Nc[p]; ++c) { 38763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %" PetscInt_FMT ":\n", p, c)); 38848a46eb9SPierre 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]])); 38984961fc4SMatthew G. Knepley } 39084961fc4SMatthew G. Knepley } 39184961fc4SMatthew G. Knepley } 39284961fc4SMatthew G. Knepley /* Symmetrize the graph */ 3939566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &G)); 3949566063dSJacob Faibussowitsch PetscCall(MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size])); 3959566063dSJacob Faibussowitsch PetscCall(MatSetUp(G)); 3969852e123SBarry Smith for (p = 0, off = 0; p < size; ++p) { 397e1d83109SMatthew G. Knepley for (c = 0; c < Nc[p]; ++c) { 398e1d83109SMatthew G. Knepley const PetscInt r = Noff[p] + c; 399e1d83109SMatthew G. Knepley PetscInt n; 400e1d83109SMatthew G. Knepley 401e1d83109SMatthew G. Knepley for (n = 0; n < N[r]; ++n, ++off) { 402e1d83109SMatthew G. Knepley const PetscInt q = Noff[adj[off].rank] + adj[off].index; 403e1d83109SMatthew G. Knepley const PetscScalar o = val[off] ? 1.0 : 0.0; 40484961fc4SMatthew G. Knepley 4059566063dSJacob Faibussowitsch PetscCall(MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES)); 4069566063dSJacob Faibussowitsch PetscCall(MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES)); 40784961fc4SMatthew G. Knepley } 40884961fc4SMatthew G. Knepley } 409e1d83109SMatthew G. Knepley } 4109566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY)); 4119566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY)); 41284961fc4SMatthew G. Knepley 4139566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(Noff[size], &seenProcs)); 4149566063dSJacob Faibussowitsch PetscCall(PetscBTMemzero(Noff[size], seenProcs)); 4159566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(Noff[size], &flippedProcs)); 4169566063dSJacob Faibussowitsch PetscCall(PetscBTMemzero(Noff[size], flippedProcs)); 4179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Noff[size], &procFIFO)); 41884961fc4SMatthew G. Knepley pTop = pBottom = 0; 4199852e123SBarry Smith for (p = 0; p < Noff[size]; ++p) { 42084961fc4SMatthew G. Knepley if (PetscBTLookup(seenProcs, p)) continue; 42184961fc4SMatthew G. Knepley /* Initialize FIFO with next proc */ 42284961fc4SMatthew G. Knepley procFIFO[pBottom++] = p; 4239566063dSJacob Faibussowitsch PetscCall(PetscBTSet(seenProcs, p)); 42484961fc4SMatthew G. Knepley /* Consider each proc in FIFO */ 42584961fc4SMatthew G. Knepley while (pTop < pBottom) { 42684961fc4SMatthew G. Knepley const PetscScalar *ornt; 42784961fc4SMatthew G. Knepley const PetscInt *neighbors; 428e1d83109SMatthew G. Knepley PetscInt proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n; 42984961fc4SMatthew G. Knepley 43084961fc4SMatthew G. Knepley proc = procFIFO[pTop++]; 43184961fc4SMatthew G. Knepley flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0; 4329566063dSJacob Faibussowitsch PetscCall(MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt)); 43384961fc4SMatthew G. Knepley /* Loop over neighboring procs */ 43484961fc4SMatthew G. Knepley for (n = 0; n < numNeighbors; ++n) { 43584961fc4SMatthew G. Knepley nproc = neighbors[n]; 43684961fc4SMatthew G. Knepley mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1; 43784961fc4SMatthew G. Knepley seen = PetscBTLookup(seenProcs, nproc); 43884961fc4SMatthew G. Knepley flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0; 43984961fc4SMatthew G. Knepley 44084961fc4SMatthew G. Knepley if (mismatch ^ (flippedA ^ flippedB)) { 44163a3b9bcSJacob 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); 44284961fc4SMatthew G. Knepley if (!flippedB) { 4439566063dSJacob Faibussowitsch PetscCall(PetscBTSet(flippedProcs, nproc)); 44484961fc4SMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable"); 4451dca8a05SBarry Smith } else PetscCheck(!mismatch || !flippedA || !flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable"); 44684961fc4SMatthew G. Knepley if (!seen) { 44784961fc4SMatthew G. Knepley procFIFO[pBottom++] = nproc; 4489566063dSJacob Faibussowitsch PetscCall(PetscBTSet(seenProcs, nproc)); 44984961fc4SMatthew G. Knepley } 45084961fc4SMatthew G. Knepley } 45184961fc4SMatthew G. Knepley } 45284961fc4SMatthew G. Knepley } 4539566063dSJacob Faibussowitsch PetscCall(PetscFree(procFIFO)); 4549566063dSJacob Faibussowitsch PetscCall(MatDestroy(&G)); 4559566063dSJacob Faibussowitsch PetscCall(PetscFree2(adj, val)); 4569566063dSJacob Faibussowitsch PetscCall(PetscBTDestroy(&seenProcs)); 45784961fc4SMatthew G. Knepley } 458e1d83109SMatthew G. Knepley /* Scatter flip flags */ 459e1d83109SMatthew G. Knepley { 460e1d83109SMatthew G. Knepley PetscBool *flips = NULL; 461e1d83109SMatthew G. Knepley 462dd400576SPatrick Sanan if (rank == 0) { 4639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Noff[size], &flips)); 4649852e123SBarry Smith for (p = 0; p < Noff[size]; ++p) { 465e1d83109SMatthew G. Knepley flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE; 4669566063dSJacob Faibussowitsch if (flg && flips[p]) PetscCall(PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p)); 467e1d83109SMatthew G. Knepley } 468ad540459SPierre Jolivet for (p = 0; p < size; ++p) displs[p + 1] = displs[p] + Nc[p]; 469e1d83109SMatthew G. Knepley } 4709566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, numComponents, MPIU_BOOL, 0, comm)); 4719566063dSJacob Faibussowitsch PetscCall(PetscFree(flips)); 47284961fc4SMatthew G. Knepley } 4739566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscBTDestroy(&flippedProcs)); 4749566063dSJacob Faibussowitsch PetscCall(PetscFree(N)); 4759566063dSJacob Faibussowitsch PetscCall(PetscFree4(recvcounts, displs, Nc, Noff)); 4769566063dSJacob Faibussowitsch PetscCall(PetscFree2(nrankComp, match)); 477e1d83109SMatthew G. Knepley 478e1d83109SMatthew G. Knepley /* Decide whether to flip cells in each component */ 4799371c9d4SSatish Balay for (c = 0; c < cEnd - cStart; ++c) { 4809371c9d4SSatish Balay if (flipped[cellComp[c]]) PetscCall(PetscBTNegate(flippedCells, c)); 4819371c9d4SSatish Balay } 4829566063dSJacob Faibussowitsch PetscCall(PetscFree(flipped)); 48384961fc4SMatthew G. Knepley } 48484961fc4SMatthew G. Knepley if (flg) { 48584961fc4SMatthew G. Knepley PetscViewer v; 48684961fc4SMatthew G. Knepley 4879566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetStdout(comm, &v)); 4889566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(v)); 4899566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank)); 4909566063dSJacob Faibussowitsch PetscCall(PetscBTView(cEnd - cStart, flippedCells, v)); 4919566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(v)); 4929566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(v)); 49384961fc4SMatthew G. Knepley } 49484961fc4SMatthew G. Knepley /* Reverse flipped cells in the mesh */ 49584961fc4SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 49648a46eb9SPierre Jolivet if (PetscBTLookup(flippedCells, c - cStart)) PetscCall(DMPlexOrientPoint(dm, c, -1)); 49784961fc4SMatthew G. Knepley } 4989566063dSJacob Faibussowitsch PetscCall(PetscBTDestroy(&seenCells)); 4999566063dSJacob Faibussowitsch PetscCall(PetscBTDestroy(&flippedCells)); 5009566063dSJacob Faibussowitsch PetscCall(PetscBTDestroy(&seenFaces)); 5019566063dSJacob Faibussowitsch PetscCall(PetscFree2(numNeighbors, neighbors)); 5023d6051bcSMatthew G. Knepley PetscCall(PetscFree3(rorntComp, lorntComp, locSupport)); 5039566063dSJacob Faibussowitsch PetscCall(PetscFree3(faceFIFO, cellComp, faceComp)); 50484961fc4SMatthew G. Knepley PetscFunctionReturn(0); 50584961fc4SMatthew G. Knepley } 506