xref: /petsc/src/dm/impls/plex/plexorient.c (revision ad227fea26420e7a31ef041a9fa1418814d171a8)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
284961fc4SMatthew G. Knepley #include <petscsf.h>
384961fc4SMatthew G. Knepley 
484961fc4SMatthew G. Knepley /*@
527d0e99aSVaclav 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.
627d0e99aSVaclav Hapla 
727d0e99aSVaclav Hapla   Not Collective
827d0e99aSVaclav Hapla 
927d0e99aSVaclav Hapla   Input Parameters:
1027d0e99aSVaclav Hapla + dm              - The DM (DMPLEX)
1127d0e99aSVaclav Hapla . p               - The DAG point whose cone is compared
129dddd249SSatish Balay . mainConeSize    - Number of the reference cone points passed (at least 2 and at most size of the cone of p)
139dddd249SSatish Balay - mainCone        - The reference cone points
1427d0e99aSVaclav Hapla 
1527d0e99aSVaclav Hapla   Output Parameters:
1627d0e99aSVaclav Hapla + start           - The new starting point within the cone of p to make it conforming with the reference cone
1727d0e99aSVaclav Hapla - reverse         - The flag whether the order of the cone points should be reversed
1827d0e99aSVaclav Hapla 
1927d0e99aSVaclav Hapla   Level: advanced
2027d0e99aSVaclav Hapla 
2127d0e99aSVaclav Hapla .seealso: DMPlexOrient(), DMPlexOrientCell()
2227d0e99aSVaclav Hapla @*/
239dddd249SSatish Balay PetscErrorCode DMPlexCompareOrientations(DM dm, PetscInt p, PetscInt mainConeSize, const PetscInt mainCone[], PetscInt *start, PetscBool *reverse)
2427d0e99aSVaclav Hapla {
2527d0e99aSVaclav Hapla   PetscInt        coneSize;
2627d0e99aSVaclav Hapla   const PetscInt *cone;
2727d0e99aSVaclav Hapla   PetscInt        i, start_;
2827d0e99aSVaclav Hapla   PetscBool       reverse_;
2927d0e99aSVaclav Hapla   PetscErrorCode  ierr;
3027d0e99aSVaclav Hapla 
3127d0e99aSVaclav Hapla   PetscFunctionBegin;
3227d0e99aSVaclav Hapla   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3327d0e99aSVaclav Hapla   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
3427d0e99aSVaclav Hapla   if (coneSize < 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %D has no cone", p);
3527d0e99aSVaclav Hapla   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
369dddd249SSatish Balay   if (mainConeSize < 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %D: mainConeSize must be at least 2", p);
379dddd249SSatish Balay   if (mainConeSize > coneSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %D: mainConeSize must be at most coneSize", p);
3827d0e99aSVaclav Hapla   start_ = 0;
3927d0e99aSVaclav Hapla   for (i=0; i<coneSize; i++) {
409dddd249SSatish Balay     if (cone[i] == mainCone[0]) {
4127d0e99aSVaclav Hapla       start_ = i;
4227d0e99aSVaclav Hapla       break;
4327d0e99aSVaclav Hapla     }
4427d0e99aSVaclav Hapla   }
459dddd249SSatish Balay   if (PetscUnlikely(i==coneSize)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Point %D: starting point of reference cone not found in worker cone", p);
4627d0e99aSVaclav Hapla   reverse_ = PETSC_FALSE;
479dddd249SSatish Balay   for (i=0; i<mainConeSize; i++) {if (cone[(start_+i)%coneSize] != mainCone[i]) break;}
489dddd249SSatish Balay   if (i != mainConeSize) {
4927d0e99aSVaclav Hapla     reverse_ = PETSC_TRUE;
509dddd249SSatish Balay     for (i=0; i<mainConeSize; i++) {if (cone[(coneSize+start_-i)%coneSize] != mainCone[i]) break;}
519dddd249SSatish Balay     if (i < mainConeSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Point %D: cone has non-conforming order of points with respect to reference cone", p);
5227d0e99aSVaclav Hapla   }
5327d0e99aSVaclav Hapla   if (start) *start = start_;
5427d0e99aSVaclav Hapla   if (reverse) *reverse = reverse_;
559dddd249SSatish Balay   if (PetscUnlikely(cone[start_] != mainCone[0])) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D: cone[%d] = %d != %d = mainCone[0]", p, start_, cone[start_], mainCone[0]);
5627d0e99aSVaclav Hapla   PetscFunctionReturn(0);
5727d0e99aSVaclav Hapla }
5827d0e99aSVaclav Hapla 
5927d0e99aSVaclav Hapla /*@
6027d0e99aSVaclav Hapla   DMPlexOrientCell - Set the desired order of cone points of this DAG point, and fix orientations accordingly.
6127d0e99aSVaclav Hapla 
6227d0e99aSVaclav Hapla   Not Collective
6327d0e99aSVaclav Hapla 
6427d0e99aSVaclav Hapla   Input Parameters:
6527d0e99aSVaclav Hapla + dm              - The DM
6627d0e99aSVaclav Hapla . p               - The DAG point (from interval given by DMPlexGetChart())
679dddd249SSatish Balay . mainConeSize    - Number of specified cone points (at least 2)
689dddd249SSatish Balay - mainCone        - Specified cone points, i.e. ordered subset of current cone in DAG numbering (not cone-local numbering)
6927d0e99aSVaclav Hapla 
7027d0e99aSVaclav Hapla   Level: intermediate
7127d0e99aSVaclav Hapla 
7227d0e99aSVaclav Hapla .seealso: DMPlexOrient(), DMPlexGetCone(), DMPlexGetConeOrientation(), DMPlexInterpolate(), DMPlexGetChart()
7327d0e99aSVaclav Hapla @*/
749dddd249SSatish Balay PetscErrorCode DMPlexOrientCell(DM dm, PetscInt p, PetscInt mainConeSize, const PetscInt mainCone[])
7527d0e99aSVaclav Hapla {
7627d0e99aSVaclav Hapla   PetscInt        coneSize;
7727d0e99aSVaclav Hapla   PetscInt        start1=0;
7827d0e99aSVaclav Hapla   PetscBool       reverse1=PETSC_FALSE;
7927d0e99aSVaclav Hapla   PetscErrorCode  ierr;
8027d0e99aSVaclav Hapla 
8127d0e99aSVaclav Hapla   PetscFunctionBegin;
8227d0e99aSVaclav Hapla   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
839dddd249SSatish Balay   if (mainConeSize) PetscValidIntPointer(mainCone,4);
849dddd249SSatish Balay   if (mainConeSize == 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "mainConeSize cannot be 1");
8527d0e99aSVaclav Hapla   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
8627d0e99aSVaclav Hapla   if (!coneSize) PetscFunctionReturn(0); /* do nothing for points with no cone */
879dddd249SSatish Balay   ierr = DMPlexCompareOrientations(dm, p, mainConeSize, mainCone, &start1, &reverse1);CHKERRQ(ierr);
8827d0e99aSVaclav Hapla   ierr = DMPlexOrientCell_Internal(dm, p, start1, reverse1);CHKERRQ(ierr);
8976bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
9027d0e99aSVaclav Hapla     PetscInt        c;
9127d0e99aSVaclav Hapla     const PetscInt *cone;
9227d0e99aSVaclav Hapla     ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
939dddd249SSatish Balay     for (c = 0; c < mainConeSize; c++) {
949dddd249SSatish Balay       if (PetscUnlikely(cone[c] != mainCone[c])) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_PLIB, "The algorithm above is wrong as cone[%d] = %d != %d = mainCone[%d]", c, cone[c], mainCone[c], c);
9527d0e99aSVaclav Hapla     }
9627d0e99aSVaclav Hapla   }
9727d0e99aSVaclav Hapla   PetscFunctionReturn(0);
9827d0e99aSVaclav Hapla }
9927d0e99aSVaclav Hapla 
10027d0e99aSVaclav Hapla PetscErrorCode DMPlexOrientCell_Internal(DM dm, PetscInt p, PetscInt start1, PetscBool reverse1)
10127d0e99aSVaclav Hapla {
10227d0e99aSVaclav Hapla   PetscInt        i, j, k, maxConeSize, coneSize, coneConeSize, supportSize, supportConeSize;
10327d0e99aSVaclav Hapla   PetscInt        start0, start;
10427d0e99aSVaclav Hapla   PetscBool       reverse0, reverse;
10527d0e99aSVaclav Hapla   PetscInt        newornt;
10627d0e99aSVaclav Hapla   const PetscInt *cone=NULL, *support=NULL, *supportCone=NULL, *ornts=NULL;
10727d0e99aSVaclav Hapla   PetscInt       *newcone=NULL, *newornts=NULL;
10827d0e99aSVaclav Hapla   PetscErrorCode  ierr;
10927d0e99aSVaclav Hapla 
11027d0e99aSVaclav Hapla   PetscFunctionBegin;
11127d0e99aSVaclav Hapla   if (!start1 && !reverse1) PetscFunctionReturn(0);
11227d0e99aSVaclav Hapla   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
11327d0e99aSVaclav Hapla   if (!coneSize) PetscFunctionReturn(0); /* do nothing for points with no cone */
11427d0e99aSVaclav Hapla   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
11527d0e99aSVaclav Hapla   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
11627d0e99aSVaclav Hapla   /* permute p's cone and orientations */
11727d0e99aSVaclav Hapla   ierr = DMPlexGetConeOrientation(dm, p, &ornts);CHKERRQ(ierr);
11827d0e99aSVaclav Hapla   ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &newcone);CHKERRQ(ierr);
11927d0e99aSVaclav Hapla   ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &newornts);CHKERRQ(ierr);
12027d0e99aSVaclav Hapla   ierr = DMPlexFixFaceOrientations_Permute_Private(coneSize, cone, start1, reverse1, newcone);CHKERRQ(ierr);
12127d0e99aSVaclav Hapla   ierr = DMPlexFixFaceOrientations_Permute_Private(coneSize, ornts, start1, reverse1, newornts);CHKERRQ(ierr);
12227d0e99aSVaclav Hapla   /* if direction of p (face) is flipped, flip also p's cone points (edges) */
12327d0e99aSVaclav Hapla   if (reverse1) {
12427d0e99aSVaclav Hapla     for (i=0; i<coneSize; i++) {
12527d0e99aSVaclav Hapla       ierr = DMPlexGetConeSize(dm, cone[i], &coneConeSize);CHKERRQ(ierr);
12627d0e99aSVaclav Hapla       ierr = DMPlexFixFaceOrientations_Translate_Private(newornts[i], &start0, &reverse0);CHKERRQ(ierr);
12727d0e99aSVaclav Hapla       ierr = DMPlexFixFaceOrientations_Combine_Private(coneConeSize, start0, reverse0, 1, PETSC_FALSE, &start, &reverse);CHKERRQ(ierr);
12827d0e99aSVaclav Hapla       ierr = DMPlexFixFaceOrientations_TranslateBack_Private(coneConeSize, start, reverse, &newornts[i]);CHKERRQ(ierr);
12927d0e99aSVaclav Hapla     }
13027d0e99aSVaclav Hapla   }
13127d0e99aSVaclav Hapla   ierr = DMPlexSetConeOrientation(dm, p, newornts);CHKERRQ(ierr);
13227d0e99aSVaclav Hapla   /* fix oriention of p within cones of p's support points */
13327d0e99aSVaclav Hapla   ierr = DMPlexGetSupport(dm, p, &support);CHKERRQ(ierr);
13427d0e99aSVaclav Hapla   ierr = DMPlexGetSupportSize(dm, p, &supportSize);CHKERRQ(ierr);
13527d0e99aSVaclav Hapla   for (j=0; j<supportSize; j++) {
13627d0e99aSVaclav Hapla     ierr = DMPlexGetCone(dm, support[j], &supportCone);CHKERRQ(ierr);
13727d0e99aSVaclav Hapla     ierr = DMPlexGetConeSize(dm, support[j], &supportConeSize);CHKERRQ(ierr);
13827d0e99aSVaclav Hapla     ierr = DMPlexGetConeOrientation(dm, support[j], &ornts);CHKERRQ(ierr);
13927d0e99aSVaclav Hapla     for (k=0; k<supportConeSize; k++) {
14027d0e99aSVaclav Hapla       if (supportCone[k] != p) continue;
14127d0e99aSVaclav Hapla       ierr = DMPlexFixFaceOrientations_Translate_Private(ornts[k], &start0, &reverse0);CHKERRQ(ierr);
14227d0e99aSVaclav Hapla       ierr = DMPlexFixFaceOrientations_Combine_Private(coneSize, start0, reverse0, start1, reverse1, &start, &reverse);CHKERRQ(ierr);
14327d0e99aSVaclav Hapla       ierr = DMPlexFixFaceOrientations_TranslateBack_Private(coneSize, start, reverse, &newornt);CHKERRQ(ierr);
14427d0e99aSVaclav Hapla       ierr = DMPlexInsertConeOrientation(dm, support[j], k, newornt);CHKERRQ(ierr);
14527d0e99aSVaclav Hapla     }
14627d0e99aSVaclav Hapla   }
14727d0e99aSVaclav Hapla   /* rewrite cone */
14827d0e99aSVaclav Hapla   ierr = DMPlexSetCone(dm, p, newcone);CHKERRQ(ierr);
14927d0e99aSVaclav Hapla   ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &newcone);CHKERRQ(ierr);
15027d0e99aSVaclav Hapla   ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &newornts);CHKERRQ(ierr);
15127d0e99aSVaclav Hapla   PetscFunctionReturn(0);
15227d0e99aSVaclav Hapla }
15327d0e99aSVaclav Hapla 
15427d0e99aSVaclav Hapla /*@
15584961fc4SMatthew G. Knepley   DMPlexReverseCell - Give a mesh cell the opposite orientation
15684961fc4SMatthew G. Knepley 
15784961fc4SMatthew G. Knepley   Input Parameters:
15884961fc4SMatthew G. Knepley + dm   - The DM
15984961fc4SMatthew G. Knepley - cell - The cell number
16084961fc4SMatthew G. Knepley 
16184961fc4SMatthew G. Knepley   Note: The modification of the DM is done in-place.
16284961fc4SMatthew G. Knepley 
16384961fc4SMatthew G. Knepley   Level: advanced
16484961fc4SMatthew G. Knepley 
16584961fc4SMatthew G. Knepley .seealso: DMPlexOrient(), DMCreate(), DMPLEX
16684961fc4SMatthew G. Knepley @*/
16784961fc4SMatthew G. Knepley PetscErrorCode DMPlexReverseCell(DM dm, PetscInt cell)
16884961fc4SMatthew G. Knepley {
16984961fc4SMatthew G. Knepley   /* Note that the reverse orientation ro of a face with orientation o is:
17084961fc4SMatthew G. Knepley 
17184961fc4SMatthew G. Knepley        ro = o >= 0 ? -(faceSize - o) : faceSize + o
17284961fc4SMatthew G. Knepley 
17384961fc4SMatthew G. Knepley      where faceSize is the size of the cone for the face.
17484961fc4SMatthew G. Knepley   */
17584961fc4SMatthew G. Knepley   const PetscInt *cone,    *coneO, *support;
17684961fc4SMatthew G. Knepley   PetscInt       *revcone, *revconeO;
17784961fc4SMatthew G. Knepley   PetscInt        maxConeSize, coneSize, supportSize, faceSize, cp, sp;
17884961fc4SMatthew G. Knepley   PetscErrorCode  ierr;
17984961fc4SMatthew G. Knepley 
18084961fc4SMatthew G. Knepley   PetscFunctionBegin;
18184961fc4SMatthew G. Knepley   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
18269291d52SBarry Smith   ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &revcone);CHKERRQ(ierr);
18369291d52SBarry Smith   ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &revconeO);CHKERRQ(ierr);
18484961fc4SMatthew G. Knepley   /* Reverse cone, and reverse orientations of faces */
18584961fc4SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
18684961fc4SMatthew G. Knepley   ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
18784961fc4SMatthew G. Knepley   ierr = DMPlexGetConeOrientation(dm, cell, &coneO);CHKERRQ(ierr);
18884961fc4SMatthew G. Knepley   for (cp = 0; cp < coneSize; ++cp) {
18984961fc4SMatthew G. Knepley     const PetscInt rcp = coneSize-cp-1;
19084961fc4SMatthew G. Knepley 
19184961fc4SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, cone[rcp], &faceSize);CHKERRQ(ierr);
19284961fc4SMatthew G. Knepley     revcone[cp]  = cone[rcp];
19384961fc4SMatthew G. Knepley     revconeO[cp] = coneO[rcp] >= 0 ? -(faceSize-coneO[rcp]) : faceSize+coneO[rcp];
19484961fc4SMatthew G. Knepley   }
19584961fc4SMatthew G. Knepley   ierr = DMPlexSetCone(dm, cell, revcone);CHKERRQ(ierr);
19684961fc4SMatthew G. Knepley   ierr = DMPlexSetConeOrientation(dm, cell, revconeO);CHKERRQ(ierr);
19784961fc4SMatthew G. Knepley   /* Reverse orientation of this cell in the support hypercells */
19884961fc4SMatthew G. Knepley   faceSize = coneSize;
19984961fc4SMatthew G. Knepley   ierr = DMPlexGetSupportSize(dm, cell, &supportSize);CHKERRQ(ierr);
20084961fc4SMatthew G. Knepley   ierr = DMPlexGetSupport(dm, cell, &support);CHKERRQ(ierr);
20184961fc4SMatthew G. Knepley   for (sp = 0; sp < supportSize; ++sp) {
20284961fc4SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, support[sp], &coneSize);CHKERRQ(ierr);
20384961fc4SMatthew G. Knepley     ierr = DMPlexGetCone(dm, support[sp], &cone);CHKERRQ(ierr);
20484961fc4SMatthew G. Knepley     ierr = DMPlexGetConeOrientation(dm, support[sp], &coneO);CHKERRQ(ierr);
20584961fc4SMatthew G. Knepley     for (cp = 0; cp < coneSize; ++cp) {
20684961fc4SMatthew G. Knepley       if (cone[cp] != cell) continue;
20784961fc4SMatthew G. Knepley       ierr = DMPlexInsertConeOrientation(dm, support[sp], cp, coneO[cp] >= 0 ? -(faceSize-coneO[cp]) : faceSize+coneO[cp]);CHKERRQ(ierr);
20884961fc4SMatthew G. Knepley     }
20984961fc4SMatthew G. Knepley   }
21069291d52SBarry Smith   ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &revcone);CHKERRQ(ierr);
21169291d52SBarry Smith   ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &revconeO);CHKERRQ(ierr);
21284961fc4SMatthew G. Knepley   PetscFunctionReturn(0);
21384961fc4SMatthew G. Knepley }
21484961fc4SMatthew G. Knepley 
2157b310ce2SMatthew G. Knepley /*
2167b310ce2SMatthew G. Knepley   - Checks face match
2177b310ce2SMatthew G. Knepley     - Flips non-matching
2187b310ce2SMatthew G. Knepley   - Inserts faces of support cells in FIFO
2197b310ce2SMatthew G. Knepley */
2207b2de0fdSMatthew 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)
2217b310ce2SMatthew G. Knepley {
2227b310ce2SMatthew G. Knepley   const PetscInt *support, *coneA, *coneB, *coneOA, *coneOB;
2237b310ce2SMatthew G. Knepley   PetscInt        supportSize, coneSizeA, coneSizeB, posA = -1, posB = -1;
2247b2de0fdSMatthew G. Knepley   PetscInt        face, dim, seenA, flippedA, seenB, flippedB, mismatch, c;
2257b310ce2SMatthew G. Knepley   PetscErrorCode  ierr;
2267b310ce2SMatthew G. Knepley 
2277b310ce2SMatthew G. Knepley   PetscFunctionBegin;
2287b310ce2SMatthew G. Knepley   face = faceFIFO[(*fTop)++];
2297b310ce2SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2307b310ce2SMatthew G. Knepley   ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
2317b310ce2SMatthew G. Knepley   ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr);
2327b310ce2SMatthew G. Knepley   if (supportSize < 2) PetscFunctionReturn(0);
2337b310ce2SMatthew G. Knepley   if (supportSize != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %d", supportSize);
2347b310ce2SMatthew G. Knepley   seenA    = PetscBTLookup(seenCells,    support[0]-cStart);
2357b310ce2SMatthew G. Knepley   flippedA = PetscBTLookup(flippedCells, support[0]-cStart) ? 1 : 0;
2367b310ce2SMatthew G. Knepley   seenB    = PetscBTLookup(seenCells,    support[1]-cStart);
2377b310ce2SMatthew G. Knepley   flippedB = PetscBTLookup(flippedCells, support[1]-cStart) ? 1 : 0;
2387b310ce2SMatthew G. Knepley 
2397b310ce2SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, support[0], &coneSizeA);CHKERRQ(ierr);
2407b310ce2SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, support[1], &coneSizeB);CHKERRQ(ierr);
2417b310ce2SMatthew G. Knepley   ierr = DMPlexGetCone(dm, support[0], &coneA);CHKERRQ(ierr);
2427b310ce2SMatthew G. Knepley   ierr = DMPlexGetCone(dm, support[1], &coneB);CHKERRQ(ierr);
2437b310ce2SMatthew G. Knepley   ierr = DMPlexGetConeOrientation(dm, support[0], &coneOA);CHKERRQ(ierr);
2447b310ce2SMatthew G. Knepley   ierr = DMPlexGetConeOrientation(dm, support[1], &coneOB);CHKERRQ(ierr);
2457b310ce2SMatthew G. Knepley   for (c = 0; c < coneSizeA; ++c) {
2467b310ce2SMatthew G. Knepley     if (!PetscBTLookup(seenFaces, coneA[c]-fStart)) {
2477b310ce2SMatthew G. Knepley       faceFIFO[(*fBottom)++] = coneA[c];
2487b310ce2SMatthew G. Knepley       ierr = PetscBTSet(seenFaces, coneA[c]-fStart);CHKERRQ(ierr);
2497b310ce2SMatthew G. Knepley     }
2507b310ce2SMatthew G. Knepley     if (coneA[c] == face) posA = c;
2517b310ce2SMatthew 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);
2527b310ce2SMatthew G. Knepley   }
2537b310ce2SMatthew 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]);
2547b310ce2SMatthew G. Knepley   for (c = 0; c < coneSizeB; ++c) {
2557b310ce2SMatthew G. Knepley     if (!PetscBTLookup(seenFaces, coneB[c]-fStart)) {
2567b310ce2SMatthew G. Knepley       faceFIFO[(*fBottom)++] = coneB[c];
2577b310ce2SMatthew G. Knepley       ierr = PetscBTSet(seenFaces, coneB[c]-fStart);CHKERRQ(ierr);
2587b310ce2SMatthew G. Knepley     }
2597b310ce2SMatthew G. Knepley     if (coneB[c] == face) posB = c;
2607b310ce2SMatthew 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);
2617b310ce2SMatthew G. Knepley   }
2627b310ce2SMatthew 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]);
2637b310ce2SMatthew G. Knepley 
2647b310ce2SMatthew G. Knepley   if (dim == 1) {
2657b310ce2SMatthew G. Knepley     mismatch = posA == posB;
2667b310ce2SMatthew G. Knepley   } else {
2677b310ce2SMatthew G. Knepley     mismatch = coneOA[posA] == coneOB[posB];
2687b310ce2SMatthew G. Knepley   }
2697b310ce2SMatthew G. Knepley 
2707b310ce2SMatthew G. Knepley   if (mismatch ^ (flippedA ^ flippedB)) {
2717b310ce2SMatthew 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]);
2727b310ce2SMatthew G. Knepley     if (!seenA && !flippedA) {
2737b310ce2SMatthew G. Knepley       ierr = PetscBTSet(flippedCells, support[0]-cStart);CHKERRQ(ierr);
2747b310ce2SMatthew G. Knepley     } else if (!seenB && !flippedB) {
2757b310ce2SMatthew G. Knepley       ierr = PetscBTSet(flippedCells, support[1]-cStart);CHKERRQ(ierr);
2767b310ce2SMatthew G. Knepley     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
2777b310ce2SMatthew 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");
2787b310ce2SMatthew G. Knepley   ierr = PetscBTSet(seenCells, support[0]-cStart);CHKERRQ(ierr);
2797b310ce2SMatthew G. Knepley   ierr = PetscBTSet(seenCells, support[1]-cStart);CHKERRQ(ierr);
2807b310ce2SMatthew G. Knepley   PetscFunctionReturn(0);
2817b310ce2SMatthew G. Knepley }
2827b310ce2SMatthew G. Knepley 
28384961fc4SMatthew G. Knepley /*@
28484961fc4SMatthew G. Knepley   DMPlexOrient - Give a consistent orientation to the input mesh
28584961fc4SMatthew G. Knepley 
28684961fc4SMatthew G. Knepley   Input Parameters:
28784961fc4SMatthew G. Knepley . dm - The DM
28884961fc4SMatthew G. Knepley 
28984961fc4SMatthew G. Knepley   Note: The orientation data for the DM are change in-place.
29084961fc4SMatthew G. Knepley $ This routine will fail for non-orientable surfaces, such as the Moebius strip.
29184961fc4SMatthew G. Knepley 
29284961fc4SMatthew G. Knepley   Level: advanced
29384961fc4SMatthew G. Knepley 
29484961fc4SMatthew G. Knepley .seealso: DMCreate(), DMPLEX
29584961fc4SMatthew G. Knepley @*/
29684961fc4SMatthew G. Knepley PetscErrorCode DMPlexOrient(DM dm)
29784961fc4SMatthew G. Knepley {
29884961fc4SMatthew G. Knepley   MPI_Comm           comm;
299e1d83109SMatthew G. Knepley   PetscSF            sf;
300e1d83109SMatthew G. Knepley   const PetscInt    *lpoints;
301e1d83109SMatthew G. Knepley   const PetscSFNode *rpoints;
302e1d83109SMatthew G. Knepley   PetscSFNode       *rorntComp = NULL, *lorntComp = NULL;
303e1d83109SMatthew G. Knepley   PetscInt          *numNeighbors, **neighbors;
304e1d83109SMatthew G. Knepley   PetscSFNode       *nrankComp;
305e1d83109SMatthew G. Knepley   PetscBool         *match, *flipped;
30684961fc4SMatthew G. Knepley   PetscBT            seenCells, flippedCells, seenFaces;
307e1d83109SMatthew G. Knepley   PetscInt          *faceFIFO, fTop, fBottom, *cellComp, *faceComp;
3087cadcfe8SMatthew G. Knepley   PetscInt           numLeaves, numRoots, dim, h, cStart, cEnd, c, cell, fStart, fEnd, face, off, totNeighbors = 0;
30939526728SToby Isaac   PetscMPIInt        rank, size, numComponents, comp = 0;
31039526728SToby Isaac   PetscBool          flg, flg2;
311a17aa47eSToby Isaac   PetscViewer        viewer = NULL, selfviewer = NULL;
31284961fc4SMatthew G. Knepley   PetscErrorCode     ierr;
31384961fc4SMatthew G. Knepley 
31484961fc4SMatthew G. Knepley   PetscFunctionBegin;
31584961fc4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
316ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
317ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
318c5929fdfSBarry Smith   ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view", &flg);CHKERRQ(ierr);
31939526728SToby Isaac   ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view_synchronized", &flg2);CHKERRQ(ierr);
320e1d83109SMatthew G. Knepley   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
321e1d83109SMatthew G. Knepley   ierr = PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints);CHKERRQ(ierr);
32284961fc4SMatthew G. Knepley   /* Truth Table
32384961fc4SMatthew G. Knepley      mismatch    flips   do action   mismatch   flipA ^ flipB   action
32484961fc4SMatthew G. Knepley          F       0 flips     no         F             F           F
32584961fc4SMatthew G. Knepley          F       1 flip      yes        F             T           T
32684961fc4SMatthew G. Knepley          F       2 flips     no         T             F           T
32784961fc4SMatthew G. Knepley          T       0 flips     yes        T             T           F
32884961fc4SMatthew G. Knepley          T       1 flip      no
32984961fc4SMatthew G. Knepley          T       2 flips     yes
33084961fc4SMatthew G. Knepley   */
33184961fc4SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
33284961fc4SMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(dm, &h);CHKERRQ(ierr);
33384961fc4SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, h,   &cStart, &cEnd);CHKERRQ(ierr);
33484961fc4SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, h+1, &fStart, &fEnd);CHKERRQ(ierr);
33584961fc4SMatthew G. Knepley   ierr = PetscBTCreate(cEnd - cStart, &seenCells);CHKERRQ(ierr);
33684961fc4SMatthew G. Knepley   ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr);
33784961fc4SMatthew G. Knepley   ierr = PetscBTCreate(cEnd - cStart, &flippedCells);CHKERRQ(ierr);
33884961fc4SMatthew G. Knepley   ierr = PetscBTMemzero(cEnd - cStart, flippedCells);CHKERRQ(ierr);
33984961fc4SMatthew G. Knepley   ierr = PetscBTCreate(fEnd - fStart, &seenFaces);CHKERRQ(ierr);
34084961fc4SMatthew G. Knepley   ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr);
341e1d83109SMatthew G. Knepley   ierr = PetscCalloc3(fEnd - fStart, &faceFIFO, cEnd-cStart, &cellComp, fEnd-fStart, &faceComp);CHKERRQ(ierr);
342e1d83109SMatthew G. Knepley   /*
343e1d83109SMatthew G. Knepley    OLD STYLE
344e1d83109SMatthew G. Knepley    - Add an integer array over cells and faces (component) for connected component number
345e1d83109SMatthew G. Knepley    Foreach component
346e1d83109SMatthew G. Knepley      - Mark the initial cell as seen
347e1d83109SMatthew G. Knepley      - Process component as usual
348e1d83109SMatthew G. Knepley      - Set component for all seenCells
349e1d83109SMatthew G. Knepley      - Wipe seenCells and seenFaces (flippedCells can stay)
350e1d83109SMatthew G. Knepley    - Generate parallel adjacency for component using SF and seenFaces
351e1d83109SMatthew G. Knepley    - Collect numComponents adj data from each proc to 0
352e1d83109SMatthew G. Knepley    - Build same serial graph
353e1d83109SMatthew G. Knepley    - Use same solver
354e1d83109SMatthew G. Knepley    - Use Scatterv to to send back flipped flags for each component
355e1d83109SMatthew G. Knepley    - Negate flippedCells by component
356e1d83109SMatthew G. Knepley 
357e1d83109SMatthew G. Knepley    NEW STYLE
358e1d83109SMatthew G. Knepley    - Create the adj on each process
359e1d83109SMatthew G. Knepley    - Bootstrap to complete graph on proc 0
360e1d83109SMatthew G. Knepley   */
361e1d83109SMatthew G. Knepley   /* Loop over components */
362e1d83109SMatthew G. Knepley   for (cell = cStart; cell < cEnd; ++cell) cellComp[cell-cStart] = -1;
363e1d83109SMatthew G. Knepley   do {
364e1d83109SMatthew G. Knepley     /* Look for first unmarked cell */
365e1d83109SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) if (cellComp[cell-cStart] < 0) break;
366e1d83109SMatthew G. Knepley     if (cell >= cEnd) break;
367e1d83109SMatthew G. Knepley     /* Initialize FIFO with first cell in component */
368e1d83109SMatthew G. Knepley     {
36984961fc4SMatthew G. Knepley       const PetscInt *cone;
37084961fc4SMatthew G. Knepley       PetscInt        coneSize;
37184961fc4SMatthew G. Knepley 
372e1d83109SMatthew G. Knepley       fTop = fBottom = 0;
373e1d83109SMatthew G. Knepley       ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
374e1d83109SMatthew G. Knepley       ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
37584961fc4SMatthew G. Knepley       for (c = 0; c < coneSize; ++c) {
37684961fc4SMatthew G. Knepley         faceFIFO[fBottom++] = cone[c];
37784961fc4SMatthew G. Knepley         ierr = PetscBTSet(seenFaces, cone[c]-fStart);CHKERRQ(ierr);
37884961fc4SMatthew G. Knepley       }
37931c8331aSMatthew G. Knepley       ierr = PetscBTSet(seenCells, cell-cStart);CHKERRQ(ierr);
38084961fc4SMatthew G. Knepley     }
38184961fc4SMatthew G. Knepley     /* Consider each face in FIFO */
38284961fc4SMatthew G. Knepley     while (fTop < fBottom) {
3837b2de0fdSMatthew G. Knepley       ierr = DMPlexCheckFace_Internal(dm, faceFIFO, &fTop, &fBottom, cStart, fStart, fEnd, seenCells, flippedCells, seenFaces);CHKERRQ(ierr);
38484961fc4SMatthew G. Knepley     }
385e1d83109SMatthew G. Knepley     /* Set component for cells and faces */
386e1d83109SMatthew G. Knepley     for (cell = 0; cell < cEnd-cStart; ++cell) {
387e1d83109SMatthew G. Knepley       if (PetscBTLookup(seenCells, cell)) cellComp[cell] = comp;
388e1d83109SMatthew G. Knepley     }
389e1d83109SMatthew G. Knepley     for (face = 0; face < fEnd-fStart; ++face) {
390e1d83109SMatthew G. Knepley       if (PetscBTLookup(seenFaces, face)) faceComp[face] = comp;
391e1d83109SMatthew G. Knepley     }
392e1d83109SMatthew G. Knepley     /* Wipe seenCells and seenFaces for next component */
393e1d83109SMatthew G. Knepley     ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr);
394e1d83109SMatthew G. Knepley     ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr);
395e1d83109SMatthew G. Knepley     ++comp;
396e1d83109SMatthew G. Knepley   } while (1);
397e1d83109SMatthew G. Knepley   numComponents = comp;
39884961fc4SMatthew G. Knepley   if (flg) {
39984961fc4SMatthew G. Knepley     PetscViewer v;
40084961fc4SMatthew G. Knepley 
4017e09831bSMatthew G. Knepley     ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr);
4021575c14dSBarry Smith     ierr = PetscViewerASCIIPushSynchronized(v);CHKERRQ(ierr);
40384961fc4SMatthew G. Knepley     ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank);CHKERRQ(ierr);
40484961fc4SMatthew G. Knepley     ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr);
4054d4c343aSBarry Smith     ierr = PetscViewerFlush(v);CHKERRQ(ierr);
4061575c14dSBarry Smith     ierr = PetscViewerASCIIPopSynchronized(v);CHKERRQ(ierr);
40784961fc4SMatthew G. Knepley   }
40884961fc4SMatthew G. Knepley   /* Now all subdomains are oriented, but we need a consistent parallel orientation */
40984961fc4SMatthew G. Knepley   if (numLeaves >= 0) {
410e1d83109SMatthew G. Knepley     /* Store orientations of boundary faces*/
411e1d83109SMatthew G. Knepley     ierr = PetscCalloc2(numRoots,&rorntComp,numRoots,&lorntComp);CHKERRQ(ierr);
41284961fc4SMatthew G. Knepley     for (face = fStart; face < fEnd; ++face) {
413e1d83109SMatthew G. Knepley       const PetscInt *cone, *support, *ornt;
414e1d83109SMatthew G. Knepley       PetscInt        coneSize, supportSize;
415e1d83109SMatthew G. Knepley 
41684961fc4SMatthew G. Knepley       ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
41784961fc4SMatthew G. Knepley       if (supportSize != 1) continue;
41884961fc4SMatthew G. Knepley       ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr);
41984961fc4SMatthew G. Knepley 
42084961fc4SMatthew G. Knepley       ierr = DMPlexGetCone(dm, support[0], &cone);CHKERRQ(ierr);
42184961fc4SMatthew G. Knepley       ierr = DMPlexGetConeSize(dm, support[0], &coneSize);CHKERRQ(ierr);
42284961fc4SMatthew G. Knepley       ierr = DMPlexGetConeOrientation(dm, support[0], &ornt);CHKERRQ(ierr);
42384961fc4SMatthew G. Knepley       for (c = 0; c < coneSize; ++c) if (cone[c] == face) break;
42484961fc4SMatthew G. Knepley       if (dim == 1) {
42584961fc4SMatthew G. Knepley         /* Use cone position instead, shifted to -1 or 1 */
4266e1eb25cSMatthew G. Knepley         if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = 1-c*2;
4276e1eb25cSMatthew G. Knepley         else                                                rorntComp[face].rank = c*2-1;
42884961fc4SMatthew G. Knepley       } else {
429e1d83109SMatthew G. Knepley         if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 :  1;
430e1d83109SMatthew G. Knepley         else                                                rorntComp[face].rank = ornt[c] < 0 ?  1 : -1;
43184961fc4SMatthew G. Knepley       }
432e1d83109SMatthew G. Knepley       rorntComp[face].index = faceComp[face-fStart];
43384961fc4SMatthew G. Knepley     }
434e1d83109SMatthew G. Knepley     /* Communicate boundary edge orientations */
435*ad227feaSJunchao Zhang     ierr = PetscSFBcastBegin(sf, MPIU_2INT, rorntComp, lorntComp,MPI_REPLACE);CHKERRQ(ierr);
436*ad227feaSJunchao Zhang     ierr = PetscSFBcastEnd(sf, MPIU_2INT, rorntComp, lorntComp,MPI_REPLACE);CHKERRQ(ierr);
437e1d83109SMatthew G. Knepley   }
438e1d83109SMatthew G. Knepley   /* Get process adjacency */
439e1d83109SMatthew G. Knepley   ierr = PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors);CHKERRQ(ierr);
440a17aa47eSToby Isaac   viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm));
441a17aa47eSToby Isaac   if (flg2) {ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);}
442a17aa47eSToby Isaac   ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&selfviewer);CHKERRQ(ierr);
443e1d83109SMatthew G. Knepley   for (comp = 0; comp < numComponents; ++comp) {
444e1d83109SMatthew G. Knepley     PetscInt  l, n;
44584961fc4SMatthew G. Knepley 
446e1d83109SMatthew G. Knepley     numNeighbors[comp] = 0;
44739ea070eSMatthew G. Knepley     ierr = PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp]);CHKERRQ(ierr);
448e1d83109SMatthew G. Knepley     /* I know this is p^2 time in general, but for bounded degree its alright */
449e1d83109SMatthew G. Knepley     for (l = 0; l < numLeaves; ++l) {
450e1d83109SMatthew G. Knepley       const PetscInt face = lpoints[l];
451e1d83109SMatthew G. Knepley 
452e1d83109SMatthew G. Knepley       /* Find a representative face (edge) separating pairs of procs */
453e1d83109SMatthew G. Knepley       if ((face >= fStart) && (face < fEnd) && (faceComp[face-fStart] == comp)) {
454269a49d5SMatthew G. Knepley         const PetscInt rrank = rpoints[l].rank;
455e1d83109SMatthew G. Knepley         const PetscInt rcomp = lorntComp[face].index;
456e1d83109SMatthew G. Knepley 
457269a49d5SMatthew G. Knepley         for (n = 0; n < numNeighbors[comp]; ++n) if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break;
458e1d83109SMatthew G. Knepley         if (n >= numNeighbors[comp]) {
459e1d83109SMatthew G. Knepley           PetscInt supportSize;
460e1d83109SMatthew G. Knepley 
461e1d83109SMatthew G. Knepley           ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
462e1d83109SMatthew G. Knepley           if (supportSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %d", supportSize);
46335041eceSToby 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);}
464e1d83109SMatthew G. Knepley           neighbors[comp][numNeighbors[comp]++] = l;
465e1d83109SMatthew G. Knepley         }
466e1d83109SMatthew G. Knepley       }
467e1d83109SMatthew G. Knepley     }
468e1d83109SMatthew G. Knepley     totNeighbors += numNeighbors[comp];
469e1d83109SMatthew G. Knepley   }
470a17aa47eSToby Isaac   ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&selfviewer);CHKERRQ(ierr);
471a17aa47eSToby Isaac   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
472a17aa47eSToby Isaac   if (flg2) {ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);}
473e1d83109SMatthew G. Knepley   ierr = PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match);CHKERRQ(ierr);
474e1d83109SMatthew G. Knepley   for (comp = 0, off = 0; comp < numComponents; ++comp) {
475e1d83109SMatthew G. Knepley     PetscInt n;
476e1d83109SMatthew G. Knepley 
477e1d83109SMatthew G. Knepley     for (n = 0; n < numNeighbors[comp]; ++n, ++off) {
478e1d83109SMatthew G. Knepley       const PetscInt face = lpoints[neighbors[comp][n]];
479e1d83109SMatthew G. Knepley       const PetscInt o    = rorntComp[face].rank*lorntComp[face].rank;
480e1d83109SMatthew G. Knepley 
481e1d83109SMatthew G. Knepley       if      (o < 0) match[off] = PETSC_TRUE;
482e1d83109SMatthew G. Knepley       else if (o > 0) match[off] = PETSC_FALSE;
483e1d83109SMatthew 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);
484e1d83109SMatthew G. Knepley       nrankComp[off].rank  = rpoints[neighbors[comp][n]].rank;
485e1d83109SMatthew G. Knepley       nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index;
486e1d83109SMatthew G. Knepley     }
487e1d83109SMatthew G. Knepley     ierr = PetscFree(neighbors[comp]);CHKERRQ(ierr);
48884961fc4SMatthew G. Knepley   }
48984961fc4SMatthew G. Knepley   /* Collect the graph on 0 */
490e1d83109SMatthew G. Knepley   if (numLeaves >= 0) {
49184961fc4SMatthew G. Knepley     Mat          G;
49284961fc4SMatthew G. Knepley     PetscBT      seenProcs, flippedProcs;
49384961fc4SMatthew G. Knepley     PetscInt    *procFIFO, pTop, pBottom;
4947cadcfe8SMatthew G. Knepley     PetscInt    *N   = NULL, *Noff;
495e1d83109SMatthew G. Knepley     PetscSFNode *adj = NULL;
49684961fc4SMatthew G. Knepley     PetscBool   *val = NULL;
4977cadcfe8SMatthew G. Knepley     PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc, p, o;
4989852e123SBarry Smith     PetscMPIInt  size = 0;
49984961fc4SMatthew G. Knepley 
500e1d83109SMatthew G. Knepley     ierr = PetscCalloc1(numComponents, &flipped);CHKERRQ(ierr);
501ffc4695bSBarry Smith     if (!rank) {ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);}
5029852e123SBarry Smith     ierr = PetscCalloc4(size, &recvcounts, size+1, &displs, size, &Nc, size+1, &Noff);CHKERRQ(ierr);
503ffc4695bSBarry Smith     ierr = MPI_Gather(&numComponents, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm);CHKERRMPI(ierr);
5049852e123SBarry Smith     for (p = 0; p < size; ++p) {
505e1d83109SMatthew G. Knepley       displs[p+1] = displs[p] + Nc[p];
506e1d83109SMatthew G. Knepley     }
5079852e123SBarry Smith     if (!rank) {ierr = PetscMalloc1(displs[size],&N);CHKERRQ(ierr);}
508ffc4695bSBarry Smith     ierr = MPI_Gatherv(numNeighbors, numComponents, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm);CHKERRMPI(ierr);
5099852e123SBarry Smith     for (p = 0, o = 0; p < size; ++p) {
510e1d83109SMatthew G. Knepley       recvcounts[p] = 0;
511e1d83109SMatthew G. Knepley       for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o];
51284961fc4SMatthew G. Knepley       displs[p+1] = displs[p] + recvcounts[p];
51384961fc4SMatthew G. Knepley     }
5149852e123SBarry Smith     if (!rank) {ierr = PetscMalloc2(displs[size], &adj, displs[size], &val);CHKERRQ(ierr);}
515ffc4695bSBarry Smith     ierr = MPI_Gatherv(nrankComp, totNeighbors, MPIU_2INT, adj, recvcounts, displs, MPIU_2INT, 0, comm);CHKERRMPI(ierr);
516ffc4695bSBarry Smith     ierr = MPI_Gatherv(match, totNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm);CHKERRMPI(ierr);
517e1d83109SMatthew G. Knepley     ierr = PetscFree2(numNeighbors, neighbors);CHKERRQ(ierr);
518e1d83109SMatthew G. Knepley     if (!rank) {
5199852e123SBarry Smith       for (p = 1; p <= size; ++p) {Noff[p] = Noff[p-1] + Nc[p-1];}
520a83982f3SMatthew G. Knepley       if (flg) {
521e1d83109SMatthew G. Knepley         PetscInt n;
522e1d83109SMatthew G. Knepley 
5239852e123SBarry Smith         for (p = 0, off = 0; p < size; ++p) {
524e1d83109SMatthew G. Knepley           for (c = 0; c < Nc[p]; ++c) {
525302440fdSBarry Smith             ierr = PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %d:\n", p, c);CHKERRQ(ierr);
526e1d83109SMatthew G. Knepley             for (n = 0; n < N[Noff[p]+c]; ++n, ++off) {
527302440fdSBarry Smith               ierr = PetscPrintf(PETSC_COMM_SELF, "  edge (%d, %d) (%d):\n", adj[off].rank, adj[off].index, val[off]);CHKERRQ(ierr);
528e1d83109SMatthew G. Knepley             }
52984961fc4SMatthew G. Knepley           }
53084961fc4SMatthew G. Knepley         }
53184961fc4SMatthew G. Knepley       }
53284961fc4SMatthew G. Knepley       /* Symmetrize the graph */
53384961fc4SMatthew G. Knepley       ierr = MatCreate(PETSC_COMM_SELF, &G);CHKERRQ(ierr);
5349852e123SBarry Smith       ierr = MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size]);CHKERRQ(ierr);
53584961fc4SMatthew G. Knepley       ierr = MatSetUp(G);CHKERRQ(ierr);
5369852e123SBarry Smith       for (p = 0, off = 0; p < size; ++p) {
537e1d83109SMatthew G. Knepley         for (c = 0; c < Nc[p]; ++c) {
538e1d83109SMatthew G. Knepley           const PetscInt r = Noff[p]+c;
539e1d83109SMatthew G. Knepley           PetscInt       n;
540e1d83109SMatthew G. Knepley 
541e1d83109SMatthew G. Knepley           for (n = 0; n < N[r]; ++n, ++off) {
542e1d83109SMatthew G. Knepley             const PetscInt    q = Noff[adj[off].rank] + adj[off].index;
543e1d83109SMatthew G. Knepley             const PetscScalar o = val[off] ? 1.0 : 0.0;
54484961fc4SMatthew G. Knepley 
54584961fc4SMatthew G. Knepley             ierr = MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES);CHKERRQ(ierr);
54684961fc4SMatthew G. Knepley             ierr = MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES);CHKERRQ(ierr);
54784961fc4SMatthew G. Knepley           }
54884961fc4SMatthew G. Knepley         }
549e1d83109SMatthew G. Knepley       }
55084961fc4SMatthew G. Knepley       ierr = MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
55184961fc4SMatthew G. Knepley       ierr = MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
55284961fc4SMatthew G. Knepley 
5539852e123SBarry Smith       ierr = PetscBTCreate(Noff[size], &seenProcs);CHKERRQ(ierr);
5549852e123SBarry Smith       ierr = PetscBTMemzero(Noff[size], seenProcs);CHKERRQ(ierr);
5559852e123SBarry Smith       ierr = PetscBTCreate(Noff[size], &flippedProcs);CHKERRQ(ierr);
5569852e123SBarry Smith       ierr = PetscBTMemzero(Noff[size], flippedProcs);CHKERRQ(ierr);
5579852e123SBarry Smith       ierr = PetscMalloc1(Noff[size], &procFIFO);CHKERRQ(ierr);
55884961fc4SMatthew G. Knepley       pTop = pBottom = 0;
5599852e123SBarry Smith       for (p = 0; p < Noff[size]; ++p) {
56084961fc4SMatthew G. Knepley         if (PetscBTLookup(seenProcs, p)) continue;
56184961fc4SMatthew G. Knepley         /* Initialize FIFO with next proc */
56284961fc4SMatthew G. Knepley         procFIFO[pBottom++] = p;
56384961fc4SMatthew G. Knepley         ierr = PetscBTSet(seenProcs, p);CHKERRQ(ierr);
56484961fc4SMatthew G. Knepley         /* Consider each proc in FIFO */
56584961fc4SMatthew G. Knepley         while (pTop < pBottom) {
56684961fc4SMatthew G. Knepley           const PetscScalar *ornt;
56784961fc4SMatthew G. Knepley           const PetscInt    *neighbors;
568e1d83109SMatthew G. Knepley           PetscInt           proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n;
56984961fc4SMatthew G. Knepley 
57084961fc4SMatthew G. Knepley           proc     = procFIFO[pTop++];
57184961fc4SMatthew G. Knepley           flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0;
57284961fc4SMatthew G. Knepley           ierr = MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt);CHKERRQ(ierr);
57384961fc4SMatthew G. Knepley           /* Loop over neighboring procs */
57484961fc4SMatthew G. Knepley           for (n = 0; n < numNeighbors; ++n) {
57584961fc4SMatthew G. Knepley             nproc    = neighbors[n];
57684961fc4SMatthew G. Knepley             mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1;
57784961fc4SMatthew G. Knepley             seen     = PetscBTLookup(seenProcs, nproc);
57884961fc4SMatthew G. Knepley             flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0;
57984961fc4SMatthew G. Knepley 
58084961fc4SMatthew G. Knepley             if (mismatch ^ (flippedA ^ flippedB)) {
58184961fc4SMatthew 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);
58284961fc4SMatthew G. Knepley               if (!flippedB) {
58384961fc4SMatthew G. Knepley                 ierr = PetscBTSet(flippedProcs, nproc);CHKERRQ(ierr);
58484961fc4SMatthew G. Knepley               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
58584961fc4SMatthew 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");
58684961fc4SMatthew G. Knepley             if (!seen) {
58784961fc4SMatthew G. Knepley               procFIFO[pBottom++] = nproc;
58884961fc4SMatthew G. Knepley               ierr = PetscBTSet(seenProcs, nproc);CHKERRQ(ierr);
58984961fc4SMatthew G. Knepley             }
59084961fc4SMatthew G. Knepley           }
59184961fc4SMatthew G. Knepley         }
59284961fc4SMatthew G. Knepley       }
59384961fc4SMatthew G. Knepley       ierr = PetscFree(procFIFO);CHKERRQ(ierr);
59484961fc4SMatthew G. Knepley       ierr = MatDestroy(&G);CHKERRQ(ierr);
59584961fc4SMatthew G. Knepley       ierr = PetscFree2(adj, val);CHKERRQ(ierr);
596e1d83109SMatthew G. Knepley       ierr = PetscBTDestroy(&seenProcs);CHKERRQ(ierr);
59784961fc4SMatthew G. Knepley     }
598e1d83109SMatthew G. Knepley     /* Scatter flip flags */
599e1d83109SMatthew G. Knepley     {
600e1d83109SMatthew G. Knepley       PetscBool *flips = NULL;
601e1d83109SMatthew G. Knepley 
602e1d83109SMatthew G. Knepley       if (!rank) {
6039852e123SBarry Smith         ierr = PetscMalloc1(Noff[size], &flips);CHKERRQ(ierr);
6049852e123SBarry Smith         for (p = 0; p < Noff[size]; ++p) {
605e1d83109SMatthew G. Knepley           flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE;
606302440fdSBarry Smith           if (flg && flips[p]) {ierr = PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p);CHKERRQ(ierr);}
607e1d83109SMatthew G. Knepley         }
6089852e123SBarry Smith         for (p = 0; p < size; ++p) {
609e1d83109SMatthew G. Knepley           displs[p+1] = displs[p] + Nc[p];
610e1d83109SMatthew G. Knepley         }
611e1d83109SMatthew G. Knepley       }
612ffc4695bSBarry Smith       ierr = MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, numComponents, MPIU_BOOL, 0, comm);CHKERRMPI(ierr);
61384961fc4SMatthew G. Knepley       ierr = PetscFree(flips);CHKERRQ(ierr);
61484961fc4SMatthew G. Knepley     }
615e1d83109SMatthew G. Knepley     if (!rank) {ierr = PetscBTDestroy(&flippedProcs);CHKERRQ(ierr);}
616e1d83109SMatthew G. Knepley     ierr = PetscFree(N);CHKERRQ(ierr);
617e1d83109SMatthew G. Knepley     ierr = PetscFree4(recvcounts, displs, Nc, Noff);CHKERRQ(ierr);
618e1d83109SMatthew G. Knepley     ierr = PetscFree2(nrankComp, match);CHKERRQ(ierr);
619e1d83109SMatthew G. Knepley 
620e1d83109SMatthew G. Knepley     /* Decide whether to flip cells in each component */
621e1d83109SMatthew G. Knepley     for (c = 0; c < cEnd-cStart; ++c) {if (flipped[cellComp[c]]) {ierr = PetscBTNegate(flippedCells, c);CHKERRQ(ierr);}}
622e1d83109SMatthew G. Knepley     ierr = PetscFree(flipped);CHKERRQ(ierr);
62384961fc4SMatthew G. Knepley   }
62484961fc4SMatthew G. Knepley   if (flg) {
62584961fc4SMatthew G. Knepley     PetscViewer v;
62684961fc4SMatthew G. Knepley 
6277e09831bSMatthew G. Knepley     ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr);
6281575c14dSBarry Smith     ierr = PetscViewerASCIIPushSynchronized(v);CHKERRQ(ierr);
62984961fc4SMatthew G. Knepley     ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank);CHKERRQ(ierr);
63084961fc4SMatthew G. Knepley     ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr);
6314d4c343aSBarry Smith     ierr = PetscViewerFlush(v);CHKERRQ(ierr);
6321575c14dSBarry Smith     ierr = PetscViewerASCIIPopSynchronized(v);CHKERRQ(ierr);
63384961fc4SMatthew G. Knepley   }
63484961fc4SMatthew G. Knepley   /* Reverse flipped cells in the mesh */
63584961fc4SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
6360d366550SMatthew G. Knepley     if (PetscBTLookup(flippedCells, c-cStart)) {
6370d366550SMatthew G. Knepley       ierr = DMPlexReverseCell(dm, c);CHKERRQ(ierr);
6380d366550SMatthew G. Knepley     }
63984961fc4SMatthew G. Knepley   }
64084961fc4SMatthew G. Knepley   ierr = PetscBTDestroy(&seenCells);CHKERRQ(ierr);
64184961fc4SMatthew G. Knepley   ierr = PetscBTDestroy(&flippedCells);CHKERRQ(ierr);
64284961fc4SMatthew G. Knepley   ierr = PetscBTDestroy(&seenFaces);CHKERRQ(ierr);
643e1d83109SMatthew G. Knepley   ierr = PetscFree2(numNeighbors, neighbors);CHKERRQ(ierr);
644e1d83109SMatthew G. Knepley   ierr = PetscFree2(rorntComp, lorntComp);CHKERRQ(ierr);
645e1d83109SMatthew G. Knepley   ierr = PetscFree3(faceFIFO, cellComp, faceComp);CHKERRQ(ierr);
64684961fc4SMatthew G. Knepley   PetscFunctionReturn(0);
64784961fc4SMatthew G. Knepley }
648