xref: /petsc/src/dm/impls/plex/plexorient.c (revision 27d0e99a0845f279d25f36eb29f4a4212ca97ece)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
284961fc4SMatthew G. Knepley #include <petscsf.h>
384961fc4SMatthew G. Knepley 
484961fc4SMatthew G. Knepley /*@
5*27d0e99aSVaclav Hapla   DMPlexCompareOrientations - Compare the cone of the given DAG point (cell) with the given reference cone (with the same cone points modulo order), and return relative orientation.
6*27d0e99aSVaclav Hapla 
7*27d0e99aSVaclav Hapla   Not Collective
8*27d0e99aSVaclav Hapla 
9*27d0e99aSVaclav Hapla   Input Parameters:
10*27d0e99aSVaclav Hapla + dm              - The DM (DMPLEX)
11*27d0e99aSVaclav Hapla . p               - The DAG point whose cone is compared
12*27d0e99aSVaclav Hapla . masterConeSize  - Number of the reference cone points passed (at least 2 and at most size of the cone of p)
13*27d0e99aSVaclav Hapla - masterCone      - The reference cone points
14*27d0e99aSVaclav Hapla 
15*27d0e99aSVaclav Hapla   Output Parameters:
16*27d0e99aSVaclav Hapla + start           - The new starting point within the cone of p to make it conforming with the reference cone
17*27d0e99aSVaclav Hapla - reverse         - The flag whether the order of the cone points should be reversed
18*27d0e99aSVaclav Hapla 
19*27d0e99aSVaclav Hapla   Level: advanced
20*27d0e99aSVaclav Hapla 
21*27d0e99aSVaclav Hapla .seealso: DMPlexOrient(), DMPlexOrientCell()
22*27d0e99aSVaclav Hapla @*/
23*27d0e99aSVaclav Hapla PetscErrorCode DMPlexCompareOrientations(DM dm, PetscInt p, PetscInt masterConeSize, const PetscInt masterCone[], PetscInt *start, PetscBool *reverse)
24*27d0e99aSVaclav Hapla {
25*27d0e99aSVaclav Hapla   PetscInt        coneSize;
26*27d0e99aSVaclav Hapla   const PetscInt *cone;
27*27d0e99aSVaclav Hapla   PetscInt        i, start_;
28*27d0e99aSVaclav Hapla   PetscBool       reverse_;
29*27d0e99aSVaclav Hapla   PetscErrorCode  ierr;
30*27d0e99aSVaclav Hapla 
31*27d0e99aSVaclav Hapla   PetscFunctionBegin;
32*27d0e99aSVaclav Hapla   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
33*27d0e99aSVaclav Hapla   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
34*27d0e99aSVaclav Hapla   if (coneSize < 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %D has no cone", p);
35*27d0e99aSVaclav Hapla   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
36*27d0e99aSVaclav Hapla   if (masterConeSize < 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %D: masterConeSize must be at least 2", p);
37*27d0e99aSVaclav Hapla   if (masterConeSize > coneSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %D: masterConeSize must be at most coneSize", p);
38*27d0e99aSVaclav Hapla   start_ = 0;
39*27d0e99aSVaclav Hapla   for (i=0; i<coneSize; i++) {
40*27d0e99aSVaclav Hapla     if (cone[i] == masterCone[0]) {
41*27d0e99aSVaclav Hapla       start_ = i;
42*27d0e99aSVaclav Hapla       break;
43*27d0e99aSVaclav Hapla     }
44*27d0e99aSVaclav Hapla   }
45*27d0e99aSVaclav Hapla   if (PetscUnlikely(i==coneSize)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Point %D: starting point of reference cone not found in slave cone", p);
46*27d0e99aSVaclav Hapla   reverse_ = PETSC_FALSE;
47*27d0e99aSVaclav Hapla   for (i=0; i<masterConeSize; i++) {if (cone[(start_+i)%coneSize] != masterCone[i]) break;}
48*27d0e99aSVaclav Hapla   if (i != masterConeSize) {
49*27d0e99aSVaclav Hapla     reverse_ = PETSC_TRUE;
50*27d0e99aSVaclav Hapla     for (i=0; i<masterConeSize; i++) {if (cone[(coneSize+start_-i)%coneSize] != masterCone[i]) break;}
51*27d0e99aSVaclav Hapla     if (i < masterConeSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Point %D: cone has non-conforming order of points with respect to reference cone", p);
52*27d0e99aSVaclav Hapla   }
53*27d0e99aSVaclav Hapla   if (start) *start = start_;
54*27d0e99aSVaclav Hapla   if (reverse) *reverse = reverse_;
55*27d0e99aSVaclav Hapla   if (PetscUnlikely(cone[start_] != masterCone[0])) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D: cone[%d] = %d != %d = masterCone[0]", p, start_, cone[start_], masterCone[0]);
56*27d0e99aSVaclav Hapla   PetscFunctionReturn(0);
57*27d0e99aSVaclav Hapla }
58*27d0e99aSVaclav Hapla 
59*27d0e99aSVaclav Hapla /*@
60*27d0e99aSVaclav Hapla   DMPlexOrientCell - Set the desired order of cone points of this DAG point, and fix orientations accordingly.
61*27d0e99aSVaclav Hapla 
62*27d0e99aSVaclav Hapla   Not Collective
63*27d0e99aSVaclav Hapla 
64*27d0e99aSVaclav Hapla   Input Parameters:
65*27d0e99aSVaclav Hapla + dm              - The DM
66*27d0e99aSVaclav Hapla . p               - The DAG point (from interval given by DMPlexGetChart())
67*27d0e99aSVaclav Hapla . masterConeSize  - Number of specified cone points (at least 2)
68*27d0e99aSVaclav Hapla - masterCone      - Specified cone points, i.e. ordered subset of current cone in DAG numbering (not cone-local numbering)
69*27d0e99aSVaclav Hapla 
70*27d0e99aSVaclav Hapla   Level: intermediate
71*27d0e99aSVaclav Hapla 
72*27d0e99aSVaclav Hapla .seealso: DMPlexOrient(), DMPlexGetCone(), DMPlexGetConeOrientation(), DMPlexInterpolate(), DMPlexGetChart()
73*27d0e99aSVaclav Hapla @*/
74*27d0e99aSVaclav Hapla PetscErrorCode DMPlexOrientCell(DM dm, PetscInt p, PetscInt masterConeSize, const PetscInt masterCone[])
75*27d0e99aSVaclav Hapla {
76*27d0e99aSVaclav Hapla   PetscInt        coneSize;
77*27d0e99aSVaclav Hapla   PetscInt        start1=0;
78*27d0e99aSVaclav Hapla   PetscBool       reverse1=PETSC_FALSE;
79*27d0e99aSVaclav Hapla   PetscErrorCode  ierr;
80*27d0e99aSVaclav Hapla 
81*27d0e99aSVaclav Hapla   PetscFunctionBegin;
82*27d0e99aSVaclav Hapla   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
83*27d0e99aSVaclav Hapla   if (masterConeSize) PetscValidIntPointer(masterCone,4);
84*27d0e99aSVaclav Hapla   if (masterConeSize == 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "masterConeSize cannot be 1");
85*27d0e99aSVaclav Hapla   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
86*27d0e99aSVaclav Hapla   if (!coneSize) PetscFunctionReturn(0); /* do nothing for points with no cone */
87*27d0e99aSVaclav Hapla   ierr = DMPlexCompareOrientations(dm, p, masterConeSize, masterCone, &start1, &reverse1);CHKERRQ(ierr);
88*27d0e99aSVaclav Hapla   ierr = DMPlexOrientCell_Internal(dm, p, start1, reverse1);CHKERRQ(ierr);
89*27d0e99aSVaclav Hapla #if defined(PETSC_USE_DEBUG)
90*27d0e99aSVaclav Hapla   {
91*27d0e99aSVaclav Hapla     PetscInt        c;
92*27d0e99aSVaclav Hapla     const PetscInt *cone;
93*27d0e99aSVaclav Hapla     ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
94*27d0e99aSVaclav Hapla     for (c = 0; c < masterConeSize; c++) {
95*27d0e99aSVaclav Hapla       if (PetscUnlikely(cone[c] != masterCone[c])) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_PLIB, "The algorithm above is wrong as cone[%d] = %d != %d = masterCone[%d]", c, cone[c], masterCone[c], c);
96*27d0e99aSVaclav Hapla     }
97*27d0e99aSVaclav Hapla   }
98*27d0e99aSVaclav Hapla #endif
99*27d0e99aSVaclav Hapla   PetscFunctionReturn(0);
100*27d0e99aSVaclav Hapla }
101*27d0e99aSVaclav Hapla 
102*27d0e99aSVaclav Hapla PetscErrorCode DMPlexOrientCell_Internal(DM dm, PetscInt p, PetscInt start1, PetscBool reverse1)
103*27d0e99aSVaclav Hapla {
104*27d0e99aSVaclav Hapla   PetscInt        i, j, k, maxConeSize, coneSize, coneConeSize, supportSize, supportConeSize;
105*27d0e99aSVaclav Hapla   PetscInt        start0, start;
106*27d0e99aSVaclav Hapla   PetscBool       reverse0, reverse;
107*27d0e99aSVaclav Hapla   PetscInt        newornt;
108*27d0e99aSVaclav Hapla   const PetscInt *cone=NULL, *support=NULL, *supportCone=NULL, *ornts=NULL;
109*27d0e99aSVaclav Hapla   PetscInt       *newcone=NULL, *newornts=NULL;
110*27d0e99aSVaclav Hapla   PetscErrorCode  ierr;
111*27d0e99aSVaclav Hapla 
112*27d0e99aSVaclav Hapla   PetscFunctionBegin;
113*27d0e99aSVaclav Hapla   if (!start1 && !reverse1) PetscFunctionReturn(0);
114*27d0e99aSVaclav Hapla   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
115*27d0e99aSVaclav Hapla   if (!coneSize) PetscFunctionReturn(0); /* do nothing for points with no cone */
116*27d0e99aSVaclav Hapla   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
117*27d0e99aSVaclav Hapla   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
118*27d0e99aSVaclav Hapla   /* permute p's cone and orientations */
119*27d0e99aSVaclav Hapla   ierr = DMPlexGetConeOrientation(dm, p, &ornts);CHKERRQ(ierr);
120*27d0e99aSVaclav Hapla   ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &newcone);CHKERRQ(ierr);
121*27d0e99aSVaclav Hapla   ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &newornts);CHKERRQ(ierr);
122*27d0e99aSVaclav Hapla   ierr = DMPlexFixFaceOrientations_Permute_Private(coneSize, cone, start1, reverse1, newcone);CHKERRQ(ierr);
123*27d0e99aSVaclav Hapla   ierr = DMPlexFixFaceOrientations_Permute_Private(coneSize, ornts, start1, reverse1, newornts);CHKERRQ(ierr);
124*27d0e99aSVaclav Hapla   /* if direction of p (face) is flipped, flip also p's cone points (edges) */
125*27d0e99aSVaclav Hapla   if (reverse1) {
126*27d0e99aSVaclav Hapla     for (i=0; i<coneSize; i++) {
127*27d0e99aSVaclav Hapla       ierr = DMPlexGetConeSize(dm, cone[i], &coneConeSize);CHKERRQ(ierr);
128*27d0e99aSVaclav Hapla       ierr = DMPlexFixFaceOrientations_Translate_Private(newornts[i], &start0, &reverse0);CHKERRQ(ierr);
129*27d0e99aSVaclav Hapla       ierr = DMPlexFixFaceOrientations_Combine_Private(coneConeSize, start0, reverse0, 1, PETSC_FALSE, &start, &reverse);CHKERRQ(ierr);
130*27d0e99aSVaclav Hapla       ierr = DMPlexFixFaceOrientations_TranslateBack_Private(coneConeSize, start, reverse, &newornts[i]);CHKERRQ(ierr);
131*27d0e99aSVaclav Hapla     }
132*27d0e99aSVaclav Hapla   }
133*27d0e99aSVaclav Hapla   ierr = DMPlexSetConeOrientation(dm, p, newornts);CHKERRQ(ierr);
134*27d0e99aSVaclav Hapla   /* fix oriention of p within cones of p's support points */
135*27d0e99aSVaclav Hapla   ierr = DMPlexGetSupport(dm, p, &support);CHKERRQ(ierr);
136*27d0e99aSVaclav Hapla   ierr = DMPlexGetSupportSize(dm, p, &supportSize);CHKERRQ(ierr);
137*27d0e99aSVaclav Hapla   for (j=0; j<supportSize; j++) {
138*27d0e99aSVaclav Hapla     ierr = DMPlexGetCone(dm, support[j], &supportCone);CHKERRQ(ierr);
139*27d0e99aSVaclav Hapla     ierr = DMPlexGetConeSize(dm, support[j], &supportConeSize);CHKERRQ(ierr);
140*27d0e99aSVaclav Hapla     ierr = DMPlexGetConeOrientation(dm, support[j], &ornts);CHKERRQ(ierr);
141*27d0e99aSVaclav Hapla     for (k=0; k<supportConeSize; k++) {
142*27d0e99aSVaclav Hapla       if (supportCone[k] != p) continue;
143*27d0e99aSVaclav Hapla       ierr = DMPlexFixFaceOrientations_Translate_Private(ornts[k], &start0, &reverse0);CHKERRQ(ierr);
144*27d0e99aSVaclav Hapla       ierr = DMPlexFixFaceOrientations_Combine_Private(coneSize, start0, reverse0, start1, reverse1, &start, &reverse);CHKERRQ(ierr);
145*27d0e99aSVaclav Hapla       ierr = DMPlexFixFaceOrientations_TranslateBack_Private(coneSize, start, reverse, &newornt);CHKERRQ(ierr);
146*27d0e99aSVaclav Hapla       ierr = DMPlexInsertConeOrientation(dm, support[j], k, newornt);CHKERRQ(ierr);
147*27d0e99aSVaclav Hapla     }
148*27d0e99aSVaclav Hapla   }
149*27d0e99aSVaclav Hapla   /* rewrite cone */
150*27d0e99aSVaclav Hapla   ierr = DMPlexSetCone(dm, p, newcone);CHKERRQ(ierr);
151*27d0e99aSVaclav Hapla   ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &newcone);CHKERRQ(ierr);
152*27d0e99aSVaclav Hapla   ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &newornts);CHKERRQ(ierr);
153*27d0e99aSVaclav Hapla   PetscFunctionReturn(0);
154*27d0e99aSVaclav Hapla }
155*27d0e99aSVaclav Hapla 
156*27d0e99aSVaclav Hapla /*@
15784961fc4SMatthew G. Knepley   DMPlexReverseCell - Give a mesh cell the opposite orientation
15884961fc4SMatthew G. Knepley 
15984961fc4SMatthew G. Knepley   Input Parameters:
16084961fc4SMatthew G. Knepley + dm   - The DM
16184961fc4SMatthew G. Knepley - cell - The cell number
16284961fc4SMatthew G. Knepley 
16384961fc4SMatthew G. Knepley   Note: The modification of the DM is done in-place.
16484961fc4SMatthew G. Knepley 
16584961fc4SMatthew G. Knepley   Level: advanced
16684961fc4SMatthew G. Knepley 
16784961fc4SMatthew G. Knepley .seealso: DMPlexOrient(), DMCreate(), DMPLEX
16884961fc4SMatthew G. Knepley @*/
16984961fc4SMatthew G. Knepley PetscErrorCode DMPlexReverseCell(DM dm, PetscInt cell)
17084961fc4SMatthew G. Knepley {
17184961fc4SMatthew G. Knepley   /* Note that the reverse orientation ro of a face with orientation o is:
17284961fc4SMatthew G. Knepley 
17384961fc4SMatthew G. Knepley        ro = o >= 0 ? -(faceSize - o) : faceSize + o
17484961fc4SMatthew G. Knepley 
17584961fc4SMatthew G. Knepley      where faceSize is the size of the cone for the face.
17684961fc4SMatthew G. Knepley   */
17784961fc4SMatthew G. Knepley   const PetscInt *cone,    *coneO, *support;
17884961fc4SMatthew G. Knepley   PetscInt       *revcone, *revconeO;
17984961fc4SMatthew G. Knepley   PetscInt        maxConeSize, coneSize, supportSize, faceSize, cp, sp;
18084961fc4SMatthew G. Knepley   PetscErrorCode  ierr;
18184961fc4SMatthew G. Knepley 
18284961fc4SMatthew G. Knepley   PetscFunctionBegin;
18384961fc4SMatthew G. Knepley   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
18469291d52SBarry Smith   ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &revcone);CHKERRQ(ierr);
18569291d52SBarry Smith   ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &revconeO);CHKERRQ(ierr);
18684961fc4SMatthew G. Knepley   /* Reverse cone, and reverse orientations of faces */
18784961fc4SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
18884961fc4SMatthew G. Knepley   ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
18984961fc4SMatthew G. Knepley   ierr = DMPlexGetConeOrientation(dm, cell, &coneO);CHKERRQ(ierr);
19084961fc4SMatthew G. Knepley   for (cp = 0; cp < coneSize; ++cp) {
19184961fc4SMatthew G. Knepley     const PetscInt rcp = coneSize-cp-1;
19284961fc4SMatthew G. Knepley 
19384961fc4SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, cone[rcp], &faceSize);CHKERRQ(ierr);
19484961fc4SMatthew G. Knepley     revcone[cp]  = cone[rcp];
19584961fc4SMatthew G. Knepley     revconeO[cp] = coneO[rcp] >= 0 ? -(faceSize-coneO[rcp]) : faceSize+coneO[rcp];
19684961fc4SMatthew G. Knepley   }
19784961fc4SMatthew G. Knepley   ierr = DMPlexSetCone(dm, cell, revcone);CHKERRQ(ierr);
19884961fc4SMatthew G. Knepley   ierr = DMPlexSetConeOrientation(dm, cell, revconeO);CHKERRQ(ierr);
19984961fc4SMatthew G. Knepley   /* Reverse orientation of this cell in the support hypercells */
20084961fc4SMatthew G. Knepley   faceSize = coneSize;
20184961fc4SMatthew G. Knepley   ierr = DMPlexGetSupportSize(dm, cell, &supportSize);CHKERRQ(ierr);
20284961fc4SMatthew G. Knepley   ierr = DMPlexGetSupport(dm, cell, &support);CHKERRQ(ierr);
20384961fc4SMatthew G. Knepley   for (sp = 0; sp < supportSize; ++sp) {
20484961fc4SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, support[sp], &coneSize);CHKERRQ(ierr);
20584961fc4SMatthew G. Knepley     ierr = DMPlexGetCone(dm, support[sp], &cone);CHKERRQ(ierr);
20684961fc4SMatthew G. Knepley     ierr = DMPlexGetConeOrientation(dm, support[sp], &coneO);CHKERRQ(ierr);
20784961fc4SMatthew G. Knepley     for (cp = 0; cp < coneSize; ++cp) {
20884961fc4SMatthew G. Knepley       if (cone[cp] != cell) continue;
20984961fc4SMatthew G. Knepley       ierr = DMPlexInsertConeOrientation(dm, support[sp], cp, coneO[cp] >= 0 ? -(faceSize-coneO[cp]) : faceSize+coneO[cp]);CHKERRQ(ierr);
21084961fc4SMatthew G. Knepley     }
21184961fc4SMatthew G. Knepley   }
21269291d52SBarry Smith   ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &revcone);CHKERRQ(ierr);
21369291d52SBarry Smith   ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &revconeO);CHKERRQ(ierr);
21484961fc4SMatthew G. Knepley   PetscFunctionReturn(0);
21584961fc4SMatthew G. Knepley }
21684961fc4SMatthew G. Knepley 
2177b310ce2SMatthew G. Knepley /*
2187b310ce2SMatthew G. Knepley   - Checks face match
2197b310ce2SMatthew G. Knepley     - Flips non-matching
2207b310ce2SMatthew G. Knepley   - Inserts faces of support cells in FIFO
2217b310ce2SMatthew G. Knepley */
2227b2de0fdSMatthew G. Knepley static PetscErrorCode DMPlexCheckFace_Internal(DM dm, PetscInt *faceFIFO, PetscInt *fTop, PetscInt *fBottom, PetscInt cStart, PetscInt fStart, PetscInt fEnd, PetscBT seenCells, PetscBT flippedCells, PetscBT seenFaces)
2237b310ce2SMatthew G. Knepley {
2247b310ce2SMatthew G. Knepley   const PetscInt *support, *coneA, *coneB, *coneOA, *coneOB;
2257b310ce2SMatthew G. Knepley   PetscInt        supportSize, coneSizeA, coneSizeB, posA = -1, posB = -1;
2267b2de0fdSMatthew G. Knepley   PetscInt        face, dim, seenA, flippedA, seenB, flippedB, mismatch, c;
2277b310ce2SMatthew G. Knepley   PetscErrorCode  ierr;
2287b310ce2SMatthew G. Knepley 
2297b310ce2SMatthew G. Knepley   PetscFunctionBegin;
2307b310ce2SMatthew G. Knepley   face = faceFIFO[(*fTop)++];
2317b310ce2SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2327b310ce2SMatthew G. Knepley   ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
2337b310ce2SMatthew G. Knepley   ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr);
2347b310ce2SMatthew G. Knepley   if (supportSize < 2) PetscFunctionReturn(0);
2357b310ce2SMatthew G. Knepley   if (supportSize != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %d", supportSize);
2367b310ce2SMatthew G. Knepley   seenA    = PetscBTLookup(seenCells,    support[0]-cStart);
2377b310ce2SMatthew G. Knepley   flippedA = PetscBTLookup(flippedCells, support[0]-cStart) ? 1 : 0;
2387b310ce2SMatthew G. Knepley   seenB    = PetscBTLookup(seenCells,    support[1]-cStart);
2397b310ce2SMatthew G. Knepley   flippedB = PetscBTLookup(flippedCells, support[1]-cStart) ? 1 : 0;
2407b310ce2SMatthew G. Knepley 
2417b310ce2SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, support[0], &coneSizeA);CHKERRQ(ierr);
2427b310ce2SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, support[1], &coneSizeB);CHKERRQ(ierr);
2437b310ce2SMatthew G. Knepley   ierr = DMPlexGetCone(dm, support[0], &coneA);CHKERRQ(ierr);
2447b310ce2SMatthew G. Knepley   ierr = DMPlexGetCone(dm, support[1], &coneB);CHKERRQ(ierr);
2457b310ce2SMatthew G. Knepley   ierr = DMPlexGetConeOrientation(dm, support[0], &coneOA);CHKERRQ(ierr);
2467b310ce2SMatthew G. Knepley   ierr = DMPlexGetConeOrientation(dm, support[1], &coneOB);CHKERRQ(ierr);
2477b310ce2SMatthew G. Knepley   for (c = 0; c < coneSizeA; ++c) {
2487b310ce2SMatthew G. Knepley     if (!PetscBTLookup(seenFaces, coneA[c]-fStart)) {
2497b310ce2SMatthew G. Knepley       faceFIFO[(*fBottom)++] = coneA[c];
2507b310ce2SMatthew G. Knepley       ierr = PetscBTSet(seenFaces, coneA[c]-fStart);CHKERRQ(ierr);
2517b310ce2SMatthew G. Knepley     }
2527b310ce2SMatthew G. Knepley     if (coneA[c] == face) posA = c;
2537b310ce2SMatthew G. Knepley     if (*fBottom > fEnd-fStart) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %d was pushed exceeding capacity %d > %d", coneA[c], *fBottom, fEnd-fStart);
2547b310ce2SMatthew G. Knepley   }
2557b310ce2SMatthew G. Knepley   if (posA < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d could not be located in cell %d", face, support[0]);
2567b310ce2SMatthew G. Knepley   for (c = 0; c < coneSizeB; ++c) {
2577b310ce2SMatthew G. Knepley     if (!PetscBTLookup(seenFaces, coneB[c]-fStart)) {
2587b310ce2SMatthew G. Knepley       faceFIFO[(*fBottom)++] = coneB[c];
2597b310ce2SMatthew G. Knepley       ierr = PetscBTSet(seenFaces, coneB[c]-fStart);CHKERRQ(ierr);
2607b310ce2SMatthew G. Knepley     }
2617b310ce2SMatthew G. Knepley     if (coneB[c] == face) posB = c;
2627b310ce2SMatthew G. Knepley     if (*fBottom > fEnd-fStart) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %d was pushed exceeding capacity %d > %d", coneA[c], *fBottom, fEnd-fStart);
2637b310ce2SMatthew G. Knepley   }
2647b310ce2SMatthew G. Knepley   if (posB < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d could not be located in cell %d", face, support[1]);
2657b310ce2SMatthew G. Knepley 
2667b310ce2SMatthew G. Knepley   if (dim == 1) {
2677b310ce2SMatthew G. Knepley     mismatch = posA == posB;
2687b310ce2SMatthew G. Knepley   } else {
2697b310ce2SMatthew G. Knepley     mismatch = coneOA[posA] == coneOB[posB];
2707b310ce2SMatthew G. Knepley   }
2717b310ce2SMatthew G. Knepley 
2727b310ce2SMatthew G. Knepley   if (mismatch ^ (flippedA ^ flippedB)) {
2737b310ce2SMatthew G. Knepley     if (seenA && seenB) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Previously seen cells %d and %d do not match: Fault mesh is non-orientable", support[0], support[1]);
2747b310ce2SMatthew G. Knepley     if (!seenA && !flippedA) {
2757b310ce2SMatthew G. Knepley       ierr = PetscBTSet(flippedCells, support[0]-cStart);CHKERRQ(ierr);
2767b310ce2SMatthew G. Knepley     } else if (!seenB && !flippedB) {
2777b310ce2SMatthew G. Knepley       ierr = PetscBTSet(flippedCells, support[1]-cStart);CHKERRQ(ierr);
2787b310ce2SMatthew G. Knepley     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
2797b310ce2SMatthew G. Knepley   } else if (mismatch && flippedA && flippedB) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
2807b310ce2SMatthew G. Knepley   ierr = PetscBTSet(seenCells, support[0]-cStart);CHKERRQ(ierr);
2817b310ce2SMatthew G. Knepley   ierr = PetscBTSet(seenCells, support[1]-cStart);CHKERRQ(ierr);
2827b310ce2SMatthew G. Knepley   PetscFunctionReturn(0);
2837b310ce2SMatthew G. Knepley }
2847b310ce2SMatthew G. Knepley 
28584961fc4SMatthew G. Knepley /*@
28684961fc4SMatthew G. Knepley   DMPlexOrient - Give a consistent orientation to the input mesh
28784961fc4SMatthew G. Knepley 
28884961fc4SMatthew G. Knepley   Input Parameters:
28984961fc4SMatthew G. Knepley . dm - The DM
29084961fc4SMatthew G. Knepley 
29184961fc4SMatthew G. Knepley   Note: The orientation data for the DM are change in-place.
29284961fc4SMatthew G. Knepley $ This routine will fail for non-orientable surfaces, such as the Moebius strip.
29384961fc4SMatthew G. Knepley 
29484961fc4SMatthew G. Knepley   Level: advanced
29584961fc4SMatthew G. Knepley 
29684961fc4SMatthew G. Knepley .seealso: DMCreate(), DMPLEX
29784961fc4SMatthew G. Knepley @*/
29884961fc4SMatthew G. Knepley PetscErrorCode DMPlexOrient(DM dm)
29984961fc4SMatthew G. Knepley {
30084961fc4SMatthew G. Knepley   MPI_Comm           comm;
301e1d83109SMatthew G. Knepley   PetscSF            sf;
302e1d83109SMatthew G. Knepley   const PetscInt    *lpoints;
303e1d83109SMatthew G. Knepley   const PetscSFNode *rpoints;
304e1d83109SMatthew G. Knepley   PetscSFNode       *rorntComp = NULL, *lorntComp = NULL;
305e1d83109SMatthew G. Knepley   PetscInt          *numNeighbors, **neighbors;
306e1d83109SMatthew G. Knepley   PetscSFNode       *nrankComp;
307e1d83109SMatthew G. Knepley   PetscBool         *match, *flipped;
30884961fc4SMatthew G. Knepley   PetscBT            seenCells, flippedCells, seenFaces;
309e1d83109SMatthew G. Knepley   PetscInt          *faceFIFO, fTop, fBottom, *cellComp, *faceComp;
3107cadcfe8SMatthew G. Knepley   PetscInt           numLeaves, numRoots, dim, h, cStart, cEnd, c, cell, fStart, fEnd, face, off, totNeighbors = 0;
31139526728SToby Isaac   PetscMPIInt        rank, size, numComponents, comp = 0;
31239526728SToby Isaac   PetscBool          flg, flg2;
313a17aa47eSToby Isaac   PetscViewer        viewer = NULL, selfviewer = NULL;
31484961fc4SMatthew G. Knepley   PetscErrorCode     ierr;
31584961fc4SMatthew G. Knepley 
31684961fc4SMatthew G. Knepley   PetscFunctionBegin;
31784961fc4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
318a83982f3SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
31939526728SToby Isaac   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
320c5929fdfSBarry Smith   ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view", &flg);CHKERRQ(ierr);
32139526728SToby Isaac   ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view_synchronized", &flg2);CHKERRQ(ierr);
322e1d83109SMatthew G. Knepley   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
323e1d83109SMatthew G. Knepley   ierr = PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints);CHKERRQ(ierr);
32484961fc4SMatthew G. Knepley   /* Truth Table
32584961fc4SMatthew G. Knepley      mismatch    flips   do action   mismatch   flipA ^ flipB   action
32684961fc4SMatthew G. Knepley          F       0 flips     no         F             F           F
32784961fc4SMatthew G. Knepley          F       1 flip      yes        F             T           T
32884961fc4SMatthew G. Knepley          F       2 flips     no         T             F           T
32984961fc4SMatthew G. Knepley          T       0 flips     yes        T             T           F
33084961fc4SMatthew G. Knepley          T       1 flip      no
33184961fc4SMatthew G. Knepley          T       2 flips     yes
33284961fc4SMatthew G. Knepley   */
33384961fc4SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
33484961fc4SMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(dm, &h);CHKERRQ(ierr);
33584961fc4SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, h,   &cStart, &cEnd);CHKERRQ(ierr);
33684961fc4SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, h+1, &fStart, &fEnd);CHKERRQ(ierr);
33784961fc4SMatthew G. Knepley   ierr = PetscBTCreate(cEnd - cStart, &seenCells);CHKERRQ(ierr);
33884961fc4SMatthew G. Knepley   ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr);
33984961fc4SMatthew G. Knepley   ierr = PetscBTCreate(cEnd - cStart, &flippedCells);CHKERRQ(ierr);
34084961fc4SMatthew G. Knepley   ierr = PetscBTMemzero(cEnd - cStart, flippedCells);CHKERRQ(ierr);
34184961fc4SMatthew G. Knepley   ierr = PetscBTCreate(fEnd - fStart, &seenFaces);CHKERRQ(ierr);
34284961fc4SMatthew G. Knepley   ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr);
343e1d83109SMatthew G. Knepley   ierr = PetscCalloc3(fEnd - fStart, &faceFIFO, cEnd-cStart, &cellComp, fEnd-fStart, &faceComp);CHKERRQ(ierr);
344e1d83109SMatthew G. Knepley   /*
345e1d83109SMatthew G. Knepley    OLD STYLE
346e1d83109SMatthew G. Knepley    - Add an integer array over cells and faces (component) for connected component number
347e1d83109SMatthew G. Knepley    Foreach component
348e1d83109SMatthew G. Knepley      - Mark the initial cell as seen
349e1d83109SMatthew G. Knepley      - Process component as usual
350e1d83109SMatthew G. Knepley      - Set component for all seenCells
351e1d83109SMatthew G. Knepley      - Wipe seenCells and seenFaces (flippedCells can stay)
352e1d83109SMatthew G. Knepley    - Generate parallel adjacency for component using SF and seenFaces
353e1d83109SMatthew G. Knepley    - Collect numComponents adj data from each proc to 0
354e1d83109SMatthew G. Knepley    - Build same serial graph
355e1d83109SMatthew G. Knepley    - Use same solver
356e1d83109SMatthew G. Knepley    - Use Scatterv to to send back flipped flags for each component
357e1d83109SMatthew G. Knepley    - Negate flippedCells by component
358e1d83109SMatthew G. Knepley 
359e1d83109SMatthew G. Knepley    NEW STYLE
360e1d83109SMatthew G. Knepley    - Create the adj on each process
361e1d83109SMatthew G. Knepley    - Bootstrap to complete graph on proc 0
362e1d83109SMatthew G. Knepley   */
363e1d83109SMatthew G. Knepley   /* Loop over components */
364e1d83109SMatthew G. Knepley   for (cell = cStart; cell < cEnd; ++cell) cellComp[cell-cStart] = -1;
365e1d83109SMatthew G. Knepley   do {
366e1d83109SMatthew G. Knepley     /* Look for first unmarked cell */
367e1d83109SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) if (cellComp[cell-cStart] < 0) break;
368e1d83109SMatthew G. Knepley     if (cell >= cEnd) break;
369e1d83109SMatthew G. Knepley     /* Initialize FIFO with first cell in component */
370e1d83109SMatthew G. Knepley     {
37184961fc4SMatthew G. Knepley       const PetscInt *cone;
37284961fc4SMatthew G. Knepley       PetscInt        coneSize;
37384961fc4SMatthew G. Knepley 
374e1d83109SMatthew G. Knepley       fTop = fBottom = 0;
375e1d83109SMatthew G. Knepley       ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
376e1d83109SMatthew G. Knepley       ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
37784961fc4SMatthew G. Knepley       for (c = 0; c < coneSize; ++c) {
37884961fc4SMatthew G. Knepley         faceFIFO[fBottom++] = cone[c];
37984961fc4SMatthew G. Knepley         ierr = PetscBTSet(seenFaces, cone[c]-fStart);CHKERRQ(ierr);
38084961fc4SMatthew G. Knepley       }
38131c8331aSMatthew G. Knepley       ierr = PetscBTSet(seenCells, cell-cStart);CHKERRQ(ierr);
38284961fc4SMatthew G. Knepley     }
38384961fc4SMatthew G. Knepley     /* Consider each face in FIFO */
38484961fc4SMatthew G. Knepley     while (fTop < fBottom) {
3857b2de0fdSMatthew G. Knepley       ierr = DMPlexCheckFace_Internal(dm, faceFIFO, &fTop, &fBottom, cStart, fStart, fEnd, seenCells, flippedCells, seenFaces);CHKERRQ(ierr);
38684961fc4SMatthew G. Knepley     }
387e1d83109SMatthew G. Knepley     /* Set component for cells and faces */
388e1d83109SMatthew G. Knepley     for (cell = 0; cell < cEnd-cStart; ++cell) {
389e1d83109SMatthew G. Knepley       if (PetscBTLookup(seenCells, cell)) cellComp[cell] = comp;
390e1d83109SMatthew G. Knepley     }
391e1d83109SMatthew G. Knepley     for (face = 0; face < fEnd-fStart; ++face) {
392e1d83109SMatthew G. Knepley       if (PetscBTLookup(seenFaces, face)) faceComp[face] = comp;
393e1d83109SMatthew G. Knepley     }
394e1d83109SMatthew G. Knepley     /* Wipe seenCells and seenFaces for next component */
395e1d83109SMatthew G. Knepley     ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr);
396e1d83109SMatthew G. Knepley     ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr);
397e1d83109SMatthew G. Knepley     ++comp;
398e1d83109SMatthew G. Knepley   } while (1);
399e1d83109SMatthew G. Knepley   numComponents = comp;
40084961fc4SMatthew G. Knepley   if (flg) {
40184961fc4SMatthew G. Knepley     PetscViewer v;
40284961fc4SMatthew G. Knepley 
4037e09831bSMatthew G. Knepley     ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr);
4041575c14dSBarry Smith     ierr = PetscViewerASCIIPushSynchronized(v);CHKERRQ(ierr);
40584961fc4SMatthew G. Knepley     ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank);CHKERRQ(ierr);
40684961fc4SMatthew G. Knepley     ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr);
4074d4c343aSBarry Smith     ierr = PetscViewerFlush(v);CHKERRQ(ierr);
4081575c14dSBarry Smith     ierr = PetscViewerASCIIPopSynchronized(v);CHKERRQ(ierr);
40984961fc4SMatthew G. Knepley   }
41084961fc4SMatthew G. Knepley   /* Now all subdomains are oriented, but we need a consistent parallel orientation */
41184961fc4SMatthew G. Knepley   if (numLeaves >= 0) {
412e1d83109SMatthew G. Knepley     /* Store orientations of boundary faces*/
413e1d83109SMatthew G. Knepley     ierr = PetscCalloc2(numRoots,&rorntComp,numRoots,&lorntComp);CHKERRQ(ierr);
41484961fc4SMatthew G. Knepley     for (face = fStart; face < fEnd; ++face) {
415e1d83109SMatthew G. Knepley       const PetscInt *cone, *support, *ornt;
416e1d83109SMatthew G. Knepley       PetscInt        coneSize, supportSize;
417e1d83109SMatthew G. Knepley 
41884961fc4SMatthew G. Knepley       ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
41984961fc4SMatthew G. Knepley       if (supportSize != 1) continue;
42084961fc4SMatthew G. Knepley       ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr);
42184961fc4SMatthew G. Knepley 
42284961fc4SMatthew G. Knepley       ierr = DMPlexGetCone(dm, support[0], &cone);CHKERRQ(ierr);
42384961fc4SMatthew G. Knepley       ierr = DMPlexGetConeSize(dm, support[0], &coneSize);CHKERRQ(ierr);
42484961fc4SMatthew G. Knepley       ierr = DMPlexGetConeOrientation(dm, support[0], &ornt);CHKERRQ(ierr);
42584961fc4SMatthew G. Knepley       for (c = 0; c < coneSize; ++c) if (cone[c] == face) break;
42684961fc4SMatthew G. Knepley       if (dim == 1) {
42784961fc4SMatthew G. Knepley         /* Use cone position instead, shifted to -1 or 1 */
4286e1eb25cSMatthew G. Knepley         if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = 1-c*2;
4296e1eb25cSMatthew G. Knepley         else                                                rorntComp[face].rank = c*2-1;
43084961fc4SMatthew G. Knepley       } else {
431e1d83109SMatthew G. Knepley         if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 :  1;
432e1d83109SMatthew G. Knepley         else                                                rorntComp[face].rank = ornt[c] < 0 ?  1 : -1;
43384961fc4SMatthew G. Knepley       }
434e1d83109SMatthew G. Knepley       rorntComp[face].index = faceComp[face-fStart];
43584961fc4SMatthew G. Knepley     }
436e1d83109SMatthew G. Knepley     /* Communicate boundary edge orientations */
437e1d83109SMatthew G. Knepley     ierr = PetscSFBcastBegin(sf, MPIU_2INT, rorntComp, lorntComp);CHKERRQ(ierr);
438e1d83109SMatthew G. Knepley     ierr = PetscSFBcastEnd(sf, MPIU_2INT, rorntComp, lorntComp);CHKERRQ(ierr);
439e1d83109SMatthew G. Knepley   }
440e1d83109SMatthew G. Knepley   /* Get process adjacency */
441e1d83109SMatthew G. Knepley   ierr = PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors);CHKERRQ(ierr);
442a17aa47eSToby Isaac   viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm));
443a17aa47eSToby Isaac   if (flg2) {ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);}
444a17aa47eSToby Isaac   ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&selfviewer);CHKERRQ(ierr);
445e1d83109SMatthew G. Knepley   for (comp = 0; comp < numComponents; ++comp) {
446e1d83109SMatthew G. Knepley     PetscInt  l, n;
44784961fc4SMatthew G. Knepley 
448e1d83109SMatthew G. Knepley     numNeighbors[comp] = 0;
44939ea070eSMatthew G. Knepley     ierr = PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp]);CHKERRQ(ierr);
450e1d83109SMatthew G. Knepley     /* I know this is p^2 time in general, but for bounded degree its alright */
451e1d83109SMatthew G. Knepley     for (l = 0; l < numLeaves; ++l) {
452e1d83109SMatthew G. Knepley       const PetscInt face = lpoints[l];
453e1d83109SMatthew G. Knepley 
454e1d83109SMatthew G. Knepley       /* Find a representative face (edge) separating pairs of procs */
455e1d83109SMatthew G. Knepley       if ((face >= fStart) && (face < fEnd) && (faceComp[face-fStart] == comp)) {
456269a49d5SMatthew G. Knepley         const PetscInt rrank = rpoints[l].rank;
457e1d83109SMatthew G. Knepley         const PetscInt rcomp = lorntComp[face].index;
458e1d83109SMatthew G. Knepley 
459269a49d5SMatthew G. Knepley         for (n = 0; n < numNeighbors[comp]; ++n) if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break;
460e1d83109SMatthew G. Knepley         if (n >= numNeighbors[comp]) {
461e1d83109SMatthew G. Knepley           PetscInt supportSize;
462e1d83109SMatthew G. Knepley 
463e1d83109SMatthew G. Knepley           ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
464e1d83109SMatthew G. Knepley           if (supportSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %d", supportSize);
46535041eceSToby Isaac           if (flg) {ierr = PetscViewerASCIIPrintf(selfviewer, "[%d]: component %d, Found representative leaf %d (face %d) connecting to face %d on (%d, %d) with orientation %d\n", rank, comp, l, face, rpoints[l].index, rrank, rcomp, lorntComp[face].rank);CHKERRQ(ierr);}
466e1d83109SMatthew G. Knepley           neighbors[comp][numNeighbors[comp]++] = l;
467e1d83109SMatthew G. Knepley         }
468e1d83109SMatthew G. Knepley       }
469e1d83109SMatthew G. Knepley     }
470e1d83109SMatthew G. Knepley     totNeighbors += numNeighbors[comp];
471e1d83109SMatthew G. Knepley   }
472a17aa47eSToby Isaac   ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&selfviewer);CHKERRQ(ierr);
473a17aa47eSToby Isaac   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
474a17aa47eSToby Isaac   if (flg2) {ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);}
475e1d83109SMatthew G. Knepley   ierr = PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match);CHKERRQ(ierr);
476e1d83109SMatthew G. Knepley   for (comp = 0, off = 0; comp < numComponents; ++comp) {
477e1d83109SMatthew G. Knepley     PetscInt n;
478e1d83109SMatthew G. Knepley 
479e1d83109SMatthew G. Knepley     for (n = 0; n < numNeighbors[comp]; ++n, ++off) {
480e1d83109SMatthew G. Knepley       const PetscInt face = lpoints[neighbors[comp][n]];
481e1d83109SMatthew G. Knepley       const PetscInt o    = rorntComp[face].rank*lorntComp[face].rank;
482e1d83109SMatthew G. Knepley 
483e1d83109SMatthew G. Knepley       if      (o < 0) match[off] = PETSC_TRUE;
484e1d83109SMatthew G. Knepley       else if (o > 0) match[off] = PETSC_FALSE;
485e1d83109SMatthew G. Knepley       else SETERRQ5(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid face %d (%d, %d) neighbor: %d comp: %d", face, rorntComp[face], lorntComp[face], neighbors[comp][n], comp);
486e1d83109SMatthew G. Knepley       nrankComp[off].rank  = rpoints[neighbors[comp][n]].rank;
487e1d83109SMatthew G. Knepley       nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index;
488e1d83109SMatthew G. Knepley     }
489e1d83109SMatthew G. Knepley     ierr = PetscFree(neighbors[comp]);CHKERRQ(ierr);
49084961fc4SMatthew G. Knepley   }
49184961fc4SMatthew G. Knepley   /* Collect the graph on 0 */
492e1d83109SMatthew G. Knepley   if (numLeaves >= 0) {
49384961fc4SMatthew G. Knepley     Mat          G;
49484961fc4SMatthew G. Knepley     PetscBT      seenProcs, flippedProcs;
49584961fc4SMatthew G. Knepley     PetscInt    *procFIFO, pTop, pBottom;
4967cadcfe8SMatthew G. Knepley     PetscInt    *N   = NULL, *Noff;
497e1d83109SMatthew G. Knepley     PetscSFNode *adj = NULL;
49884961fc4SMatthew G. Knepley     PetscBool   *val = NULL;
4997cadcfe8SMatthew G. Knepley     PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc, p, o;
5009852e123SBarry Smith     PetscMPIInt  size = 0;
50184961fc4SMatthew G. Knepley 
502e1d83109SMatthew G. Knepley     ierr = PetscCalloc1(numComponents, &flipped);CHKERRQ(ierr);
5039852e123SBarry Smith     if (!rank) {ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);}
5049852e123SBarry Smith     ierr = PetscCalloc4(size, &recvcounts, size+1, &displs, size, &Nc, size+1, &Noff);CHKERRQ(ierr);
505e1d83109SMatthew G. Knepley     ierr = MPI_Gather(&numComponents, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm);CHKERRQ(ierr);
5069852e123SBarry Smith     for (p = 0; p < size; ++p) {
507e1d83109SMatthew G. Knepley       displs[p+1] = displs[p] + Nc[p];
508e1d83109SMatthew G. Knepley     }
5099852e123SBarry Smith     if (!rank) {ierr = PetscMalloc1(displs[size],&N);CHKERRQ(ierr);}
5107cadcfe8SMatthew G. Knepley     ierr = MPI_Gatherv(numNeighbors, numComponents, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm);CHKERRQ(ierr);
5119852e123SBarry Smith     for (p = 0, o = 0; p < size; ++p) {
512e1d83109SMatthew G. Knepley       recvcounts[p] = 0;
513e1d83109SMatthew G. Knepley       for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o];
51484961fc4SMatthew G. Knepley       displs[p+1] = displs[p] + recvcounts[p];
51584961fc4SMatthew G. Knepley     }
5169852e123SBarry Smith     if (!rank) {ierr = PetscMalloc2(displs[size], &adj, displs[size], &val);CHKERRQ(ierr);}
517e1d83109SMatthew G. Knepley     ierr = MPI_Gatherv(nrankComp, totNeighbors, MPIU_2INT, adj, recvcounts, displs, MPIU_2INT, 0, comm);CHKERRQ(ierr);
518e1d83109SMatthew G. Knepley     ierr = MPI_Gatherv(match, totNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
519e1d83109SMatthew G. Knepley     ierr = PetscFree2(numNeighbors, neighbors);CHKERRQ(ierr);
520e1d83109SMatthew G. Knepley     if (!rank) {
5219852e123SBarry Smith       for (p = 1; p <= size; ++p) {Noff[p] = Noff[p-1] + Nc[p-1];}
522a83982f3SMatthew G. Knepley       if (flg) {
523e1d83109SMatthew G. Knepley         PetscInt n;
524e1d83109SMatthew G. Knepley 
5259852e123SBarry Smith         for (p = 0, off = 0; p < size; ++p) {
526e1d83109SMatthew G. Knepley           for (c = 0; c < Nc[p]; ++c) {
527302440fdSBarry Smith             ierr = PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %d:\n", p, c);CHKERRQ(ierr);
528e1d83109SMatthew G. Knepley             for (n = 0; n < N[Noff[p]+c]; ++n, ++off) {
529302440fdSBarry Smith               ierr = PetscPrintf(PETSC_COMM_SELF, "  edge (%d, %d) (%d):\n", adj[off].rank, adj[off].index, val[off]);CHKERRQ(ierr);
530e1d83109SMatthew G. Knepley             }
53184961fc4SMatthew G. Knepley           }
53284961fc4SMatthew G. Knepley         }
53384961fc4SMatthew G. Knepley       }
53484961fc4SMatthew G. Knepley       /* Symmetrize the graph */
53584961fc4SMatthew G. Knepley       ierr = MatCreate(PETSC_COMM_SELF, &G);CHKERRQ(ierr);
5369852e123SBarry Smith       ierr = MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size]);CHKERRQ(ierr);
53784961fc4SMatthew G. Knepley       ierr = MatSetUp(G);CHKERRQ(ierr);
5389852e123SBarry Smith       for (p = 0, off = 0; p < size; ++p) {
539e1d83109SMatthew G. Knepley         for (c = 0; c < Nc[p]; ++c) {
540e1d83109SMatthew G. Knepley           const PetscInt r = Noff[p]+c;
541e1d83109SMatthew G. Knepley           PetscInt       n;
542e1d83109SMatthew G. Knepley 
543e1d83109SMatthew G. Knepley           for (n = 0; n < N[r]; ++n, ++off) {
544e1d83109SMatthew G. Knepley             const PetscInt    q = Noff[adj[off].rank] + adj[off].index;
545e1d83109SMatthew G. Knepley             const PetscScalar o = val[off] ? 1.0 : 0.0;
54684961fc4SMatthew G. Knepley 
54784961fc4SMatthew G. Knepley             ierr = MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES);CHKERRQ(ierr);
54884961fc4SMatthew G. Knepley             ierr = MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES);CHKERRQ(ierr);
54984961fc4SMatthew G. Knepley           }
55084961fc4SMatthew G. Knepley         }
551e1d83109SMatthew G. Knepley       }
55284961fc4SMatthew G. Knepley       ierr = MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
55384961fc4SMatthew G. Knepley       ierr = MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
55484961fc4SMatthew G. Knepley 
5559852e123SBarry Smith       ierr = PetscBTCreate(Noff[size], &seenProcs);CHKERRQ(ierr);
5569852e123SBarry Smith       ierr = PetscBTMemzero(Noff[size], seenProcs);CHKERRQ(ierr);
5579852e123SBarry Smith       ierr = PetscBTCreate(Noff[size], &flippedProcs);CHKERRQ(ierr);
5589852e123SBarry Smith       ierr = PetscBTMemzero(Noff[size], flippedProcs);CHKERRQ(ierr);
5599852e123SBarry Smith       ierr = PetscMalloc1(Noff[size], &procFIFO);CHKERRQ(ierr);
56084961fc4SMatthew G. Knepley       pTop = pBottom = 0;
5619852e123SBarry Smith       for (p = 0; p < Noff[size]; ++p) {
56284961fc4SMatthew G. Knepley         if (PetscBTLookup(seenProcs, p)) continue;
56384961fc4SMatthew G. Knepley         /* Initialize FIFO with next proc */
56484961fc4SMatthew G. Knepley         procFIFO[pBottom++] = p;
56584961fc4SMatthew G. Knepley         ierr = PetscBTSet(seenProcs, p);CHKERRQ(ierr);
56684961fc4SMatthew G. Knepley         /* Consider each proc in FIFO */
56784961fc4SMatthew G. Knepley         while (pTop < pBottom) {
56884961fc4SMatthew G. Knepley           const PetscScalar *ornt;
56984961fc4SMatthew G. Knepley           const PetscInt    *neighbors;
570e1d83109SMatthew G. Knepley           PetscInt           proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n;
57184961fc4SMatthew G. Knepley 
57284961fc4SMatthew G. Knepley           proc     = procFIFO[pTop++];
57384961fc4SMatthew G. Knepley           flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0;
57484961fc4SMatthew G. Knepley           ierr = MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt);CHKERRQ(ierr);
57584961fc4SMatthew G. Knepley           /* Loop over neighboring procs */
57684961fc4SMatthew G. Knepley           for (n = 0; n < numNeighbors; ++n) {
57784961fc4SMatthew G. Knepley             nproc    = neighbors[n];
57884961fc4SMatthew G. Knepley             mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1;
57984961fc4SMatthew G. Knepley             seen     = PetscBTLookup(seenProcs, nproc);
58084961fc4SMatthew G. Knepley             flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0;
58184961fc4SMatthew G. Knepley 
58284961fc4SMatthew G. Knepley             if (mismatch ^ (flippedA ^ flippedB)) {
58384961fc4SMatthew G. Knepley               if (seen) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Previously seen procs %d and %d do not match: Fault mesh is non-orientable", proc, nproc);
58484961fc4SMatthew G. Knepley               if (!flippedB) {
58584961fc4SMatthew G. Knepley                 ierr = PetscBTSet(flippedProcs, nproc);CHKERRQ(ierr);
58684961fc4SMatthew G. Knepley               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
58784961fc4SMatthew G. Knepley             } else if (mismatch && flippedA && flippedB) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
58884961fc4SMatthew G. Knepley             if (!seen) {
58984961fc4SMatthew G. Knepley               procFIFO[pBottom++] = nproc;
59084961fc4SMatthew G. Knepley               ierr = PetscBTSet(seenProcs, nproc);CHKERRQ(ierr);
59184961fc4SMatthew G. Knepley             }
59284961fc4SMatthew G. Knepley           }
59384961fc4SMatthew G. Knepley         }
59484961fc4SMatthew G. Knepley       }
59584961fc4SMatthew G. Knepley       ierr = PetscFree(procFIFO);CHKERRQ(ierr);
59684961fc4SMatthew G. Knepley       ierr = MatDestroy(&G);CHKERRQ(ierr);
59784961fc4SMatthew G. Knepley       ierr = PetscFree2(adj, val);CHKERRQ(ierr);
598e1d83109SMatthew G. Knepley       ierr = PetscBTDestroy(&seenProcs);CHKERRQ(ierr);
59984961fc4SMatthew G. Knepley     }
600e1d83109SMatthew G. Knepley     /* Scatter flip flags */
601e1d83109SMatthew G. Knepley     {
602e1d83109SMatthew G. Knepley       PetscBool *flips = NULL;
603e1d83109SMatthew G. Knepley 
604e1d83109SMatthew G. Knepley       if (!rank) {
6059852e123SBarry Smith         ierr = PetscMalloc1(Noff[size], &flips);CHKERRQ(ierr);
6069852e123SBarry Smith         for (p = 0; p < Noff[size]; ++p) {
607e1d83109SMatthew G. Knepley           flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE;
608302440fdSBarry Smith           if (flg && flips[p]) {ierr = PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p);CHKERRQ(ierr);}
609e1d83109SMatthew G. Knepley         }
6109852e123SBarry Smith         for (p = 0; p < size; ++p) {
611e1d83109SMatthew G. Knepley           displs[p+1] = displs[p] + Nc[p];
612e1d83109SMatthew G. Knepley         }
613e1d83109SMatthew G. Knepley       }
614e1d83109SMatthew G. Knepley       ierr = MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, numComponents, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
61584961fc4SMatthew G. Knepley       ierr = PetscFree(flips);CHKERRQ(ierr);
61684961fc4SMatthew G. Knepley     }
617e1d83109SMatthew G. Knepley     if (!rank) {ierr = PetscBTDestroy(&flippedProcs);CHKERRQ(ierr);}
618e1d83109SMatthew G. Knepley     ierr = PetscFree(N);CHKERRQ(ierr);
619e1d83109SMatthew G. Knepley     ierr = PetscFree4(recvcounts, displs, Nc, Noff);CHKERRQ(ierr);
620e1d83109SMatthew G. Knepley     ierr = PetscFree2(nrankComp, match);CHKERRQ(ierr);
621e1d83109SMatthew G. Knepley 
622e1d83109SMatthew G. Knepley     /* Decide whether to flip cells in each component */
623e1d83109SMatthew G. Knepley     for (c = 0; c < cEnd-cStart; ++c) {if (flipped[cellComp[c]]) {ierr = PetscBTNegate(flippedCells, c);CHKERRQ(ierr);}}
624e1d83109SMatthew G. Knepley     ierr = PetscFree(flipped);CHKERRQ(ierr);
62584961fc4SMatthew G. Knepley   }
62684961fc4SMatthew G. Knepley   if (flg) {
62784961fc4SMatthew G. Knepley     PetscViewer v;
62884961fc4SMatthew G. Knepley 
6297e09831bSMatthew G. Knepley     ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr);
6301575c14dSBarry Smith     ierr = PetscViewerASCIIPushSynchronized(v);CHKERRQ(ierr);
63184961fc4SMatthew G. Knepley     ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank);CHKERRQ(ierr);
63284961fc4SMatthew G. Knepley     ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr);
6334d4c343aSBarry Smith     ierr = PetscViewerFlush(v);CHKERRQ(ierr);
6341575c14dSBarry Smith     ierr = PetscViewerASCIIPopSynchronized(v);CHKERRQ(ierr);
63584961fc4SMatthew G. Knepley   }
63684961fc4SMatthew G. Knepley   /* Reverse flipped cells in the mesh */
63784961fc4SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
6380d366550SMatthew G. Knepley     if (PetscBTLookup(flippedCells, c-cStart)) {
6390d366550SMatthew G. Knepley       ierr = DMPlexReverseCell(dm, c);CHKERRQ(ierr);
6400d366550SMatthew G. Knepley     }
64184961fc4SMatthew G. Knepley   }
64284961fc4SMatthew G. Knepley   ierr = PetscBTDestroy(&seenCells);CHKERRQ(ierr);
64384961fc4SMatthew G. Knepley   ierr = PetscBTDestroy(&flippedCells);CHKERRQ(ierr);
64484961fc4SMatthew G. Knepley   ierr = PetscBTDestroy(&seenFaces);CHKERRQ(ierr);
645e1d83109SMatthew G. Knepley   ierr = PetscFree2(numNeighbors, neighbors);CHKERRQ(ierr);
646e1d83109SMatthew G. Knepley   ierr = PetscFree2(rorntComp, lorntComp);CHKERRQ(ierr);
647e1d83109SMatthew G. Knepley   ierr = PetscFree3(faceFIFO, cellComp, faceComp);CHKERRQ(ierr);
64884961fc4SMatthew G. Knepley   PetscFunctionReturn(0);
64984961fc4SMatthew G. Knepley }
650