xref: /petsc/src/dm/impls/plex/plexorient.c (revision 8112c1cbf372cb53bf7c5aca994a84a6a303db4d)
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:
10a1cb98faSBarry 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 
161cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexOrient()`, `DMPlexGetCone()`, `DMPlexGetConeOrientation()`, `DMPlexInterpolate()`, `DMPlexGetChart()`
1727d0e99aSVaclav Hapla @*/
DMPlexOrientPoint(DM dm,PetscInt p,PetscInt o)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));
2885036b15SMatthew G. Knepley   arr = DMPolytopeTypeGetArrangement(ct, o);
290851c46eSMatthew G. Knepley   if (!arr) PetscFunctionReturn(PETSC_SUCCESS);
309566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
319566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCone(dm, p, &cone));
329566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeOrientation(dm, p, &ornt));
339566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, coneSize, MPIU_INT, &newcone));
349566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, coneSize, MPIU_INT, &newornt));
35b5a892a1SMatthew G. Knepley   for (c = 0; c < coneSize; ++c) {
36b5a892a1SMatthew G. Knepley     DMPolytopeType ft;
37b5a892a1SMatthew G. Knepley     PetscInt       nO;
3827d0e99aSVaclav Hapla 
399566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCellType(dm, cone[c], &ft));
4085036b15SMatthew G. Knepley     nO         = DMPolytopeTypeGetNumArrangements(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]]);
4363a3b9bcSJacob 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]);
4427d0e99aSVaclav Hapla   }
459566063dSJacob Faibussowitsch   PetscCall(DMPlexSetCone(dm, p, newcone));
469566063dSJacob Faibussowitsch   PetscCall(DMPlexSetConeOrientation(dm, p, newornt));
479566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm, coneSize, MPIU_INT, &newcone));
489566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm, coneSize, MPIU_INT, &newornt));
49b5a892a1SMatthew G. Knepley   /* Update orientation of this point in the support points */
509566063dSJacob Faibussowitsch   PetscCall(DMPlexGetSupportSize(dm, p, &supportSize));
519566063dSJacob Faibussowitsch   PetscCall(DMPlexGetSupport(dm, p, &support));
52b5a892a1SMatthew G. Knepley   for (s = 0; s < supportSize; ++s) {
539566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeSize(dm, support[s], &coneSize));
549566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCone(dm, support[s], &cone));
559566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeOrientation(dm, support[s], &ornt));
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);
629566063dSJacob Faibussowitsch       PetscCall(DMPlexInsertConeOrientation(dm, support[s], c, po));
6384961fc4SMatthew G. Knepley     }
6484961fc4SMatthew G. Knepley   }
653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6684961fc4SMatthew G. Knepley }
6784961fc4SMatthew G. Knepley 
GetPointIndex(PetscInt point,PetscInt pStart,PetscInt pEnd,const PetscInt points[])680851c46eSMatthew G. Knepley static PetscInt GetPointIndex(PetscInt point, PetscInt pStart, PetscInt pEnd, const PetscInt points[])
690851c46eSMatthew G. Knepley {
700851c46eSMatthew G. Knepley   if (points) {
710851c46eSMatthew G. Knepley     PetscInt loc;
720851c46eSMatthew G. Knepley 
733387f462SBarry Smith     PetscCallAbort(PETSC_COMM_SELF, PetscFindInt(point, pEnd - pStart, points, &loc));
740851c46eSMatthew G. Knepley     if (loc >= 0) return loc;
750851c46eSMatthew G. Knepley   } else {
760851c46eSMatthew G. Knepley     if (point >= pStart && point < pEnd) return point - pStart;
770851c46eSMatthew G. Knepley   }
780851c46eSMatthew G. Knepley   return -1;
790851c46eSMatthew G. Knepley }
800851c46eSMatthew G. Knepley 
817b310ce2SMatthew G. Knepley /*
827b310ce2SMatthew G. Knepley   - Checks face match
837b310ce2SMatthew G. Knepley     - Flips non-matching
847b310ce2SMatthew G. Knepley   - Inserts faces of support cells in FIFO
857b310ce2SMatthew G. Knepley */
DMPlexCheckFace_Internal(DM dm,PetscInt * faceFIFO,PetscInt * fTop,PetscInt * fBottom,IS cellIS,IS faceIS,PetscBT seenCells,PetscBT flippedCells,PetscBT seenFaces)860851c46eSMatthew G. Knepley static PetscErrorCode DMPlexCheckFace_Internal(DM dm, PetscInt *faceFIFO, PetscInt *fTop, PetscInt *fBottom, IS cellIS, IS faceIS, PetscBT seenCells, PetscBT flippedCells, PetscBT seenFaces)
870851c46eSMatthew G. Knepley {
880851c46eSMatthew G. Knepley   const PetscInt *supp, *coneA, *coneB, *coneOA, *coneOB;
890851c46eSMatthew G. Knepley   PetscInt        suppSize, Ns = 0, coneSizeA, coneSizeB, posA = -1, posB = -1;
900851c46eSMatthew G. Knepley   PetscInt        face, dim, indC[3], indS[3], seenA, flippedA, seenB, flippedB, mismatch;
910851c46eSMatthew G. Knepley   const PetscInt *cells, *faces;
920851c46eSMatthew G. Knepley   PetscInt        cStart, cEnd, fStart, fEnd;
930851c46eSMatthew G. Knepley 
940851c46eSMatthew G. Knepley   PetscFunctionBegin;
950851c46eSMatthew G. Knepley   face = faceFIFO[(*fTop)++];
960851c46eSMatthew G. Knepley   PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
970851c46eSMatthew G. Knepley   PetscCall(ISGetPointRange(faceIS, &fStart, &fEnd, &faces));
980851c46eSMatthew G. Knepley   PetscCall(DMPlexGetPointDepth(dm, cells ? cells[cStart] : cStart, &dim));
990851c46eSMatthew G. Knepley   PetscCall(DMPlexGetSupportSize(dm, face, &suppSize));
1000851c46eSMatthew G. Knepley   PetscCall(DMPlexGetSupport(dm, face, &supp));
1010851c46eSMatthew G. Knepley   // Filter the support
1020851c46eSMatthew G. Knepley   for (PetscInt s = 0; s < suppSize; ++s) {
1030851c46eSMatthew G. Knepley     // Filter support
1040851c46eSMatthew G. Knepley     indC[Ns] = GetPointIndex(supp[s], cStart, cEnd, cells);
1050851c46eSMatthew G. Knepley     indS[Ns] = s;
1060851c46eSMatthew G. Knepley     if (indC[Ns] >= 0) ++Ns;
1070851c46eSMatthew G. Knepley   }
1080851c46eSMatthew G. Knepley   if (Ns < 2) PetscFunctionReturn(PETSC_SUCCESS);
1090851c46eSMatthew G. Knepley   PetscCheck(Ns == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %" PetscInt_FMT, Ns);
1100851c46eSMatthew G. Knepley   PetscCheck(indC[0] >= 0 && indC[1] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Support cells %" PetscInt_FMT " (%" PetscInt_FMT ") and %" PetscInt_FMT " (%" PetscInt_FMT ") are not both valid", supp[0], indC[0], supp[1], indC[1]);
1110851c46eSMatthew G. Knepley   seenA    = PetscBTLookup(seenCells, indC[0]);
1120851c46eSMatthew G. Knepley   flippedA = PetscBTLookup(flippedCells, indC[0]) ? 1 : 0;
1130851c46eSMatthew G. Knepley   seenB    = PetscBTLookup(seenCells, indC[1]);
1140851c46eSMatthew G. Knepley   flippedB = PetscBTLookup(flippedCells, indC[1]) ? 1 : 0;
1150851c46eSMatthew G. Knepley 
1160851c46eSMatthew G. Knepley   PetscCall(DMPlexGetConeSize(dm, supp[indS[0]], &coneSizeA));
1170851c46eSMatthew G. Knepley   PetscCall(DMPlexGetConeSize(dm, supp[indS[1]], &coneSizeB));
1180851c46eSMatthew G. Knepley   PetscCall(DMPlexGetCone(dm, supp[indS[0]], &coneA));
1190851c46eSMatthew G. Knepley   PetscCall(DMPlexGetCone(dm, supp[indS[1]], &coneB));
1200851c46eSMatthew G. Knepley   PetscCall(DMPlexGetConeOrientation(dm, supp[indS[0]], &coneOA));
1210851c46eSMatthew G. Knepley   PetscCall(DMPlexGetConeOrientation(dm, supp[indS[1]], &coneOB));
1220851c46eSMatthew G. Knepley   for (PetscInt c = 0; c < coneSizeA; ++c) {
1230851c46eSMatthew G. Knepley     const PetscInt indF = GetPointIndex(coneA[c], fStart, fEnd, faces);
1240851c46eSMatthew G. Knepley 
1250851c46eSMatthew G. Knepley     // Filter cone
1260851c46eSMatthew G. Knepley     if (indF < 0) continue;
1270851c46eSMatthew G. Knepley     if (!PetscBTLookup(seenFaces, indF)) {
1280851c46eSMatthew G. Knepley       faceFIFO[(*fBottom)++] = coneA[c];
1290851c46eSMatthew G. Knepley       PetscCall(PetscBTSet(seenFaces, indF));
1300851c46eSMatthew G. Knepley     }
1310851c46eSMatthew G. Knepley     if (coneA[c] == face) posA = c;
1320851c46eSMatthew G. Knepley     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);
1330851c46eSMatthew G. Knepley   }
1340851c46eSMatthew G. Knepley   PetscCheck(posA >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " could not be located in cell %" PetscInt_FMT, face, supp[indS[0]]);
1350851c46eSMatthew G. Knepley   for (PetscInt c = 0; c < coneSizeB; ++c) {
1360851c46eSMatthew G. Knepley     const PetscInt indF = GetPointIndex(coneB[c], fStart, fEnd, faces);
1370851c46eSMatthew G. Knepley 
1380851c46eSMatthew G. Knepley     // Filter cone
1390851c46eSMatthew G. Knepley     if (indF < 0) continue;
1400851c46eSMatthew G. Knepley     if (!PetscBTLookup(seenFaces, indF)) {
1410851c46eSMatthew G. Knepley       faceFIFO[(*fBottom)++] = coneB[c];
1420851c46eSMatthew G. Knepley       PetscCall(PetscBTSet(seenFaces, indF));
1430851c46eSMatthew G. Knepley     }
1440851c46eSMatthew G. Knepley     if (coneB[c] == face) posB = c;
1450851c46eSMatthew G. Knepley     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);
1460851c46eSMatthew G. Knepley   }
1470851c46eSMatthew G. Knepley   PetscCheck(posB >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " could not be located in cell %" PetscInt_FMT, face, supp[indS[1]]);
1480851c46eSMatthew G. Knepley 
1490851c46eSMatthew G. Knepley   if (dim == 1) {
1500851c46eSMatthew G. Knepley     mismatch = posA == posB;
1510851c46eSMatthew G. Knepley   } else {
1520851c46eSMatthew G. Knepley     mismatch = coneOA[posA] == coneOB[posB];
1530851c46eSMatthew G. Knepley   }
1540851c46eSMatthew G. Knepley 
1550851c46eSMatthew G. Knepley   if (mismatch ^ (flippedA ^ flippedB)) {
1560851c46eSMatthew G. Knepley     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", supp[indS[0]], supp[indS[1]]);
1570851c46eSMatthew G. Knepley     if (!seenA && !flippedA) {
1580851c46eSMatthew G. Knepley       PetscCall(PetscBTSet(flippedCells, indC[0]));
1590851c46eSMatthew G. Knepley     } else if (!seenB && !flippedB) {
1600851c46eSMatthew G. Knepley       PetscCall(PetscBTSet(flippedCells, indC[1]));
1610851c46eSMatthew G. Knepley     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
1620851c46eSMatthew G. Knepley   } else PetscCheck(!mismatch || !flippedA || !flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
1630851c46eSMatthew G. Knepley   PetscCall(PetscBTSet(seenCells, indC[0]));
1640851c46eSMatthew G. Knepley   PetscCall(PetscBTSet(seenCells, indC[1]));
1650851c46eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1660851c46eSMatthew G. Knepley }
1670851c46eSMatthew G. Knepley 
DMPlexCheckFace_Old_Internal(DM dm,PetscInt * faceFIFO,PetscInt * fTop,PetscInt * fBottom,PetscInt cStart,PetscInt fStart,PetscInt fEnd,PetscBT seenCells,PetscBT flippedCells,PetscBT seenFaces)1680851c46eSMatthew G. Knepley static PetscErrorCode DMPlexCheckFace_Old_Internal(DM dm, PetscInt *faceFIFO, PetscInt *fTop, PetscInt *fBottom, PetscInt cStart, PetscInt fStart, PetscInt fEnd, PetscBT seenCells, PetscBT flippedCells, PetscBT seenFaces)
169d71ae5a4SJacob Faibussowitsch {
1707b310ce2SMatthew G. Knepley   const PetscInt *support, *coneA, *coneB, *coneOA, *coneOB;
1717b310ce2SMatthew G. Knepley   PetscInt        supportSize, coneSizeA, coneSizeB, posA = -1, posB = -1;
1727b2de0fdSMatthew G. Knepley   PetscInt        face, dim, seenA, flippedA, seenB, flippedB, mismatch, c;
1737b310ce2SMatthew G. Knepley 
1747b310ce2SMatthew G. Knepley   PetscFunctionBegin;
1757b310ce2SMatthew G. Knepley   face = faceFIFO[(*fTop)++];
1769566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1779566063dSJacob Faibussowitsch   PetscCall(DMPlexGetSupportSize(dm, face, &supportSize));
1789566063dSJacob Faibussowitsch   PetscCall(DMPlexGetSupport(dm, face, &support));
1793ba16761SJacob Faibussowitsch   if (supportSize < 2) PetscFunctionReturn(PETSC_SUCCESS);
18063a3b9bcSJacob Faibussowitsch   PetscCheck(supportSize == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %" PetscInt_FMT, supportSize);
1817b310ce2SMatthew G. Knepley   seenA    = PetscBTLookup(seenCells, support[0] - cStart);
1827b310ce2SMatthew G. Knepley   flippedA = PetscBTLookup(flippedCells, support[0] - cStart) ? 1 : 0;
1837b310ce2SMatthew G. Knepley   seenB    = PetscBTLookup(seenCells, support[1] - cStart);
1847b310ce2SMatthew G. Knepley   flippedB = PetscBTLookup(flippedCells, support[1] - cStart) ? 1 : 0;
1857b310ce2SMatthew G. Knepley 
1869566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeSize(dm, support[0], &coneSizeA));
1879566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeSize(dm, support[1], &coneSizeB));
1889566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCone(dm, support[0], &coneA));
1899566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCone(dm, support[1], &coneB));
1909566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeOrientation(dm, support[0], &coneOA));
1919566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeOrientation(dm, support[1], &coneOB));
1927b310ce2SMatthew G. Knepley   for (c = 0; c < coneSizeA; ++c) {
1937b310ce2SMatthew G. Knepley     if (!PetscBTLookup(seenFaces, coneA[c] - fStart)) {
1947b310ce2SMatthew G. Knepley       faceFIFO[(*fBottom)++] = coneA[c];
1959566063dSJacob Faibussowitsch       PetscCall(PetscBTSet(seenFaces, coneA[c] - fStart));
1967b310ce2SMatthew G. Knepley     }
1977b310ce2SMatthew G. Knepley     if (coneA[c] == face) posA = c;
19863a3b9bcSJacob 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);
1997b310ce2SMatthew G. Knepley   }
20063a3b9bcSJacob 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]);
2017b310ce2SMatthew G. Knepley   for (c = 0; c < coneSizeB; ++c) {
2027b310ce2SMatthew G. Knepley     if (!PetscBTLookup(seenFaces, coneB[c] - fStart)) {
2037b310ce2SMatthew G. Knepley       faceFIFO[(*fBottom)++] = coneB[c];
2049566063dSJacob Faibussowitsch       PetscCall(PetscBTSet(seenFaces, coneB[c] - fStart));
2057b310ce2SMatthew G. Knepley     }
2067b310ce2SMatthew G. Knepley     if (coneB[c] == face) posB = c;
20763a3b9bcSJacob 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);
2087b310ce2SMatthew G. Knepley   }
20963a3b9bcSJacob 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]);
2107b310ce2SMatthew G. Knepley 
2117b310ce2SMatthew G. Knepley   if (dim == 1) {
2127b310ce2SMatthew G. Knepley     mismatch = posA == posB;
2137b310ce2SMatthew G. Knepley   } else {
2147b310ce2SMatthew G. Knepley     mismatch = coneOA[posA] == coneOB[posB];
2157b310ce2SMatthew G. Knepley   }
2167b310ce2SMatthew G. Knepley 
2177b310ce2SMatthew G. Knepley   if (mismatch ^ (flippedA ^ flippedB)) {
21863a3b9bcSJacob 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]);
2197b310ce2SMatthew G. Knepley     if (!seenA && !flippedA) {
2209566063dSJacob Faibussowitsch       PetscCall(PetscBTSet(flippedCells, support[0] - cStart));
2217b310ce2SMatthew G. Knepley     } else if (!seenB && !flippedB) {
2229566063dSJacob Faibussowitsch       PetscCall(PetscBTSet(flippedCells, support[1] - cStart));
2237b310ce2SMatthew G. Knepley     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
2241dca8a05SBarry Smith   } else PetscCheck(!mismatch || !flippedA || !flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
2259566063dSJacob Faibussowitsch   PetscCall(PetscBTSet(seenCells, support[0] - cStart));
2269566063dSJacob Faibussowitsch   PetscCall(PetscBTSet(seenCells, support[1] - cStart));
2273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2287b310ce2SMatthew G. Knepley }
2297b310ce2SMatthew G. Knepley 
2300851c46eSMatthew G. Knepley /*
2310851c46eSMatthew G. Knepley   DMPlexOrient_Serial - Compute valid orientation for local connected components
2320851c46eSMatthew G. Knepley 
2330851c46eSMatthew G. Knepley   Not collective
2340851c46eSMatthew G. Knepley 
2350851c46eSMatthew G. Knepley   Input Parameters:
2360851c46eSMatthew G. Knepley   + dm         - The `DM`
2370851c46eSMatthew G. Knepley   - cellHeight - The height of k-cells to be oriented
2380851c46eSMatthew G. Knepley 
2390851c46eSMatthew G. Knepley   Output Parameters:
2400851c46eSMatthew G. Knepley   + Ncomp        - The number of connected component
2410851c46eSMatthew G. Knepley   . cellComp     - The connected component for each local cell
2420851c46eSMatthew G. Knepley   . faceComp     - The connected component for each local face
2430851c46eSMatthew G. Knepley   - flippedCells - Marked cells should be inverted
2440851c46eSMatthew G. Knepley 
2450851c46eSMatthew G. Knepley   Level: developer
2460851c46eSMatthew G. Knepley 
2475ef53592SPierre Jolivet .seealso: `DMPlexOrient()`
2480851c46eSMatthew G. Knepley */
DMPlexOrient_Serial(DM dm,IS cellIS,IS faceIS,PetscInt * Ncomp,PetscInt cellComp[],PetscInt faceComp[],PetscBT flippedCells)2490851c46eSMatthew G. Knepley static PetscErrorCode DMPlexOrient_Serial(DM dm, IS cellIS, IS faceIS, PetscInt *Ncomp, PetscInt cellComp[], PetscInt faceComp[], PetscBT flippedCells)
2500851c46eSMatthew G. Knepley {
2510851c46eSMatthew G. Knepley   PetscBT         seenCells, seenFaces;
2520851c46eSMatthew G. Knepley   PetscInt       *faceFIFO;
2530851c46eSMatthew G. Knepley   const PetscInt *cells = NULL, *faces = NULL;
2540851c46eSMatthew G. Knepley   PetscInt        cStart = 0, cEnd = 0, fStart = 0, fEnd = 0;
2550851c46eSMatthew G. Knepley 
2560851c46eSMatthew G. Knepley   PetscFunctionBegin;
2570851c46eSMatthew G. Knepley   if (cellIS) PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
2580851c46eSMatthew G. Knepley   if (faceIS) PetscCall(ISGetPointRange(faceIS, &fStart, &fEnd, &faces));
2590851c46eSMatthew G. Knepley   PetscCall(PetscBTCreate(cEnd - cStart, &seenCells));
2600851c46eSMatthew G. Knepley   PetscCall(PetscBTMemzero(cEnd - cStart, seenCells));
2610851c46eSMatthew G. Knepley   PetscCall(PetscBTCreate(fEnd - fStart, &seenFaces));
2620851c46eSMatthew G. Knepley   PetscCall(PetscBTMemzero(fEnd - fStart, seenFaces));
2630851c46eSMatthew G. Knepley   PetscCall(PetscMalloc1(fEnd - fStart, &faceFIFO));
2640851c46eSMatthew G. Knepley   *Ncomp = 0;
2650851c46eSMatthew G. Knepley   for (PetscInt c = 0; c < cEnd - cStart; ++c) cellComp[c] = -1;
2660851c46eSMatthew G. Knepley   do {
2670851c46eSMatthew G. Knepley     PetscInt cc, fTop, fBottom;
2680851c46eSMatthew G. Knepley 
2690851c46eSMatthew G. Knepley     // Look for first unmarked cell
2700851c46eSMatthew G. Knepley     for (cc = cStart; cc < cEnd; ++cc)
2710851c46eSMatthew G. Knepley       if (cellComp[cc - cStart] < 0) break;
2720851c46eSMatthew G. Knepley     if (cc >= cEnd) break;
2730851c46eSMatthew G. Knepley     // Initialize FIFO with first cell in component
2740851c46eSMatthew G. Knepley     {
2750851c46eSMatthew G. Knepley       const PetscInt  cell = cells ? cells[cc] : cc;
2760851c46eSMatthew G. Knepley       const PetscInt *cone;
2770851c46eSMatthew G. Knepley       PetscInt        coneSize;
2780851c46eSMatthew G. Knepley 
2790851c46eSMatthew G. Knepley       fTop = fBottom = 0;
2800851c46eSMatthew G. Knepley       PetscCall(DMPlexGetConeSize(dm, cell, &coneSize));
2810851c46eSMatthew G. Knepley       PetscCall(DMPlexGetCone(dm, cell, &cone));
2820851c46eSMatthew G. Knepley       for (PetscInt c = 0; c < coneSize; ++c) {
2836de65179SMatthew G. Knepley         const PetscInt idx = GetPointIndex(cone[c], fStart, fEnd, faces);
2846de65179SMatthew G. Knepley 
2850851c46eSMatthew G. Knepley         // Cell faces are guaranteed to be in the face set
2866de65179SMatthew G. Knepley         PetscCheck(idx >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " of cell %" PetscInt_FMT " is not present in the label", cone[c], cell);
2870851c46eSMatthew G. Knepley         faceFIFO[fBottom++] = cone[c];
2886de65179SMatthew G. Knepley         PetscCall(PetscBTSet(seenFaces, idx));
2890851c46eSMatthew G. Knepley       }
2900851c46eSMatthew G. Knepley       PetscCall(PetscBTSet(seenCells, cc - cStart));
2910851c46eSMatthew G. Knepley     }
2920851c46eSMatthew G. Knepley     // Consider each face in FIFO
2930851c46eSMatthew G. Knepley     while (fTop < fBottom) PetscCall(DMPlexCheckFace_Internal(dm, faceFIFO, &fTop, &fBottom, cellIS, faceIS, seenCells, flippedCells, seenFaces));
2940851c46eSMatthew G. Knepley     // Set component for cells and faces
2950851c46eSMatthew G. Knepley     for (PetscInt c = 0; c < cEnd - cStart; ++c) {
2960851c46eSMatthew G. Knepley       if (PetscBTLookup(seenCells, c)) cellComp[c] = *Ncomp;
2970851c46eSMatthew G. Knepley     }
2980851c46eSMatthew G. Knepley     for (PetscInt f = 0; f < fEnd - fStart; ++f) {
2990851c46eSMatthew G. Knepley       if (PetscBTLookup(seenFaces, f)) faceComp[f] = *Ncomp;
3000851c46eSMatthew G. Knepley     }
3010851c46eSMatthew G. Knepley     // Wipe seenCells and seenFaces for next component
3020851c46eSMatthew G. Knepley     PetscCall(PetscBTMemzero(fEnd - fStart, seenFaces));
3030851c46eSMatthew G. Knepley     PetscCall(PetscBTMemzero(cEnd - cStart, seenCells));
3040851c46eSMatthew G. Knepley     ++(*Ncomp);
3050851c46eSMatthew G. Knepley   } while (1);
3060851c46eSMatthew G. Knepley   PetscCall(PetscBTDestroy(&seenCells));
3070851c46eSMatthew G. Knepley   PetscCall(PetscBTDestroy(&seenFaces));
3080851c46eSMatthew G. Knepley   PetscCall(PetscFree(faceFIFO));
3090851c46eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
3100851c46eSMatthew G. Knepley }
3110851c46eSMatthew G. Knepley 
31284961fc4SMatthew G. Knepley /*@
31384961fc4SMatthew G. Knepley   DMPlexOrient - Give a consistent orientation to the input mesh
31484961fc4SMatthew G. Knepley 
3152fe279fdSBarry Smith   Input Parameter:
316a1cb98faSBarry Smith . dm - The `DM`
31784961fc4SMatthew G. Knepley 
318a1cb98faSBarry Smith   Note:
319a1cb98faSBarry Smith   The orientation data for the `DM` are change in-place.
320a1cb98faSBarry Smith 
321a1cb98faSBarry Smith   This routine will fail for non-orientable surfaces, such as the Moebius strip.
32284961fc4SMatthew G. Knepley 
32384961fc4SMatthew G. Knepley   Level: advanced
32484961fc4SMatthew G. Knepley 
32560225df5SJacob Faibussowitsch .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMCreate()`
32684961fc4SMatthew G. Knepley @*/
DMPlexOrient(DM dm)327d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexOrient(DM dm)
328d71ae5a4SJacob Faibussowitsch {
3290851c46eSMatthew G. Knepley #if 0
3300851c46eSMatthew G. Knepley   IS cellIS, faceIS;
3310851c46eSMatthew G. Knepley 
3320851c46eSMatthew G. Knepley   PetscFunctionBegin;
3330851c46eSMatthew G. Knepley   PetscCall(DMPlexGetAllCells_Internal(dm, &cellIS));
3340851c46eSMatthew G. Knepley   PetscCall(DMPlexGetAllFaces_Internal(dm, &faceIS));
3350851c46eSMatthew G. Knepley   PetscCall(DMPlexOrientCells_Internal(dm, cellIS, faceIS));
3360851c46eSMatthew G. Knepley   PetscCall(ISDestroy(&cellIS));
3370851c46eSMatthew G. Knepley   PetscCall(ISDestroy(&faceIS));
3380851c46eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
3390851c46eSMatthew G. Knepley #else
34084961fc4SMatthew G. Knepley   MPI_Comm           comm;
341e1d83109SMatthew G. Knepley   PetscSF            sf;
342e1d83109SMatthew G. Knepley   const PetscInt    *lpoints;
343e1d83109SMatthew G. Knepley   const PetscSFNode *rpoints;
344e1d83109SMatthew G. Knepley   PetscSFNode       *rorntComp = NULL, *lorntComp = NULL;
3453d6051bcSMatthew G. Knepley   PetscInt          *numNeighbors, **neighbors, *locSupport = NULL;
346e1d83109SMatthew G. Knepley   PetscSFNode       *nrankComp;
347e1d83109SMatthew G. Knepley   PetscBool         *match, *flipped;
34884961fc4SMatthew G. Knepley   PetscBT            seenCells, flippedCells, seenFaces;
349e1d83109SMatthew G. Knepley   PetscInt          *faceFIFO, fTop, fBottom, *cellComp, *faceComp;
3507cadcfe8SMatthew G. Knepley   PetscInt           numLeaves, numRoots, dim, h, cStart, cEnd, c, cell, fStart, fEnd, face, off, totNeighbors = 0;
35139526728SToby Isaac   PetscMPIInt        rank, size, numComponents, comp = 0;
35239526728SToby Isaac   PetscBool          flg, flg2;
353a17aa47eSToby Isaac   PetscViewer        viewer = NULL, selfviewer = NULL;
35484961fc4SMatthew G. Knepley 
35584961fc4SMatthew G. Knepley   PetscFunctionBegin;
3569566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
3579566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
3589566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
3599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-orientation_view", &flg));
3609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-orientation_view_synchronized", &flg2));
3619566063dSJacob Faibussowitsch   PetscCall(DMGetPointSF(dm, &sf));
3629566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints));
36384961fc4SMatthew G. Knepley   /* Truth Table
36484961fc4SMatthew G. Knepley      mismatch    flips   do action   mismatch   flipA ^ flipB   action
36584961fc4SMatthew G. Knepley          F       0 flips     no         F             F           F
36684961fc4SMatthew G. Knepley          F       1 flip      yes        F             T           T
36784961fc4SMatthew G. Knepley          F       2 flips     no         T             F           T
36884961fc4SMatthew G. Knepley          T       0 flips     yes        T             T           F
36984961fc4SMatthew G. Knepley          T       1 flip      no
37084961fc4SMatthew G. Knepley          T       2 flips     yes
37184961fc4SMatthew G. Knepley   */
3729566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
3739566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(dm, &h));
3749566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, h, &cStart, &cEnd));
3759566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, h + 1, &fStart, &fEnd));
3769566063dSJacob Faibussowitsch   PetscCall(PetscBTCreate(cEnd - cStart, &seenCells));
3779566063dSJacob Faibussowitsch   PetscCall(PetscBTMemzero(cEnd - cStart, seenCells));
3789566063dSJacob Faibussowitsch   PetscCall(PetscBTCreate(cEnd - cStart, &flippedCells));
3799566063dSJacob Faibussowitsch   PetscCall(PetscBTMemzero(cEnd - cStart, flippedCells));
3809566063dSJacob Faibussowitsch   PetscCall(PetscBTCreate(fEnd - fStart, &seenFaces));
3819566063dSJacob Faibussowitsch   PetscCall(PetscBTMemzero(fEnd - fStart, seenFaces));
3829566063dSJacob Faibussowitsch   PetscCall(PetscCalloc3(fEnd - fStart, &faceFIFO, cEnd - cStart, &cellComp, fEnd - fStart, &faceComp));
383e1d83109SMatthew G. Knepley   /*
384e1d83109SMatthew G. Knepley    OLD STYLE
385e1d83109SMatthew G. Knepley    - Add an integer array over cells and faces (component) for connected component number
386e1d83109SMatthew G. Knepley    Foreach component
387e1d83109SMatthew G. Knepley      - Mark the initial cell as seen
388e1d83109SMatthew G. Knepley      - Process component as usual
389e1d83109SMatthew G. Knepley      - Set component for all seenCells
390e1d83109SMatthew G. Knepley      - Wipe seenCells and seenFaces (flippedCells can stay)
391e1d83109SMatthew G. Knepley    - Generate parallel adjacency for component using SF and seenFaces
392e1d83109SMatthew G. Knepley    - Collect numComponents adj data from each proc to 0
393e1d83109SMatthew G. Knepley    - Build same serial graph
394e1d83109SMatthew G. Knepley    - Use same solver
39515229ffcSPierre Jolivet    - Use Scatterv to send back flipped flags for each component
396e1d83109SMatthew G. Knepley    - Negate flippedCells by component
397e1d83109SMatthew G. Knepley 
398e1d83109SMatthew G. Knepley    NEW STYLE
399e1d83109SMatthew G. Knepley    - Create the adj on each process
400e1d83109SMatthew G. Knepley    - Bootstrap to complete graph on proc 0
401e1d83109SMatthew G. Knepley   */
402e1d83109SMatthew G. Knepley   /* Loop over components */
403e1d83109SMatthew G. Knepley   for (cell = cStart; cell < cEnd; ++cell) cellComp[cell - cStart] = -1;
404e1d83109SMatthew G. Knepley   do {
405e1d83109SMatthew G. Knepley     /* Look for first unmarked cell */
4069371c9d4SSatish Balay     for (cell = cStart; cell < cEnd; ++cell)
4079371c9d4SSatish Balay       if (cellComp[cell - cStart] < 0) break;
408e1d83109SMatthew G. Knepley     if (cell >= cEnd) break;
409e1d83109SMatthew G. Knepley     /* Initialize FIFO with first cell in component */
410e1d83109SMatthew G. Knepley     {
41184961fc4SMatthew G. Knepley       const PetscInt *cone;
41284961fc4SMatthew G. Knepley       PetscInt        coneSize;
41384961fc4SMatthew G. Knepley 
414e1d83109SMatthew G. Knepley       fTop = fBottom = 0;
4159566063dSJacob Faibussowitsch       PetscCall(DMPlexGetConeSize(dm, cell, &coneSize));
4169566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(dm, cell, &cone));
41784961fc4SMatthew G. Knepley       for (c = 0; c < coneSize; ++c) {
41884961fc4SMatthew G. Knepley         faceFIFO[fBottom++] = cone[c];
4199566063dSJacob Faibussowitsch         PetscCall(PetscBTSet(seenFaces, cone[c] - fStart));
42084961fc4SMatthew G. Knepley       }
4219566063dSJacob Faibussowitsch       PetscCall(PetscBTSet(seenCells, cell - cStart));
42284961fc4SMatthew G. Knepley     }
42384961fc4SMatthew G. Knepley     /* Consider each face in FIFO */
4240851c46eSMatthew G. Knepley     while (fTop < fBottom) PetscCall(DMPlexCheckFace_Old_Internal(dm, faceFIFO, &fTop, &fBottom, cStart, fStart, fEnd, seenCells, flippedCells, seenFaces));
425e1d83109SMatthew G. Knepley     /* Set component for cells and faces */
426e1d83109SMatthew G. Knepley     for (cell = 0; cell < cEnd - cStart; ++cell) {
427e1d83109SMatthew G. Knepley       if (PetscBTLookup(seenCells, cell)) cellComp[cell] = comp;
428e1d83109SMatthew G. Knepley     }
429e1d83109SMatthew G. Knepley     for (face = 0; face < fEnd - fStart; ++face) {
430e1d83109SMatthew G. Knepley       if (PetscBTLookup(seenFaces, face)) faceComp[face] = comp;
431e1d83109SMatthew G. Knepley     }
432e1d83109SMatthew G. Knepley     /* Wipe seenCells and seenFaces for next component */
4339566063dSJacob Faibussowitsch     PetscCall(PetscBTMemzero(fEnd - fStart, seenFaces));
4349566063dSJacob Faibussowitsch     PetscCall(PetscBTMemzero(cEnd - cStart, seenCells));
435e1d83109SMatthew G. Knepley     ++comp;
436e1d83109SMatthew G. Knepley   } while (1);
437e1d83109SMatthew G. Knepley   numComponents = comp;
43884961fc4SMatthew G. Knepley   if (flg) {
43984961fc4SMatthew G. Knepley     PetscViewer v;
44084961fc4SMatthew G. Knepley 
4419566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIGetStdout(comm, &v));
4429566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(v));
4439566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank));
4449566063dSJacob Faibussowitsch     PetscCall(PetscBTView(cEnd - cStart, flippedCells, v));
4459566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(v));
4469566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(v));
44784961fc4SMatthew G. Knepley   }
44884961fc4SMatthew G. Knepley   /* Now all subdomains are oriented, but we need a consistent parallel orientation */
44984961fc4SMatthew G. Knepley   if (numLeaves >= 0) {
4503d6051bcSMatthew G. Knepley     PetscInt maxSupportSize, neighbor;
4513d6051bcSMatthew G. Knepley 
452e1d83109SMatthew G. Knepley     /* Store orientations of boundary faces*/
4533d6051bcSMatthew G. Knepley     PetscCall(DMPlexGetMaxSizes(dm, NULL, &maxSupportSize));
4543d6051bcSMatthew G. Knepley     PetscCall(PetscCalloc3(numRoots, &rorntComp, numRoots, &lorntComp, maxSupportSize, &locSupport));
45584961fc4SMatthew G. Knepley     for (face = fStart; face < fEnd; ++face) {
456e1d83109SMatthew G. Knepley       const PetscInt *cone, *support, *ornt;
4573d6051bcSMatthew G. Knepley       PetscInt        coneSize, supportSize, Ns = 0, s, l;
458e1d83109SMatthew G. Knepley 
4599566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupportSize(dm, face, &supportSize));
4603d6051bcSMatthew G. Knepley       /* Ignore overlapping cells */
4619566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupport(dm, face, &support));
4623d6051bcSMatthew G. Knepley       for (s = 0; s < supportSize; ++s) {
463658a0971SMatthew G. Knepley         if (lpoints) PetscCall(PetscFindInt(support[s], numLeaves, lpoints, &l));
464658a0971SMatthew G. Knepley         else {
465658a0971SMatthew G. Knepley           if (support[s] >= 0 && support[s] < numLeaves) l = support[s];
466658a0971SMatthew G. Knepley           else l = -1;
467658a0971SMatthew G. Knepley         }
4683d6051bcSMatthew G. Knepley         if (l >= 0) continue;
4693d6051bcSMatthew G. Knepley         locSupport[Ns++] = support[s];
4703d6051bcSMatthew G. Knepley       }
4713d6051bcSMatthew G. Knepley       if (Ns != 1) continue;
4723d6051bcSMatthew G. Knepley       neighbor = locSupport[0];
4733d6051bcSMatthew G. Knepley       PetscCall(DMPlexGetCone(dm, neighbor, &cone));
4743d6051bcSMatthew G. Knepley       PetscCall(DMPlexGetConeSize(dm, neighbor, &coneSize));
4753d6051bcSMatthew G. Knepley       PetscCall(DMPlexGetConeOrientation(dm, neighbor, &ornt));
4769371c9d4SSatish Balay       for (c = 0; c < coneSize; ++c)
4779371c9d4SSatish Balay         if (cone[c] == face) break;
47884961fc4SMatthew G. Knepley       if (dim == 1) {
47984961fc4SMatthew G. Knepley         /* Use cone position instead, shifted to -1 or 1 */
480835f2295SStefano Zampini         if (PetscBTLookup(flippedCells, neighbor - cStart)) rorntComp[face].rank = 1 - c * 2;
481835f2295SStefano Zampini         else rorntComp[face].rank = c * 2 - 1;
48284961fc4SMatthew G. Knepley       } else {
4833d6051bcSMatthew G. Knepley         if (PetscBTLookup(flippedCells, neighbor - cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 : 1;
484e1d83109SMatthew G. Knepley         else rorntComp[face].rank = ornt[c] < 0 ? 1 : -1;
48584961fc4SMatthew G. Knepley       }
486e1d83109SMatthew G. Knepley       rorntComp[face].index = faceComp[face - fStart];
48784961fc4SMatthew G. Knepley     }
488e1d83109SMatthew G. Knepley     /* Communicate boundary edge orientations */
4896497c311SBarry Smith     PetscCall(PetscSFBcastBegin(sf, MPIU_SF_NODE, rorntComp, lorntComp, MPI_REPLACE));
4906497c311SBarry Smith     PetscCall(PetscSFBcastEnd(sf, MPIU_SF_NODE, rorntComp, lorntComp, MPI_REPLACE));
491e1d83109SMatthew G. Knepley   }
492e1d83109SMatthew G. Knepley   /* Get process adjacency */
4939566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors));
494a17aa47eSToby Isaac   viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm));
4959566063dSJacob Faibussowitsch   if (flg2) PetscCall(PetscViewerASCIIPushSynchronized(viewer));
4969566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &selfviewer));
497e1d83109SMatthew G. Knepley   for (comp = 0; comp < numComponents; ++comp) {
498e1d83109SMatthew G. Knepley     PetscInt l, n;
49984961fc4SMatthew G. Knepley 
500e1d83109SMatthew G. Knepley     numNeighbors[comp] = 0;
5019566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp]));
502e1d83109SMatthew G. Knepley     /* I know this is p^2 time in general, but for bounded degree its alright */
503e1d83109SMatthew G. Knepley     for (l = 0; l < numLeaves; ++l) {
504658a0971SMatthew G. Knepley       const PetscInt face = lpoints ? lpoints[l] : l;
505e1d83109SMatthew G. Knepley 
506e1d83109SMatthew G. Knepley       /* Find a representative face (edge) separating pairs of procs */
5073d6051bcSMatthew G. Knepley       if ((face >= fStart) && (face < fEnd) && (faceComp[face - fStart] == comp) && rorntComp[face].rank) {
508835f2295SStefano Zampini         const PetscInt rrank = rpoints[l].rank;
509e1d83109SMatthew G. Knepley         const PetscInt rcomp = lorntComp[face].index;
510e1d83109SMatthew G. Knepley 
5119371c9d4SSatish Balay         for (n = 0; n < numNeighbors[comp]; ++n)
5129371c9d4SSatish Balay           if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break;
513e1d83109SMatthew G. Knepley         if (n >= numNeighbors[comp]) {
514e1d83109SMatthew G. Knepley           PetscInt supportSize;
515e1d83109SMatthew G. Knepley 
5169566063dSJacob Faibussowitsch           PetscCall(DMPlexGetSupportSize(dm, face, &supportSize));
5174bb7faabSMatthew G. Knepley           // We can have internal faces in the SF if we have cells in the SF
5184bb7faabSMatthew G. Knepley           if (supportSize > 1) continue;
5199371c9d4SSatish Balay           if (flg)
520835f2295SStefano Zampini             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,
521835f2295SStefano Zampini                                              rpoints[l].index, rrank, rcomp, lorntComp[face].rank));
522e1d83109SMatthew G. Knepley           neighbors[comp][numNeighbors[comp]++] = l;
523e1d83109SMatthew G. Knepley         }
524e1d83109SMatthew G. Knepley       }
525e1d83109SMatthew G. Knepley     }
526e1d83109SMatthew G. Knepley     totNeighbors += numNeighbors[comp];
527e1d83109SMatthew G. Knepley   }
5289566063dSJacob Faibussowitsch   PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &selfviewer));
5299566063dSJacob Faibussowitsch   if (flg2) PetscCall(PetscViewerASCIIPopSynchronized(viewer));
5309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match));
531e1d83109SMatthew G. Knepley   for (comp = 0, off = 0; comp < numComponents; ++comp) {
532e1d83109SMatthew G. Knepley     PetscInt n;
533e1d83109SMatthew G. Knepley 
534e1d83109SMatthew G. Knepley     for (n = 0; n < numNeighbors[comp]; ++n, ++off) {
535658a0971SMatthew G. Knepley       const PetscInt face = lpoints ? lpoints[neighbors[comp][n]] : neighbors[comp][n];
536e1d83109SMatthew G. Knepley       const PetscInt o    = rorntComp[face].rank * lorntComp[face].rank;
537e1d83109SMatthew G. Knepley 
538e1d83109SMatthew G. Knepley       if (o < 0) match[off] = PETSC_TRUE;
539e1d83109SMatthew G. Knepley       else if (o > 0) match[off] = PETSC_FALSE;
540835f2295SStefano Zampini       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);
541e1d83109SMatthew G. Knepley       nrankComp[off].rank  = rpoints[neighbors[comp][n]].rank;
542658a0971SMatthew G. Knepley       nrankComp[off].index = lorntComp[lpoints ? lpoints[neighbors[comp][n]] : neighbors[comp][n]].index;
543e1d83109SMatthew G. Knepley     }
5449566063dSJacob Faibussowitsch     PetscCall(PetscFree(neighbors[comp]));
54584961fc4SMatthew G. Knepley   }
54684961fc4SMatthew G. Knepley   /* Collect the graph on 0 */
547e1d83109SMatthew G. Knepley   if (numLeaves >= 0) {
54884961fc4SMatthew G. Knepley     Mat          G;
54984961fc4SMatthew G. Knepley     PetscBT      seenProcs, flippedProcs;
55084961fc4SMatthew G. Knepley     PetscInt    *procFIFO, pTop, pBottom;
5517cadcfe8SMatthew G. Knepley     PetscInt    *N          = NULL, *Noff;
552e1d83109SMatthew G. Knepley     PetscSFNode *adj        = NULL;
55384961fc4SMatthew G. Knepley     PetscBool   *val        = NULL;
5546497c311SBarry Smith     PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc, p, o, itotNeighbors;
5559852e123SBarry Smith     PetscMPIInt  size = 0;
55684961fc4SMatthew G. Knepley 
5579566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(numComponents, &flipped));
5589566063dSJacob Faibussowitsch     if (rank == 0) PetscCallMPI(MPI_Comm_size(comm, &size));
5599566063dSJacob Faibussowitsch     PetscCall(PetscCalloc4(size, &recvcounts, size + 1, &displs, size, &Nc, size + 1, &Noff));
5609566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Gather(&numComponents, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm));
561ad540459SPierre Jolivet     for (p = 0; p < size; ++p) displs[p + 1] = displs[p] + Nc[p];
5629566063dSJacob Faibussowitsch     if (rank == 0) PetscCall(PetscMalloc1(displs[size], &N));
5639566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Gatherv(numNeighbors, numComponents, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm));
5649852e123SBarry Smith     for (p = 0, o = 0; p < size; ++p) {
565e1d83109SMatthew G. Knepley       recvcounts[p] = 0;
566e1d83109SMatthew G. Knepley       for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o];
56784961fc4SMatthew G. Knepley       displs[p + 1] = displs[p] + recvcounts[p];
56884961fc4SMatthew G. Knepley     }
5699566063dSJacob Faibussowitsch     if (rank == 0) PetscCall(PetscMalloc2(displs[size], &adj, displs[size], &val));
5706497c311SBarry Smith     PetscCall(PetscMPIIntCast(totNeighbors, &itotNeighbors));
5716497c311SBarry Smith     PetscCallMPI(MPI_Gatherv(nrankComp, itotNeighbors, MPIU_SF_NODE, adj, recvcounts, displs, MPIU_SF_NODE, 0, comm));
572*5440e5dcSBarry Smith     PetscCallMPI(MPI_Gatherv(match, itotNeighbors, MPI_C_BOOL, val, recvcounts, displs, MPI_C_BOOL, 0, comm));
5739566063dSJacob Faibussowitsch     PetscCall(PetscFree2(numNeighbors, neighbors));
574dd400576SPatrick Sanan     if (rank == 0) {
575ad540459SPierre Jolivet       for (p = 1; p <= size; ++p) Noff[p] = Noff[p - 1] + Nc[p - 1];
576a83982f3SMatthew G. Knepley       if (flg) {
577e1d83109SMatthew G. Knepley         PetscInt n;
578e1d83109SMatthew G. Knepley 
5799852e123SBarry Smith         for (p = 0, off = 0; p < size; ++p) {
580e1d83109SMatthew G. Knepley           for (c = 0; c < Nc[p]; ++c) {
58163a3b9bcSJacob Faibussowitsch             PetscCall(PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %" PetscInt_FMT ":\n", p, c));
582835f2295SStefano Zampini             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]]));
58384961fc4SMatthew G. Knepley           }
58484961fc4SMatthew G. Knepley         }
58584961fc4SMatthew G. Knepley       }
58684961fc4SMatthew G. Knepley       /* Symmetrize the graph */
5879566063dSJacob Faibussowitsch       PetscCall(MatCreate(PETSC_COMM_SELF, &G));
5889566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size]));
5899566063dSJacob Faibussowitsch       PetscCall(MatSetUp(G));
5909852e123SBarry Smith       for (p = 0, off = 0; p < size; ++p) {
591e1d83109SMatthew G. Knepley         for (c = 0; c < Nc[p]; ++c) {
592e1d83109SMatthew G. Knepley           const PetscInt r = Noff[p] + c;
593e1d83109SMatthew G. Knepley           PetscInt       n;
594e1d83109SMatthew G. Knepley 
595e1d83109SMatthew G. Knepley           for (n = 0; n < N[r]; ++n, ++off) {
596e1d83109SMatthew G. Knepley             const PetscInt    q = Noff[adj[off].rank] + adj[off].index;
597e1d83109SMatthew G. Knepley             const PetscScalar o = val[off] ? 1.0 : 0.0;
59884961fc4SMatthew G. Knepley 
5999566063dSJacob Faibussowitsch             PetscCall(MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES));
6009566063dSJacob Faibussowitsch             PetscCall(MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES));
60184961fc4SMatthew G. Knepley           }
60284961fc4SMatthew G. Knepley         }
603e1d83109SMatthew G. Knepley       }
6049566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY));
6059566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY));
60684961fc4SMatthew G. Knepley 
6079566063dSJacob Faibussowitsch       PetscCall(PetscBTCreate(Noff[size], &seenProcs));
6089566063dSJacob Faibussowitsch       PetscCall(PetscBTMemzero(Noff[size], seenProcs));
6099566063dSJacob Faibussowitsch       PetscCall(PetscBTCreate(Noff[size], &flippedProcs));
6109566063dSJacob Faibussowitsch       PetscCall(PetscBTMemzero(Noff[size], flippedProcs));
6119566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(Noff[size], &procFIFO));
61284961fc4SMatthew G. Knepley       pTop = pBottom = 0;
6139852e123SBarry Smith       for (p = 0; p < Noff[size]; ++p) {
61484961fc4SMatthew G. Knepley         if (PetscBTLookup(seenProcs, p)) continue;
61584961fc4SMatthew G. Knepley         /* Initialize FIFO with next proc */
61684961fc4SMatthew G. Knepley         procFIFO[pBottom++] = p;
6179566063dSJacob Faibussowitsch         PetscCall(PetscBTSet(seenProcs, p));
61884961fc4SMatthew G. Knepley         /* Consider each proc in FIFO */
61984961fc4SMatthew G. Knepley         while (pTop < pBottom) {
62084961fc4SMatthew G. Knepley           const PetscScalar *ornt;
62184961fc4SMatthew G. Knepley           const PetscInt    *neighbors;
622e1d83109SMatthew G. Knepley           PetscInt           proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n;
62384961fc4SMatthew G. Knepley 
62484961fc4SMatthew G. Knepley           proc     = procFIFO[pTop++];
62584961fc4SMatthew G. Knepley           flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0;
6269566063dSJacob Faibussowitsch           PetscCall(MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt));
62784961fc4SMatthew G. Knepley           /* Loop over neighboring procs */
62884961fc4SMatthew G. Knepley           for (n = 0; n < numNeighbors; ++n) {
62984961fc4SMatthew G. Knepley             nproc    = neighbors[n];
63084961fc4SMatthew G. Knepley             mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1;
63184961fc4SMatthew G. Knepley             seen     = PetscBTLookup(seenProcs, nproc);
63284961fc4SMatthew G. Knepley             flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0;
63384961fc4SMatthew G. Knepley 
63484961fc4SMatthew G. Knepley             if (mismatch ^ (flippedA ^ flippedB)) {
63563a3b9bcSJacob 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);
636966bd95aSPierre Jolivet               PetscCheck(!flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
6379566063dSJacob Faibussowitsch               PetscCall(PetscBTSet(flippedProcs, nproc));
6381dca8a05SBarry Smith             } else PetscCheck(!mismatch || !flippedA || !flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
63984961fc4SMatthew G. Knepley             if (!seen) {
64084961fc4SMatthew G. Knepley               procFIFO[pBottom++] = nproc;
6419566063dSJacob Faibussowitsch               PetscCall(PetscBTSet(seenProcs, nproc));
64284961fc4SMatthew G. Knepley             }
64384961fc4SMatthew G. Knepley           }
64484961fc4SMatthew G. Knepley         }
64584961fc4SMatthew G. Knepley       }
6469566063dSJacob Faibussowitsch       PetscCall(PetscFree(procFIFO));
6479566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&G));
6489566063dSJacob Faibussowitsch       PetscCall(PetscFree2(adj, val));
6499566063dSJacob Faibussowitsch       PetscCall(PetscBTDestroy(&seenProcs));
65084961fc4SMatthew G. Knepley     }
651e1d83109SMatthew G. Knepley     /* Scatter flip flags */
652e1d83109SMatthew G. Knepley     {
653e1d83109SMatthew G. Knepley       PetscBool *flips = NULL;
654e1d83109SMatthew G. Knepley 
655dd400576SPatrick Sanan       if (rank == 0) {
6569566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(Noff[size], &flips));
6579852e123SBarry Smith         for (p = 0; p < Noff[size]; ++p) {
658e1d83109SMatthew G. Knepley           flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE;
6599566063dSJacob Faibussowitsch           if (flg && flips[p]) PetscCall(PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p));
660e1d83109SMatthew G. Knepley         }
661ad540459SPierre Jolivet         for (p = 0; p < size; ++p) displs[p + 1] = displs[p] + Nc[p];
662e1d83109SMatthew G. Knepley       }
663*5440e5dcSBarry Smith       PetscCallMPI(MPI_Scatterv(flips, Nc, displs, MPI_C_BOOL, flipped, numComponents, MPI_C_BOOL, 0, comm));
6649566063dSJacob Faibussowitsch       PetscCall(PetscFree(flips));
66584961fc4SMatthew G. Knepley     }
6669566063dSJacob Faibussowitsch     if (rank == 0) PetscCall(PetscBTDestroy(&flippedProcs));
6679566063dSJacob Faibussowitsch     PetscCall(PetscFree(N));
6689566063dSJacob Faibussowitsch     PetscCall(PetscFree4(recvcounts, displs, Nc, Noff));
6699566063dSJacob Faibussowitsch     PetscCall(PetscFree2(nrankComp, match));
670e1d83109SMatthew G. Knepley 
671e1d83109SMatthew G. Knepley     /* Decide whether to flip cells in each component */
6729371c9d4SSatish Balay     for (c = 0; c < cEnd - cStart; ++c) {
6739371c9d4SSatish Balay       if (flipped[cellComp[c]]) PetscCall(PetscBTNegate(flippedCells, c));
6749371c9d4SSatish Balay     }
6759566063dSJacob Faibussowitsch     PetscCall(PetscFree(flipped));
67684961fc4SMatthew G. Knepley   }
67784961fc4SMatthew G. Knepley   if (flg) {
67884961fc4SMatthew G. Knepley     PetscViewer v;
67984961fc4SMatthew G. Knepley 
6809566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIGetStdout(comm, &v));
6819566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(v));
6829566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank));
6839566063dSJacob Faibussowitsch     PetscCall(PetscBTView(cEnd - cStart, flippedCells, v));
6849566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(v));
6859566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(v));
68684961fc4SMatthew G. Knepley   }
68784961fc4SMatthew G. Knepley   /* Reverse flipped cells in the mesh */
68884961fc4SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
68948a46eb9SPierre Jolivet     if (PetscBTLookup(flippedCells, c - cStart)) PetscCall(DMPlexOrientPoint(dm, c, -1));
69084961fc4SMatthew G. Knepley   }
6919566063dSJacob Faibussowitsch   PetscCall(PetscBTDestroy(&seenCells));
6929566063dSJacob Faibussowitsch   PetscCall(PetscBTDestroy(&flippedCells));
6939566063dSJacob Faibussowitsch   PetscCall(PetscBTDestroy(&seenFaces));
6949566063dSJacob Faibussowitsch   PetscCall(PetscFree2(numNeighbors, neighbors));
6953d6051bcSMatthew G. Knepley   PetscCall(PetscFree3(rorntComp, lorntComp, locSupport));
6969566063dSJacob Faibussowitsch   PetscCall(PetscFree3(faceFIFO, cellComp, faceComp));
6973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6980851c46eSMatthew G. Knepley #endif
6990851c46eSMatthew G. Knepley }
7000851c46eSMatthew G. Knepley 
CreateCellAndFaceIS_Private(DM dm,DMLabel label,IS * cellIS,IS * faceIS)7010851c46eSMatthew G. Knepley static PetscErrorCode CreateCellAndFaceIS_Private(DM dm, DMLabel label, IS *cellIS, IS *faceIS)
7020851c46eSMatthew G. Knepley {
7030851c46eSMatthew G. Knepley   IS              valueIS;
7040851c46eSMatthew G. Knepley   const PetscInt *values;
7050851c46eSMatthew G. Knepley   PetscInt        Nv, depth = 0;
7060851c46eSMatthew G. Knepley 
7070851c46eSMatthew G. Knepley   PetscFunctionBegin;
7080851c46eSMatthew G. Knepley   PetscCall(DMLabelGetValueIS(label, &valueIS));
7090851c46eSMatthew G. Knepley   PetscCall(ISGetLocalSize(valueIS, &Nv));
7100851c46eSMatthew G. Knepley   PetscCall(ISGetIndices(valueIS, &values));
7110851c46eSMatthew G. Knepley   for (PetscInt v = 0; v < Nv; ++v) {
7120851c46eSMatthew G. Knepley     const PetscInt val = values[v] < 0 || values[v] >= 100 ? 0 : values[v];
713f369e074SMatthew G. Knepley     PetscInt       n;
7140851c46eSMatthew G. Knepley 
715f369e074SMatthew G. Knepley     PetscCall(DMLabelGetStratumSize(label, val, &n));
716f369e074SMatthew G. Knepley     if (!n) continue;
7170851c46eSMatthew G. Knepley     depth = PetscMax(val, depth);
7180851c46eSMatthew G. Knepley   }
7190851c46eSMatthew G. Knepley   PetscCall(ISDestroy(&valueIS));
7200851c46eSMatthew G. Knepley   PetscCheck(depth >= 1 || !Nv, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Depth for interface must be at least 1, not %" PetscInt_FMT, depth);
7210851c46eSMatthew G. Knepley   PetscCall(DMLabelGetStratumIS(label, depth, cellIS));
7220851c46eSMatthew G. Knepley   PetscCall(DMLabelGetStratumIS(label, depth - 1, faceIS));
72357508eceSPierre Jolivet   if (!*cellIS) PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, cellIS));
72457508eceSPierre Jolivet   if (!*faceIS) PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, faceIS));
7250851c46eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
7260851c46eSMatthew G. Knepley }
7270851c46eSMatthew G. Knepley 
DMPlexOrientLabel(DM dm,DMLabel label)7280851c46eSMatthew G. Knepley PetscErrorCode DMPlexOrientLabel(DM dm, DMLabel label)
7290851c46eSMatthew G. Knepley {
7300851c46eSMatthew G. Knepley   IS cellIS, faceIS;
7310851c46eSMatthew G. Knepley 
7320851c46eSMatthew G. Knepley   PetscFunctionBegin;
7330851c46eSMatthew G. Knepley   PetscCall(CreateCellAndFaceIS_Private(dm, label, &cellIS, &faceIS));
7340851c46eSMatthew G. Knepley   PetscCall(DMPlexOrientCells_Internal(dm, cellIS, faceIS));
7350851c46eSMatthew G. Knepley   PetscCall(ISDestroy(&cellIS));
7360851c46eSMatthew G. Knepley   PetscCall(ISDestroy(&faceIS));
7370851c46eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
7380851c46eSMatthew G. Knepley }
7390851c46eSMatthew G. Knepley 
DMPlexOrientCells_Internal(DM dm,IS cellIS,IS faceIS)7400851c46eSMatthew G. Knepley PetscErrorCode DMPlexOrientCells_Internal(DM dm, IS cellIS, IS faceIS)
7410851c46eSMatthew G. Knepley {
7420851c46eSMatthew G. Knepley   MPI_Comm           comm;
7430851c46eSMatthew G. Knepley   PetscSF            sf;
7440851c46eSMatthew G. Knepley   const PetscInt    *lpoints;
7450851c46eSMatthew G. Knepley   const PetscSFNode *rpoints;
7460851c46eSMatthew G. Knepley   PetscSFNode       *rorntComp = NULL, *lorntComp = NULL;
7470851c46eSMatthew G. Knepley   PetscInt          *numNeighbors, **neighbors, *locSupp = NULL;
7480851c46eSMatthew G. Knepley   PetscSFNode       *nrankComp;
7490851c46eSMatthew G. Knepley   PetscBool         *match, *flipped;
7500851c46eSMatthew G. Knepley   PetscBT            flippedCells;
7510851c46eSMatthew G. Knepley   PetscInt          *cellComp, *faceComp;
7520851c46eSMatthew G. Knepley   const PetscInt    *cells = NULL, *faces = NULL;
7530851c46eSMatthew G. Knepley   PetscInt           cStart = 0, cEnd = 0, fStart = 0, fEnd = 0;
7540851c46eSMatthew G. Knepley   PetscInt           numLeaves, numRoots, dim, Ncomp, totNeighbors = 0;
7550851c46eSMatthew G. Knepley   PetscMPIInt        rank, size;
7560851c46eSMatthew G. Knepley   PetscBool          view, viewSync;
7570851c46eSMatthew G. Knepley   PetscViewer        viewer = NULL, selfviewer = NULL;
7580851c46eSMatthew G. Knepley 
7590851c46eSMatthew G. Knepley   PetscFunctionBegin;
7600851c46eSMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
7610851c46eSMatthew G. Knepley   PetscCallMPI(MPI_Comm_rank(comm, &rank));
7620851c46eSMatthew G. Knepley   PetscCallMPI(MPI_Comm_size(comm, &size));
7630851c46eSMatthew G. Knepley   PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-orientation_view", &view));
7640851c46eSMatthew G. Knepley   PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-orientation_view_synchronized", &viewSync));
7650851c46eSMatthew G. Knepley 
7660851c46eSMatthew G. Knepley   if (cellIS) PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
7670851c46eSMatthew G. Knepley   if (faceIS) PetscCall(ISGetPointRange(faceIS, &fStart, &fEnd, &faces));
7680851c46eSMatthew G. Knepley   PetscCall(DMGetPointSF(dm, &sf));
7690851c46eSMatthew G. Knepley   PetscCall(PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints));
7700851c46eSMatthew G. Knepley   /* Truth Table
7710851c46eSMatthew G. Knepley      mismatch    flips   do action   mismatch   flipA ^ flipB   action
7720851c46eSMatthew G. Knepley          F       0 flips     no         F             F           F
7730851c46eSMatthew G. Knepley          F       1 flip      yes        F             T           T
7740851c46eSMatthew G. Knepley          F       2 flips     no         T             F           T
7750851c46eSMatthew G. Knepley          T       0 flips     yes        T             T           F
7760851c46eSMatthew G. Knepley          T       1 flip      no
7770851c46eSMatthew G. Knepley          T       2 flips     yes
7780851c46eSMatthew G. Knepley   */
7790851c46eSMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
7800851c46eSMatthew G. Knepley   PetscCall(PetscBTCreate(cEnd - cStart, &flippedCells));
7810851c46eSMatthew G. Knepley   PetscCall(PetscBTMemzero(cEnd - cStart, flippedCells));
7820851c46eSMatthew G. Knepley   PetscCall(PetscCalloc2(cEnd - cStart, &cellComp, fEnd - fStart, &faceComp));
7830851c46eSMatthew G. Knepley   /*
7840851c46eSMatthew G. Knepley    OLD STYLE
7850851c46eSMatthew G. Knepley    - Add an integer array over cells and faces (component) for connected component number
7860851c46eSMatthew G. Knepley    Foreach component
7870851c46eSMatthew G. Knepley      - Mark the initial cell as seen
7880851c46eSMatthew G. Knepley      - Process component as usual
7890851c46eSMatthew G. Knepley      - Set component for all seenCells
7900851c46eSMatthew G. Knepley      - Wipe seenCells and seenFaces (flippedCells can stay)
7910851c46eSMatthew G. Knepley    - Generate parallel adjacency for component using SF and seenFaces
7920851c46eSMatthew G. Knepley    - Collect Ncomp adj data from each proc to 0
7930851c46eSMatthew G. Knepley    - Build same serial graph
7940851c46eSMatthew G. Knepley    - Use same solver
7950851c46eSMatthew G. Knepley    - Use Scatterv to send back flipped flags for each component
7960851c46eSMatthew G. Knepley    - Negate flippedCells by component
7970851c46eSMatthew G. Knepley 
7980851c46eSMatthew G. Knepley    NEW STYLE
7990851c46eSMatthew G. Knepley    - Create the adj on each process
8000851c46eSMatthew G. Knepley    - Bootstrap to complete graph on proc 0
8010851c46eSMatthew G. Knepley   */
8020851c46eSMatthew G. Knepley   PetscCall(DMPlexOrient_Serial(dm, cellIS, faceIS, &Ncomp, cellComp, faceComp, flippedCells));
8030851c46eSMatthew G. Knepley   if (view) {
8040851c46eSMatthew G. Knepley     PetscViewer v;
80576dc4c7dSMatthew G. Knepley     PetscInt    cdepth = -1;
8060851c46eSMatthew G. Knepley 
8070851c46eSMatthew G. Knepley     PetscCall(PetscViewerASCIIGetStdout(comm, &v));
8080851c46eSMatthew G. Knepley     PetscCall(PetscViewerASCIIPushSynchronized(v));
80976dc4c7dSMatthew G. Knepley     if (cEnd > cStart) PetscCall(DMPlexGetPointDepth(dm, cells ? cells[cStart] : cStart, &cdepth));
81076dc4c7dSMatthew G. Knepley     PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]New Orientation %" PetscInt_FMT " cells (depth %" PetscInt_FMT ") and %" PetscInt_FMT " faces\n", rank, cEnd - cStart, cdepth, fEnd - fStart));
8110851c46eSMatthew G. Knepley     PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank));
8120851c46eSMatthew G. Knepley     PetscCall(PetscBTView(cEnd - cStart, flippedCells, v));
8130851c46eSMatthew G. Knepley     PetscCall(PetscViewerFlush(v));
8140851c46eSMatthew G. Knepley     PetscCall(PetscViewerASCIIPopSynchronized(v));
8150851c46eSMatthew G. Knepley   }
8160851c46eSMatthew G. Knepley   /* Now all subdomains are oriented, but we need a consistent parallel orientation */
8170851c46eSMatthew G. Knepley   // TODO: This all has to be rewritten to filter cones/supports to the ISes
8180851c46eSMatthew G. Knepley   if (numLeaves >= 0) {
8190851c46eSMatthew G. Knepley     PetscInt maxSuppSize, neighbor;
8200851c46eSMatthew G. Knepley 
8210851c46eSMatthew G. Knepley     // Store orientations of boundary faces
8220851c46eSMatthew G. Knepley     PetscCall(DMPlexGetMaxSizes(dm, NULL, &maxSuppSize));
8230851c46eSMatthew G. Knepley     PetscCall(PetscCalloc3(numRoots, &rorntComp, numRoots, &lorntComp, maxSuppSize, &locSupp));
8240851c46eSMatthew G. Knepley     for (PetscInt f = fStart; f < fEnd; ++f) {
8250851c46eSMatthew G. Knepley       const PetscInt  face = faces ? faces[f] : f;
8260851c46eSMatthew G. Knepley       const PetscInt *cone, *supp, *ornt;
8270851c46eSMatthew G. Knepley       PetscInt        coneSize, suppSize, nind, c, Ns = 0;
8280851c46eSMatthew G. Knepley 
8290851c46eSMatthew G. Knepley       PetscCall(DMPlexGetSupportSize(dm, face, &suppSize));
8300851c46eSMatthew G. Knepley       PetscCall(DMPlexGetSupport(dm, face, &supp));
8310851c46eSMatthew G. Knepley       for (PetscInt s = 0; s < suppSize; ++s) {
8320851c46eSMatthew G. Knepley         PetscInt ind, l;
8330851c46eSMatthew G. Knepley 
8340851c46eSMatthew G. Knepley         // Filter support
8350851c46eSMatthew G. Knepley         ind = GetPointIndex(supp[s], cStart, cEnd, cells);
8360851c46eSMatthew G. Knepley         if (ind < 0) continue;
8370851c46eSMatthew G. Knepley         // Ignore overlapping cells
8380851c46eSMatthew G. Knepley         PetscCall(PetscFindInt(supp[s], numLeaves, lpoints, &l));
8390851c46eSMatthew G. Knepley         if (l >= 0) continue;
8400851c46eSMatthew G. Knepley         locSupp[Ns++] = supp[s];
8410851c46eSMatthew G. Knepley       }
84217d25f2dSMatthew G. Knepley       PetscCheck(Ns < maxSuppSize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Index %" PetscInt_FMT " exceeds array size %" PetscInt_FMT, Ns, maxSuppSize);
8430851c46eSMatthew G. Knepley       if (Ns != 1) continue;
8440851c46eSMatthew G. Knepley       neighbor = locSupp[0];
8450851c46eSMatthew G. Knepley       nind     = GetPointIndex(neighbor, cStart, cEnd, cells);
8460851c46eSMatthew G. Knepley       PetscCall(DMPlexGetCone(dm, neighbor, &cone));
8470851c46eSMatthew G. Knepley       PetscCall(DMPlexGetConeSize(dm, neighbor, &coneSize));
8480851c46eSMatthew G. Knepley       PetscCall(DMPlexGetConeOrientation(dm, neighbor, &ornt));
8490851c46eSMatthew G. Knepley       for (c = 0; c < coneSize; ++c)
8500851c46eSMatthew G. Knepley         if (cone[c] == face) break;
8510851c46eSMatthew G. Knepley       if (dim == 1) {
8520851c46eSMatthew G. Knepley         /* Use cone position instead, shifted to -1 or 1 */
853835f2295SStefano Zampini         if (PetscBTLookup(flippedCells, nind)) rorntComp[face].rank = 1 - c * 2;
854835f2295SStefano Zampini         else rorntComp[face].rank = c * 2 - 1;
8550851c46eSMatthew G. Knepley       } else {
8560851c46eSMatthew G. Knepley         if (PetscBTLookup(flippedCells, nind)) rorntComp[face].rank = ornt[c] < 0 ? -1 : 1;
8570851c46eSMatthew G. Knepley         else rorntComp[face].rank = ornt[c] < 0 ? 1 : -1;
8580851c46eSMatthew G. Knepley       }
8590851c46eSMatthew G. Knepley       rorntComp[face].index = faceComp[GetPointIndex(face, fStart, fEnd, faces)];
8600851c46eSMatthew G. Knepley     }
8610851c46eSMatthew G. Knepley     // Communicate boundary edge orientations
8626497c311SBarry Smith     PetscCall(PetscSFBcastBegin(sf, MPIU_SF_NODE, rorntComp, lorntComp, MPI_REPLACE));
8636497c311SBarry Smith     PetscCall(PetscSFBcastEnd(sf, MPIU_SF_NODE, rorntComp, lorntComp, MPI_REPLACE));
8640851c46eSMatthew G. Knepley   }
8650851c46eSMatthew G. Knepley   /* Get process adjacency */
8660851c46eSMatthew G. Knepley   PetscCall(PetscMalloc2(Ncomp, &numNeighbors, Ncomp, &neighbors));
8670851c46eSMatthew G. Knepley   viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm));
8680851c46eSMatthew G. Knepley   if (viewSync) PetscCall(PetscViewerASCIIPushSynchronized(viewer));
8690851c46eSMatthew G. Knepley   PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &selfviewer));
8700851c46eSMatthew G. Knepley   for (PetscInt comp = 0; comp < Ncomp; ++comp) {
8710851c46eSMatthew G. Knepley     PetscInt n;
8720851c46eSMatthew G. Knepley 
8730851c46eSMatthew G. Knepley     numNeighbors[comp] = 0;
8740851c46eSMatthew G. Knepley     PetscCall(PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp]));
8750851c46eSMatthew G. Knepley     /* I know this is p^2 time in general, but for bounded degree its alright */
8760851c46eSMatthew G. Knepley     for (PetscInt l = 0; l < numLeaves; ++l) {
8770851c46eSMatthew G. Knepley       const PetscInt face = lpoints[l];
8780851c46eSMatthew G. Knepley       PetscInt       find;
8790851c46eSMatthew G. Knepley 
8800851c46eSMatthew G. Knepley       /* Find a representative face (edge) separating pairs of procs */
8810851c46eSMatthew G. Knepley       find = GetPointIndex(face, fStart, fEnd, faces);
8820851c46eSMatthew G. Knepley       if ((find >= 0) && (faceComp[find] == comp) && rorntComp[face].rank) {
8830851c46eSMatthew G. Knepley         const PetscInt rrank = rpoints[l].rank;
8840851c46eSMatthew G. Knepley         const PetscInt rcomp = lorntComp[face].index;
8850851c46eSMatthew G. Knepley 
8860851c46eSMatthew G. Knepley         for (n = 0; n < numNeighbors[comp]; ++n)
8870851c46eSMatthew G. Knepley           if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break;
8880851c46eSMatthew G. Knepley         if (n >= numNeighbors[comp]) {
8890851c46eSMatthew G. Knepley           const PetscInt *supp;
8900851c46eSMatthew G. Knepley           PetscInt        suppSize, Ns = 0;
8910851c46eSMatthew G. Knepley 
8920851c46eSMatthew G. Knepley           PetscCall(DMPlexGetSupport(dm, face, &supp));
8930851c46eSMatthew G. Knepley           PetscCall(DMPlexGetSupportSize(dm, face, &suppSize));
8940851c46eSMatthew G. Knepley           for (PetscInt s = 0; s < suppSize; ++s) {
8950851c46eSMatthew G. Knepley             // Filter support
8960851c46eSMatthew G. Knepley             if (GetPointIndex(supp[s], cStart, cEnd, cells) >= 0) ++Ns;
8970851c46eSMatthew G. Knepley           }
8980851c46eSMatthew G. Knepley           PetscCheck(Ns == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary face %" PetscInt_FMT " should see one cell, not %" PetscInt_FMT, face, Ns);
8990851c46eSMatthew G. Knepley           if (view)
900835f2295SStefano Zampini             PetscCall(PetscViewerASCIIPrintf(selfviewer, "[%d]: component %" PetscInt_FMT ", 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,
901835f2295SStefano Zampini                                              rpoints[l].index, rrank, rcomp, lorntComp[face].rank));
9020851c46eSMatthew G. Knepley           neighbors[comp][numNeighbors[comp]++] = l;
9030851c46eSMatthew G. Knepley         }
9040851c46eSMatthew G. Knepley       }
9050851c46eSMatthew G. Knepley     }
9060851c46eSMatthew G. Knepley     totNeighbors += numNeighbors[comp];
9070851c46eSMatthew G. Knepley   }
9080851c46eSMatthew G. Knepley   PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &selfviewer));
9090851c46eSMatthew G. Knepley   if (viewSync) PetscCall(PetscViewerASCIIPopSynchronized(viewer));
9100851c46eSMatthew G. Knepley   PetscCall(PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match));
9110851c46eSMatthew G. Knepley   for (PetscInt comp = 0, off = 0; comp < Ncomp; ++comp) {
9120851c46eSMatthew G. Knepley     for (PetscInt n = 0; n < numNeighbors[comp]; ++n, ++off) {
9130851c46eSMatthew G. Knepley       const PetscInt face = lpoints[neighbors[comp][n]];
9140851c46eSMatthew G. Knepley       const PetscInt o    = rorntComp[face].rank * lorntComp[face].rank;
9150851c46eSMatthew G. Knepley 
9160851c46eSMatthew G. Knepley       if (o < 0) match[off] = PETSC_TRUE;
9170851c46eSMatthew G. Knepley       else if (o > 0) match[off] = PETSC_FALSE;
9180851c46eSMatthew G. Knepley       else
919835f2295SStefano Zampini         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid face %" PetscInt_FMT " (%" PetscInt_FMT ", %" PetscInt_FMT ") neighbor: %" PetscInt_FMT " comp: %" PetscInt_FMT, face, rorntComp[face].rank, lorntComp[face].rank, neighbors[comp][n], comp);
9200851c46eSMatthew G. Knepley       nrankComp[off].rank  = rpoints[neighbors[comp][n]].rank;
9210851c46eSMatthew G. Knepley       nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index;
9220851c46eSMatthew G. Knepley     }
9230851c46eSMatthew G. Knepley     PetscCall(PetscFree(neighbors[comp]));
9240851c46eSMatthew G. Knepley   }
9250851c46eSMatthew G. Knepley   /* Collect the graph on 0 */
9260851c46eSMatthew G. Knepley   if (numLeaves >= 0) {
9270851c46eSMatthew G. Knepley     Mat          G;
9280851c46eSMatthew G. Knepley     PetscBT      seenProcs, flippedProcs;
9290851c46eSMatthew G. Knepley     PetscInt    *procFIFO, pTop, pBottom;
9300851c46eSMatthew G. Knepley     PetscInt    *N          = NULL, *Noff;
9310851c46eSMatthew G. Knepley     PetscSFNode *adj        = NULL;
9320851c46eSMatthew G. Knepley     PetscBool   *val        = NULL;
9330851c46eSMatthew G. Knepley     PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc;
9346497c311SBarry Smith     PetscMPIInt  size = 0, iNcomp, itotNeighbors;
9350851c46eSMatthew G. Knepley 
9360851c46eSMatthew G. Knepley     PetscCall(PetscCalloc1(Ncomp, &flipped));
9370851c46eSMatthew G. Knepley     if (rank == 0) PetscCallMPI(MPI_Comm_size(comm, &size));
9380851c46eSMatthew G. Knepley     PetscCall(PetscCalloc4(size, &recvcounts, size + 1, &displs, size, &Nc, size + 1, &Noff));
9390851c46eSMatthew G. Knepley     PetscCallMPI(MPI_Gather(&Ncomp, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm));
9400851c46eSMatthew G. Knepley     for (PetscInt p = 0; p < size; ++p) displs[p + 1] = displs[p] + Nc[p];
9410851c46eSMatthew G. Knepley     if (rank == 0) PetscCall(PetscMalloc1(displs[size], &N));
9426497c311SBarry Smith     PetscCall(PetscMPIIntCast(Ncomp, &iNcomp));
9436497c311SBarry Smith     PetscCallMPI(MPI_Gatherv(numNeighbors, iNcomp, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm));
9440851c46eSMatthew G. Knepley     for (PetscInt p = 0, o = 0; p < size; ++p) {
9450851c46eSMatthew G. Knepley       recvcounts[p] = 0;
9460851c46eSMatthew G. Knepley       for (PetscInt c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o];
9470851c46eSMatthew G. Knepley       displs[p + 1] = displs[p] + recvcounts[p];
9480851c46eSMatthew G. Knepley     }
9490851c46eSMatthew G. Knepley     if (rank == 0) PetscCall(PetscMalloc2(displs[size], &adj, displs[size], &val));
9506497c311SBarry Smith     PetscCall(PetscMPIIntCast(totNeighbors, &itotNeighbors));
9516497c311SBarry Smith     PetscCallMPI(MPI_Gatherv(nrankComp, itotNeighbors, MPIU_SF_NODE, adj, recvcounts, displs, MPIU_SF_NODE, 0, comm));
952*5440e5dcSBarry Smith     PetscCallMPI(MPI_Gatherv(match, itotNeighbors, MPI_C_BOOL, val, recvcounts, displs, MPI_C_BOOL, 0, comm));
9530851c46eSMatthew G. Knepley     PetscCall(PetscFree2(numNeighbors, neighbors));
9540851c46eSMatthew G. Knepley     if (rank == 0) {
9550851c46eSMatthew G. Knepley       for (PetscInt p = 1; p <= size; ++p) Noff[p] = Noff[p - 1] + Nc[p - 1];
9560851c46eSMatthew G. Knepley       if (view) {
9570851c46eSMatthew G. Knepley         for (PetscInt p = 0, off = 0; p < size; ++p) {
9580851c46eSMatthew G. Knepley           for (PetscInt c = 0; c < Nc[p]; ++c) {
9590851c46eSMatthew G. Knepley             PetscCall(PetscPrintf(PETSC_COMM_SELF, "Proc %" PetscInt_FMT " Comp %" PetscInt_FMT ":\n", p, c));
960835f2295SStefano Zampini             for (PetscInt 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]]));
9610851c46eSMatthew G. Knepley           }
9620851c46eSMatthew G. Knepley         }
9630851c46eSMatthew G. Knepley       }
9640851c46eSMatthew G. Knepley       /* Symmetrize the graph */
9650851c46eSMatthew G. Knepley       PetscCall(MatCreate(PETSC_COMM_SELF, &G));
9660851c46eSMatthew G. Knepley       PetscCall(MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size]));
9670851c46eSMatthew G. Knepley       PetscCall(MatSetUp(G));
9680851c46eSMatthew G. Knepley       for (PetscInt p = 0, off = 0; p < size; ++p) {
9690851c46eSMatthew G. Knepley         for (PetscInt c = 0; c < Nc[p]; ++c) {
9700851c46eSMatthew G. Knepley           const PetscInt r = Noff[p] + c;
9710851c46eSMatthew G. Knepley 
9720851c46eSMatthew G. Knepley           for (PetscInt n = 0; n < N[r]; ++n, ++off) {
9730851c46eSMatthew G. Knepley             const PetscInt    q = Noff[adj[off].rank] + adj[off].index;
9740851c46eSMatthew G. Knepley             const PetscScalar o = val[off] ? 1.0 : 0.0;
9750851c46eSMatthew G. Knepley 
9760851c46eSMatthew G. Knepley             PetscCall(MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES));
9770851c46eSMatthew G. Knepley             PetscCall(MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES));
9780851c46eSMatthew G. Knepley           }
9790851c46eSMatthew G. Knepley         }
9800851c46eSMatthew G. Knepley       }
9810851c46eSMatthew G. Knepley       PetscCall(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY));
9820851c46eSMatthew G. Knepley       PetscCall(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY));
9830851c46eSMatthew G. Knepley 
9840851c46eSMatthew G. Knepley       PetscCall(PetscBTCreate(Noff[size], &seenProcs));
9850851c46eSMatthew G. Knepley       PetscCall(PetscBTMemzero(Noff[size], seenProcs));
9860851c46eSMatthew G. Knepley       PetscCall(PetscBTCreate(Noff[size], &flippedProcs));
9870851c46eSMatthew G. Knepley       PetscCall(PetscBTMemzero(Noff[size], flippedProcs));
9880851c46eSMatthew G. Knepley       PetscCall(PetscMalloc1(Noff[size], &procFIFO));
9890851c46eSMatthew G. Knepley       pTop = pBottom = 0;
9900851c46eSMatthew G. Knepley       for (PetscInt p = 0; p < Noff[size]; ++p) {
9910851c46eSMatthew G. Knepley         if (PetscBTLookup(seenProcs, p)) continue;
9920851c46eSMatthew G. Knepley         /* Initialize FIFO with next proc */
9930851c46eSMatthew G. Knepley         procFIFO[pBottom++] = p;
9940851c46eSMatthew G. Knepley         PetscCall(PetscBTSet(seenProcs, p));
9950851c46eSMatthew G. Knepley         /* Consider each proc in FIFO */
9960851c46eSMatthew G. Knepley         while (pTop < pBottom) {
9970851c46eSMatthew G. Knepley           const PetscScalar *ornt;
9980851c46eSMatthew G. Knepley           const PetscInt    *neighbors;
9990851c46eSMatthew G. Knepley           PetscInt           proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors;
10000851c46eSMatthew G. Knepley 
10010851c46eSMatthew G. Knepley           proc     = procFIFO[pTop++];
10020851c46eSMatthew G. Knepley           flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0;
10030851c46eSMatthew G. Knepley           PetscCall(MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt));
10040851c46eSMatthew G. Knepley           /* Loop over neighboring procs */
10050851c46eSMatthew G. Knepley           for (PetscInt n = 0; n < numNeighbors; ++n) {
10060851c46eSMatthew G. Knepley             nproc    = neighbors[n];
10070851c46eSMatthew G. Knepley             mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1;
10080851c46eSMatthew G. Knepley             seen     = PetscBTLookup(seenProcs, nproc);
10090851c46eSMatthew G. Knepley             flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0;
10100851c46eSMatthew G. Knepley 
10110851c46eSMatthew G. Knepley             if (mismatch ^ (flippedA ^ flippedB)) {
10120851c46eSMatthew G. Knepley               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);
1013966bd95aSPierre Jolivet               PetscCheck(!flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
10140851c46eSMatthew G. Knepley               PetscCall(PetscBTSet(flippedProcs, nproc));
10150851c46eSMatthew G. Knepley             } else PetscCheck(!mismatch || !flippedA || !flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
10160851c46eSMatthew G. Knepley             if (!seen) {
10170851c46eSMatthew G. Knepley               procFIFO[pBottom++] = nproc;
10180851c46eSMatthew G. Knepley               PetscCall(PetscBTSet(seenProcs, nproc));
10190851c46eSMatthew G. Knepley             }
10200851c46eSMatthew G. Knepley           }
10210851c46eSMatthew G. Knepley         }
10220851c46eSMatthew G. Knepley       }
10230851c46eSMatthew G. Knepley       PetscCall(PetscFree(procFIFO));
10240851c46eSMatthew G. Knepley       PetscCall(MatDestroy(&G));
10250851c46eSMatthew G. Knepley       PetscCall(PetscFree2(adj, val));
10260851c46eSMatthew G. Knepley       PetscCall(PetscBTDestroy(&seenProcs));
10270851c46eSMatthew G. Knepley     }
10280851c46eSMatthew G. Knepley     /* Scatter flip flags */
10290851c46eSMatthew G. Knepley     {
10300851c46eSMatthew G. Knepley       PetscBool *flips = NULL;
10310851c46eSMatthew G. Knepley 
10320851c46eSMatthew G. Knepley       if (rank == 0) {
10330851c46eSMatthew G. Knepley         PetscCall(PetscMalloc1(Noff[size], &flips));
10340851c46eSMatthew G. Knepley         for (PetscInt p = 0; p < Noff[size]; ++p) {
10350851c46eSMatthew G. Knepley           flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE;
10360851c46eSMatthew G. Knepley           if (view && flips[p]) PetscCall(PetscPrintf(comm, "Flipping Proc+Comp %" PetscInt_FMT ":\n", p));
10370851c46eSMatthew G. Knepley         }
10380851c46eSMatthew G. Knepley         for (PetscInt p = 0; p < size; ++p) displs[p + 1] = displs[p] + Nc[p];
10390851c46eSMatthew G. Knepley       }
10406497c311SBarry Smith       PetscCall(PetscMPIIntCast(Ncomp, &iNcomp));
1041*5440e5dcSBarry Smith       PetscCallMPI(MPI_Scatterv(flips, Nc, displs, MPI_C_BOOL, flipped, iNcomp, MPI_C_BOOL, 0, comm));
10420851c46eSMatthew G. Knepley       PetscCall(PetscFree(flips));
10430851c46eSMatthew G. Knepley     }
10440851c46eSMatthew G. Knepley     if (rank == 0) PetscCall(PetscBTDestroy(&flippedProcs));
10450851c46eSMatthew G. Knepley     PetscCall(PetscFree(N));
10460851c46eSMatthew G. Knepley     PetscCall(PetscFree4(recvcounts, displs, Nc, Noff));
10470851c46eSMatthew G. Knepley     PetscCall(PetscFree2(nrankComp, match));
10480851c46eSMatthew G. Knepley 
10490851c46eSMatthew G. Knepley     /* Decide whether to flip cells in each component */
10500851c46eSMatthew G. Knepley     for (PetscInt c = 0; c < cEnd - cStart; ++c) {
10510851c46eSMatthew G. Knepley       if (flipped[cellComp[c]]) PetscCall(PetscBTNegate(flippedCells, c));
10520851c46eSMatthew G. Knepley     }
10530851c46eSMatthew G. Knepley     PetscCall(PetscFree(flipped));
10540851c46eSMatthew G. Knepley   }
10550851c46eSMatthew G. Knepley   if (view) {
10560851c46eSMatthew G. Knepley     PetscViewer v;
10570851c46eSMatthew G. Knepley 
10580851c46eSMatthew G. Knepley     PetscCall(PetscViewerASCIIGetStdout(comm, &v));
10590851c46eSMatthew G. Knepley     PetscCall(PetscViewerASCIIPushSynchronized(v));
10600851c46eSMatthew G. Knepley     PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank));
10610851c46eSMatthew G. Knepley     PetscCall(PetscBTView(cEnd - cStart, flippedCells, v));
10620851c46eSMatthew G. Knepley     PetscCall(PetscViewerFlush(v));
10630851c46eSMatthew G. Knepley     PetscCall(PetscViewerASCIIPopSynchronized(v));
10640851c46eSMatthew G. Knepley   }
10650851c46eSMatthew G. Knepley   // Reverse flipped cells in the mesh
10660851c46eSMatthew G. Knepley   PetscViewer     v;
10670851c46eSMatthew G. Knepley   const PetscInt *degree = NULL;
10680851c46eSMatthew G. Knepley   PetscInt       *points;
10690851c46eSMatthew G. Knepley   PetscInt        pStart, pEnd;
10700851c46eSMatthew G. Knepley 
10710851c46eSMatthew G. Knepley   if (view) {
10720851c46eSMatthew G. Knepley     PetscCall(PetscViewerASCIIGetStdout(comm, &v));
10730851c46eSMatthew G. Knepley     PetscCall(PetscViewerASCIIPushSynchronized(v));
10740851c46eSMatthew G. Knepley   }
10750851c46eSMatthew G. Knepley   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
10760851c46eSMatthew G. Knepley   if (numRoots >= 0) {
10770851c46eSMatthew G. Knepley     PetscCall(PetscSFComputeDegreeBegin(sf, &degree));
10780851c46eSMatthew G. Knepley     PetscCall(PetscSFComputeDegreeEnd(sf, &degree));
10790851c46eSMatthew G. Knepley   }
10800851c46eSMatthew G. Knepley   PetscCall(PetscCalloc1(pEnd - pStart, &points));
10810851c46eSMatthew G. Knepley   for (PetscInt c = cStart; c < cEnd; ++c) {
10820851c46eSMatthew G. Knepley     if (PetscBTLookup(flippedCells, c - cStart)) {
10830851c46eSMatthew G. Knepley       const PetscInt cell = cells ? cells[c] : c;
10840851c46eSMatthew G. Knepley 
10850851c46eSMatthew G. Knepley       PetscCall(DMPlexOrientPoint(dm, cell, -1));
10860851c46eSMatthew G. Knepley       if (degree && degree[cell]) points[cell] = 1;
10870851c46eSMatthew G. Knepley       if (view) PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]Flipping cell %" PetscInt_FMT "%s\n", rank, cell, degree && degree[cell] ? " and sending to overlap" : ""));
10880851c46eSMatthew G. Knepley     }
10890851c46eSMatthew G. Knepley   }
10900851c46eSMatthew G. Knepley   // Must propagate flips for cells in the overlap
10910851c46eSMatthew G. Knepley   if (numRoots >= 0) {
10920851c46eSMatthew G. Knepley     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, points, points, MPI_SUM));
10930851c46eSMatthew G. Knepley     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, points, points, MPI_SUM));
10940851c46eSMatthew G. Knepley   }
10950851c46eSMatthew G. Knepley   for (PetscInt c = cStart; c < cEnd; ++c) {
10960851c46eSMatthew G. Knepley     const PetscInt cell = cells ? cells[c] : c;
10970851c46eSMatthew G. Knepley 
10980851c46eSMatthew G. Knepley     if (points[cell] && !PetscBTLookup(flippedCells, c - cStart)) {
10990851c46eSMatthew G. Knepley       PetscCall(DMPlexOrientPoint(dm, cell, -1));
11000851c46eSMatthew G. Knepley       if (view) PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]Flipping cell %" PetscInt_FMT " through overlap\n", rank, cell));
11010851c46eSMatthew G. Knepley     }
11020851c46eSMatthew G. Knepley   }
11030851c46eSMatthew G. Knepley   if (view) {
11040851c46eSMatthew G. Knepley     PetscCall(PetscViewerFlush(v));
11050851c46eSMatthew G. Knepley     PetscCall(PetscViewerASCIIPopSynchronized(v));
11060851c46eSMatthew G. Knepley   }
11070851c46eSMatthew G. Knepley   PetscCall(PetscFree(points));
11080851c46eSMatthew G. Knepley   PetscCall(PetscBTDestroy(&flippedCells));
11090851c46eSMatthew G. Knepley   PetscCall(PetscFree2(numNeighbors, neighbors));
11100851c46eSMatthew G. Knepley   PetscCall(PetscFree3(rorntComp, lorntComp, locSupp));
11110851c46eSMatthew G. Knepley   PetscCall(PetscFree2(cellComp, faceComp));
11120851c46eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
111384961fc4SMatthew G. Knepley }
1114