xref: /petsc/src/dm/impls/plex/plexorient.c (revision 7e09831bc92b7f5435999991a3b282b8ca6d1e95)
184961fc4SMatthew G. Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
284961fc4SMatthew G. Knepley #include <petscsf.h>
384961fc4SMatthew G. Knepley 
484961fc4SMatthew G. Knepley #undef __FUNCT__
584961fc4SMatthew G. Knepley #define __FUNCT__ "DMPlexReverseCell"
684961fc4SMatthew G. Knepley /*@
784961fc4SMatthew G. Knepley   DMPlexReverseCell - Give a mesh cell the opposite orientation
884961fc4SMatthew G. Knepley 
984961fc4SMatthew G. Knepley   Input Parameters:
1084961fc4SMatthew G. Knepley + dm   - The DM
1184961fc4SMatthew G. Knepley - cell - The cell number
1284961fc4SMatthew G. Knepley 
1384961fc4SMatthew G. Knepley   Note: The modification of the DM is done in-place.
1484961fc4SMatthew G. Knepley 
1584961fc4SMatthew G. Knepley   Level: advanced
1684961fc4SMatthew G. Knepley 
1784961fc4SMatthew G. Knepley .seealso: DMPlexOrient(), DMCreate(), DMPLEX
1884961fc4SMatthew G. Knepley @*/
1984961fc4SMatthew G. Knepley PetscErrorCode DMPlexReverseCell(DM dm, PetscInt cell)
2084961fc4SMatthew G. Knepley {
2184961fc4SMatthew G. Knepley   /* Note that the reverse orientation ro of a face with orientation o is:
2284961fc4SMatthew G. Knepley 
2384961fc4SMatthew G. Knepley        ro = o >= 0 ? -(faceSize - o) : faceSize + o
2484961fc4SMatthew G. Knepley 
2584961fc4SMatthew G. Knepley      where faceSize is the size of the cone for the face.
2684961fc4SMatthew G. Knepley   */
2784961fc4SMatthew G. Knepley   const PetscInt *cone,    *coneO, *support;
2884961fc4SMatthew G. Knepley   PetscInt       *revcone, *revconeO;
2984961fc4SMatthew G. Knepley   PetscInt        maxConeSize, coneSize, supportSize, faceSize, cp, sp;
3084961fc4SMatthew G. Knepley   PetscErrorCode  ierr;
3184961fc4SMatthew G. Knepley 
3284961fc4SMatthew G. Knepley   PetscFunctionBegin;
3384961fc4SMatthew G. Knepley   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
3484961fc4SMatthew G. Knepley   ierr = DMGetWorkArray(dm, maxConeSize, PETSC_INT, &revcone);CHKERRQ(ierr);
3584961fc4SMatthew G. Knepley   ierr = DMGetWorkArray(dm, maxConeSize, PETSC_INT, &revconeO);CHKERRQ(ierr);
3684961fc4SMatthew G. Knepley   /* Reverse cone, and reverse orientations of faces */
3784961fc4SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
3884961fc4SMatthew G. Knepley   ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
3984961fc4SMatthew G. Knepley   ierr = DMPlexGetConeOrientation(dm, cell, &coneO);CHKERRQ(ierr);
4084961fc4SMatthew G. Knepley   for (cp = 0; cp < coneSize; ++cp) {
4184961fc4SMatthew G. Knepley     const PetscInt rcp = coneSize-cp-1;
4284961fc4SMatthew G. Knepley 
4384961fc4SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, cone[rcp], &faceSize);CHKERRQ(ierr);
4484961fc4SMatthew G. Knepley     revcone[cp]  = cone[rcp];
4584961fc4SMatthew G. Knepley     revconeO[cp] = coneO[rcp] >= 0 ? -(faceSize-coneO[rcp]) : faceSize+coneO[rcp];
4684961fc4SMatthew G. Knepley   }
4784961fc4SMatthew G. Knepley   ierr = DMPlexSetCone(dm, cell, revcone);CHKERRQ(ierr);
4884961fc4SMatthew G. Knepley   ierr = DMPlexSetConeOrientation(dm, cell, revconeO);CHKERRQ(ierr);
4984961fc4SMatthew G. Knepley   /* Reverse orientation of this cell in the support hypercells */
5084961fc4SMatthew G. Knepley   faceSize = coneSize;
5184961fc4SMatthew G. Knepley   ierr = DMPlexGetSupportSize(dm, cell, &supportSize);CHKERRQ(ierr);
5284961fc4SMatthew G. Knepley   ierr = DMPlexGetSupport(dm, cell, &support);CHKERRQ(ierr);
5384961fc4SMatthew G. Knepley   for (sp = 0; sp < supportSize; ++sp) {
5484961fc4SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, support[sp], &coneSize);CHKERRQ(ierr);
5584961fc4SMatthew G. Knepley     ierr = DMPlexGetCone(dm, support[sp], &cone);CHKERRQ(ierr);
5684961fc4SMatthew G. Knepley     ierr = DMPlexGetConeOrientation(dm, support[sp], &coneO);CHKERRQ(ierr);
5784961fc4SMatthew G. Knepley     for (cp = 0; cp < coneSize; ++cp) {
5884961fc4SMatthew G. Knepley       if (cone[cp] != cell) continue;
5984961fc4SMatthew G. Knepley       ierr = DMPlexInsertConeOrientation(dm, support[sp], cp, coneO[cp] >= 0 ? -(faceSize-coneO[cp]) : faceSize+coneO[cp]);CHKERRQ(ierr);
6084961fc4SMatthew G. Knepley     }
6184961fc4SMatthew G. Knepley   }
6284961fc4SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, maxConeSize, PETSC_INT, &revcone);CHKERRQ(ierr);
6384961fc4SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, maxConeSize, PETSC_INT, &revconeO);CHKERRQ(ierr);
6484961fc4SMatthew G. Knepley   PetscFunctionReturn(0);
6584961fc4SMatthew G. Knepley }
6684961fc4SMatthew G. Knepley 
6784961fc4SMatthew G. Knepley #undef __FUNCT__
6884961fc4SMatthew G. Knepley #define __FUNCT__ "DMPlexOrient"
6984961fc4SMatthew G. Knepley /*@
7084961fc4SMatthew G. Knepley   DMPlexOrient - Give a consistent orientation to the input mesh
7184961fc4SMatthew G. Knepley 
7284961fc4SMatthew G. Knepley   Input Parameters:
7384961fc4SMatthew G. Knepley . dm - The DM
7484961fc4SMatthew G. Knepley 
7584961fc4SMatthew G. Knepley   Note: The orientation data for the DM are change in-place.
7684961fc4SMatthew G. Knepley $ This routine will fail for non-orientable surfaces, such as the Moebius strip.
7784961fc4SMatthew G. Knepley 
7884961fc4SMatthew G. Knepley   Level: advanced
7984961fc4SMatthew G. Knepley 
8084961fc4SMatthew G. Knepley .seealso: DMCreate(), DMPLEX
8184961fc4SMatthew G. Knepley @*/
8284961fc4SMatthew G. Knepley PetscErrorCode DMPlexOrient(DM dm)
8384961fc4SMatthew G. Knepley {
8484961fc4SMatthew G. Knepley   MPI_Comm       comm;
8584961fc4SMatthew G. Knepley   PetscBT        seenCells, flippedCells, seenFaces;
8684961fc4SMatthew G. Knepley   PetscInt      *faceFIFO, fTop, fBottom;
8784961fc4SMatthew G. Knepley   PetscInt       dim, h, cStart, cEnd, c, fStart, fEnd, face;
8884961fc4SMatthew G. Knepley   PetscBool      flg;
8984961fc4SMatthew G. Knepley   PetscErrorCode ierr;
9084961fc4SMatthew G. Knepley 
9184961fc4SMatthew G. Knepley   PetscFunctionBegin;
9284961fc4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
9384961fc4SMatthew G. Knepley   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-orientation_view", &flg);CHKERRQ(ierr);
9484961fc4SMatthew G. Knepley   /* Truth Table
9584961fc4SMatthew G. Knepley      mismatch    flips   do action   mismatch   flipA ^ flipB   action
9684961fc4SMatthew G. Knepley          F       0 flips     no         F             F           F
9784961fc4SMatthew G. Knepley          F       1 flip      yes        F             T           T
9884961fc4SMatthew G. Knepley          F       2 flips     no         T             F           T
9984961fc4SMatthew G. Knepley          T       0 flips     yes        T             T           F
10084961fc4SMatthew G. Knepley          T       1 flip      no
10184961fc4SMatthew G. Knepley          T       2 flips     yes
10284961fc4SMatthew G. Knepley   */
10384961fc4SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
10484961fc4SMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(dm, &h);CHKERRQ(ierr);
10584961fc4SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, h,   &cStart, &cEnd);CHKERRQ(ierr);
10684961fc4SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, h+1, &fStart, &fEnd);CHKERRQ(ierr);
10784961fc4SMatthew G. Knepley   ierr = PetscBTCreate(cEnd - cStart, &seenCells);CHKERRQ(ierr);
10884961fc4SMatthew G. Knepley   ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr);
10984961fc4SMatthew G. Knepley   ierr = PetscBTCreate(cEnd - cStart, &flippedCells);CHKERRQ(ierr);
11084961fc4SMatthew G. Knepley   ierr = PetscBTMemzero(cEnd - cStart, flippedCells);CHKERRQ(ierr);
11184961fc4SMatthew G. Knepley   ierr = PetscBTCreate(fEnd - fStart, &seenFaces);CHKERRQ(ierr);
11284961fc4SMatthew G. Knepley   ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr);
11384961fc4SMatthew G. Knepley   ierr = PetscMalloc1((fEnd - fStart), &faceFIFO);CHKERRQ(ierr);
11484961fc4SMatthew G. Knepley   fTop = fBottom = 0;
11584961fc4SMatthew G. Knepley   /* Initialize FIFO with first cell */
11684961fc4SMatthew G. Knepley   if (cEnd > cStart) {
11784961fc4SMatthew G. Knepley     const PetscInt *cone;
11884961fc4SMatthew G. Knepley     PetscInt        coneSize;
11984961fc4SMatthew G. Knepley 
12084961fc4SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, cStart, &coneSize);CHKERRQ(ierr);
12184961fc4SMatthew G. Knepley     ierr = DMPlexGetCone(dm, cStart, &cone);CHKERRQ(ierr);
12284961fc4SMatthew G. Knepley     for (c = 0; c < coneSize; ++c) {
12384961fc4SMatthew G. Knepley       faceFIFO[fBottom++] = cone[c];
12484961fc4SMatthew G. Knepley       ierr = PetscBTSet(seenFaces, cone[c]-fStart);CHKERRQ(ierr);
12584961fc4SMatthew G. Knepley     }
12684961fc4SMatthew G. Knepley   }
12784961fc4SMatthew G. Knepley   /* Consider each face in FIFO */
12884961fc4SMatthew G. Knepley   while (fTop < fBottom) {
12984961fc4SMatthew G. Knepley     const PetscInt *support, *coneA, *coneB, *coneOA, *coneOB;
13084961fc4SMatthew G. Knepley     PetscInt        supportSize, coneSizeA, coneSizeB, posA = -1, posB = -1;
13184961fc4SMatthew G. Knepley     PetscInt        seenA, flippedA, seenB, flippedB, mismatch;
13284961fc4SMatthew G. Knepley 
13384961fc4SMatthew G. Knepley     face = faceFIFO[fTop++];
13484961fc4SMatthew G. Knepley     ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
13584961fc4SMatthew G. Knepley     ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr);
13684961fc4SMatthew G. Knepley     if (supportSize < 2) continue;
13784961fc4SMatthew G. Knepley     if (supportSize != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %d", supportSize);
13884961fc4SMatthew G. Knepley     seenA    = PetscBTLookup(seenCells,    support[0]-cStart);
13984961fc4SMatthew G. Knepley     flippedA = PetscBTLookup(flippedCells, support[0]-cStart) ? 1 : 0;
14084961fc4SMatthew G. Knepley     seenB    = PetscBTLookup(seenCells,    support[1]-cStart);
14184961fc4SMatthew G. Knepley     flippedB = PetscBTLookup(flippedCells, support[1]-cStart) ? 1 : 0;
14284961fc4SMatthew G. Knepley 
14384961fc4SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, support[0], &coneSizeA);CHKERRQ(ierr);
14484961fc4SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, support[1], &coneSizeB);CHKERRQ(ierr);
14584961fc4SMatthew G. Knepley     ierr = DMPlexGetCone(dm, support[0], &coneA);CHKERRQ(ierr);
14684961fc4SMatthew G. Knepley     ierr = DMPlexGetCone(dm, support[1], &coneB);CHKERRQ(ierr);
14784961fc4SMatthew G. Knepley     ierr = DMPlexGetConeOrientation(dm, support[0], &coneOA);CHKERRQ(ierr);
14884961fc4SMatthew G. Knepley     ierr = DMPlexGetConeOrientation(dm, support[1], &coneOB);CHKERRQ(ierr);
14984961fc4SMatthew G. Knepley     for (c = 0; c < coneSizeA; ++c) {
15084961fc4SMatthew G. Knepley       if (!PetscBTLookup(seenFaces, coneA[c]-fStart)) {
15184961fc4SMatthew G. Knepley         faceFIFO[fBottom++] = coneA[c];
15284961fc4SMatthew G. Knepley         ierr = PetscBTSet(seenFaces, coneA[c]-fStart);CHKERRQ(ierr);
15384961fc4SMatthew G. Knepley       }
15484961fc4SMatthew G. Knepley       if (coneA[c] == face) posA = c;
15584961fc4SMatthew 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);
15684961fc4SMatthew G. Knepley     }
15784961fc4SMatthew 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]);
15884961fc4SMatthew G. Knepley     for (c = 0; c < coneSizeB; ++c) {
15984961fc4SMatthew G. Knepley       if (!PetscBTLookup(seenFaces, coneB[c]-fStart)) {
16084961fc4SMatthew G. Knepley         faceFIFO[fBottom++] = coneB[c];
16184961fc4SMatthew G. Knepley         ierr = PetscBTSet(seenFaces, coneB[c]-fStart);CHKERRQ(ierr);
16284961fc4SMatthew G. Knepley       }
16384961fc4SMatthew G. Knepley       if (coneB[c] == face) posB = c;
16484961fc4SMatthew 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);
16584961fc4SMatthew G. Knepley     }
16684961fc4SMatthew 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]);
16784961fc4SMatthew G. Knepley 
16884961fc4SMatthew G. Knepley     if (dim == 1) {
16984961fc4SMatthew G. Knepley       mismatch = posA == posB;
17084961fc4SMatthew G. Knepley     } else {
17184961fc4SMatthew G. Knepley       mismatch = coneOA[posA] == coneOB[posB];
17284961fc4SMatthew G. Knepley     }
17384961fc4SMatthew G. Knepley 
17484961fc4SMatthew G. Knepley     if (mismatch ^ (flippedA ^ flippedB)) {
17584961fc4SMatthew 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]);
17684961fc4SMatthew G. Knepley       if (!seenA && !flippedA) {
17784961fc4SMatthew G. Knepley         ierr = PetscBTSet(flippedCells, support[0]-cStart);CHKERRQ(ierr);
17884961fc4SMatthew G. Knepley       } else if (!seenB && !flippedB) {
17984961fc4SMatthew G. Knepley         ierr = PetscBTSet(flippedCells, support[1]-cStart);CHKERRQ(ierr);
18084961fc4SMatthew G. Knepley       } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
18184961fc4SMatthew 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");
18284961fc4SMatthew G. Knepley     ierr = PetscBTSet(seenCells, support[0]-cStart);CHKERRQ(ierr);
18384961fc4SMatthew G. Knepley     ierr = PetscBTSet(seenCells, support[1]-cStart);CHKERRQ(ierr);
18484961fc4SMatthew G. Knepley   }
18584961fc4SMatthew G. Knepley   if (flg) {
18684961fc4SMatthew G. Knepley     PetscViewer v;
18784961fc4SMatthew G. Knepley     PetscMPIInt rank;
18884961fc4SMatthew G. Knepley 
18984961fc4SMatthew G. Knepley     ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
190*7e09831bSMatthew G. Knepley     ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr);
191*7e09831bSMatthew G. Knepley     ierr = PetscViewerASCIISynchronizedAllow(v, PETSC_TRUE);CHKERRQ(ierr);
192*7e09831bSMatthew G. Knepley     ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial seen faces:\n", rank);CHKERRQ(ierr);
193*7e09831bSMatthew G. Knepley     ierr = PetscBTView(fEnd-fStart, seenFaces,    v);CHKERRQ(ierr);
19484961fc4SMatthew G. Knepley     ierr = PetscViewerASCIISynchronizedAllow(v, PETSC_TRUE);CHKERRQ(ierr);
19584961fc4SMatthew G. Knepley     ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank);CHKERRQ(ierr);
19684961fc4SMatthew G. Knepley     ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr);
19784961fc4SMatthew G. Knepley   }
19884961fc4SMatthew G. Knepley   /* Now all subdomains are oriented, but we need a consistent parallel orientation */
19984961fc4SMatthew G. Knepley   {
20084961fc4SMatthew G. Knepley     /* Find a representative face (edge) separating pairs of procs */
20184961fc4SMatthew G. Knepley     PetscSF            sf;
20284961fc4SMatthew G. Knepley     const PetscInt    *lpoints;
20384961fc4SMatthew G. Knepley     const PetscSFNode *rpoints;
20484961fc4SMatthew G. Knepley     PetscInt          *neighbors, *nranks;
20584961fc4SMatthew G. Knepley     PetscInt           numLeaves, numRoots, numNeighbors = 0, l, n;
20684961fc4SMatthew G. Knepley 
20784961fc4SMatthew G. Knepley     ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
20884961fc4SMatthew G. Knepley     ierr = PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints);CHKERRQ(ierr);
20984961fc4SMatthew G. Knepley     if (numLeaves >= 0) {
21084961fc4SMatthew G. Knepley       const PetscInt *cone, *ornt, *support;
21184961fc4SMatthew G. Knepley       PetscInt        coneSize, supportSize;
21284961fc4SMatthew G. Knepley       int            *rornt, *lornt; /* PetscSF cannot handle smaller than int */
21384961fc4SMatthew G. Knepley       PetscBool      *match, flipped = PETSC_FALSE;
21484961fc4SMatthew G. Knepley 
21584961fc4SMatthew G. Knepley       ierr = PetscMalloc1(numLeaves,&neighbors);CHKERRQ(ierr);
21684961fc4SMatthew G. Knepley       /* I know this is p^2 time in general, but for bounded degree its alright */
21784961fc4SMatthew G. Knepley       for (l = 0; l < numLeaves; ++l) {
21884961fc4SMatthew G. Knepley         const PetscInt face = lpoints[l];
21984961fc4SMatthew G. Knepley         if ((face >= fStart) && (face < fEnd)) {
22084961fc4SMatthew G. Knepley           const PetscInt rank = rpoints[l].rank;
22184961fc4SMatthew G. Knepley           for (n = 0; n < numNeighbors; ++n) if (rank == rpoints[neighbors[n]].rank) break;
22284961fc4SMatthew G. Knepley           if (n >= numNeighbors) {
22384961fc4SMatthew G. Knepley             PetscInt supportSize;
22484961fc4SMatthew G. Knepley             ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
22584961fc4SMatthew G. Knepley             if (supportSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %d", supportSize);
22684961fc4SMatthew G. Knepley             neighbors[numNeighbors++] = l;
22784961fc4SMatthew G. Knepley           }
22884961fc4SMatthew G. Knepley         }
22984961fc4SMatthew G. Knepley       }
23084961fc4SMatthew G. Knepley       ierr = PetscCalloc4(numNeighbors,&match,numNeighbors,&nranks,numRoots,&rornt,numRoots,&lornt);CHKERRQ(ierr);
23184961fc4SMatthew G. Knepley       for (face = fStart; face < fEnd; ++face) {
23284961fc4SMatthew G. Knepley         ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
23384961fc4SMatthew G. Knepley         if (supportSize != 1) continue;
23484961fc4SMatthew G. Knepley         ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr);
23584961fc4SMatthew G. Knepley 
23684961fc4SMatthew G. Knepley         ierr = DMPlexGetCone(dm, support[0], &cone);CHKERRQ(ierr);
23784961fc4SMatthew G. Knepley         ierr = DMPlexGetConeSize(dm, support[0], &coneSize);CHKERRQ(ierr);
23884961fc4SMatthew G. Knepley         ierr = DMPlexGetConeOrientation(dm, support[0], &ornt);CHKERRQ(ierr);
23984961fc4SMatthew G. Knepley         for (c = 0; c < coneSize; ++c) if (cone[c] == face) break;
24084961fc4SMatthew G. Knepley         if (dim == 1) {
24184961fc4SMatthew G. Knepley           /* Use cone position instead, shifted to -1 or 1 */
24284961fc4SMatthew G. Knepley           rornt[face] = c*2-1;
24384961fc4SMatthew G. Knepley         } else {
24484961fc4SMatthew G. Knepley           if (PetscBTLookup(flippedCells, support[0]-cStart)) rornt[face] = ornt[c] < 0 ? -1 :  1;
24584961fc4SMatthew G. Knepley           else                                                rornt[face] = ornt[c] < 0 ?  1 : -1;
24684961fc4SMatthew G. Knepley         }
24784961fc4SMatthew G. Knepley       }
24884961fc4SMatthew G. Knepley       /* Mark each edge with match or nomatch */
24984961fc4SMatthew G. Knepley       ierr = PetscSFBcastBegin(sf, MPI_INT, rornt, lornt);CHKERRQ(ierr);
25084961fc4SMatthew G. Knepley       ierr = PetscSFBcastEnd(sf, MPI_INT, rornt, lornt);CHKERRQ(ierr);
25184961fc4SMatthew G. Knepley       for (n = 0; n < numNeighbors; ++n) {
25284961fc4SMatthew G. Knepley         const PetscInt face = lpoints[neighbors[n]];
25384961fc4SMatthew G. Knepley 
25484961fc4SMatthew G. Knepley         if (rornt[face]*lornt[face] < 0) match[n] = PETSC_TRUE;
25584961fc4SMatthew G. Knepley         else                             match[n] = PETSC_FALSE;
25684961fc4SMatthew G. Knepley         nranks[n] = rpoints[neighbors[n]].rank;
25784961fc4SMatthew G. Knepley       }
25884961fc4SMatthew G. Knepley       /* Collect the graph on 0 */
25984961fc4SMatthew G. Knepley       {
26084961fc4SMatthew G. Knepley         Mat          G;
26184961fc4SMatthew G. Knepley         PetscBT      seenProcs, flippedProcs;
26284961fc4SMatthew G. Knepley         PetscInt    *procFIFO, pTop, pBottom;
26384961fc4SMatthew G. Knepley         PetscInt    *adj = NULL;
26484961fc4SMatthew G. Knepley         PetscBool   *val = NULL;
26584961fc4SMatthew G. Knepley         PetscMPIInt *recvcounts = NULL, *displs = NULL, p;
26684961fc4SMatthew G. Knepley         PetscMPIInt  N = numNeighbors, numProcs = 0, rank;
26784961fc4SMatthew G. Knepley         PetscInt     debug = 0;
26884961fc4SMatthew G. Knepley 
26984961fc4SMatthew G. Knepley         ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
27084961fc4SMatthew G. Knepley         if (!rank) {ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);}
27184961fc4SMatthew G. Knepley         ierr = PetscCalloc2(numProcs,&recvcounts,numProcs+1,&displs);CHKERRQ(ierr);
27284961fc4SMatthew G. Knepley         ierr = MPI_Gather(&N, 1, MPI_INT, recvcounts, 1, MPI_INT, 0, comm);CHKERRQ(ierr);
27384961fc4SMatthew G. Knepley         for (p = 0; p < numProcs; ++p) {
27484961fc4SMatthew G. Knepley           displs[p+1] = displs[p] + recvcounts[p];
27584961fc4SMatthew G. Knepley         }
27684961fc4SMatthew G. Knepley         if (!rank) {ierr = PetscMalloc2(displs[numProcs],&adj,displs[numProcs],&val);CHKERRQ(ierr);}
27784961fc4SMatthew G. Knepley         ierr = MPI_Gatherv(nranks, numNeighbors, MPIU_INT, adj, recvcounts, displs, MPIU_INT, 0, comm);CHKERRQ(ierr);
27884961fc4SMatthew G. Knepley         ierr = MPI_Gatherv(match, numNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
27984961fc4SMatthew G. Knepley         if (debug) {
28084961fc4SMatthew G. Knepley           for (p = 0; p < numProcs; ++p) {
28184961fc4SMatthew G. Knepley             ierr = PetscPrintf(comm, "Proc %d:\n", p);
28284961fc4SMatthew G. Knepley             for (n = 0; n < recvcounts[p]; ++n) {
28384961fc4SMatthew G. Knepley               ierr = PetscPrintf(comm, "  edge %d (%d):\n", adj[displs[p]+n], val[displs[p]+n]);
28484961fc4SMatthew G. Knepley             }
28584961fc4SMatthew G. Knepley           }
28684961fc4SMatthew G. Knepley         }
28784961fc4SMatthew G. Knepley         /* Symmetrize the graph */
28884961fc4SMatthew G. Knepley         ierr = MatCreate(PETSC_COMM_SELF, &G);CHKERRQ(ierr);
28984961fc4SMatthew G. Knepley         ierr = MatSetSizes(G, numProcs, numProcs, numProcs, numProcs);CHKERRQ(ierr);
29084961fc4SMatthew G. Knepley         ierr = MatSetUp(G);CHKERRQ(ierr);
29184961fc4SMatthew G. Knepley         for (p = 0; p < numProcs; ++p) {
29284961fc4SMatthew G. Knepley           for (n = 0; n < recvcounts[p]; ++n) {
29384961fc4SMatthew G. Knepley             const PetscInt    r = p;
29484961fc4SMatthew G. Knepley             const PetscInt    q = adj[displs[p]+n];
29584961fc4SMatthew G. Knepley             const PetscScalar o = val[displs[p]+n] ? 1.0 : 0.0;
29684961fc4SMatthew G. Knepley 
29784961fc4SMatthew G. Knepley             ierr = MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES);CHKERRQ(ierr);
29884961fc4SMatthew G. Knepley             ierr = MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES);CHKERRQ(ierr);
29984961fc4SMatthew G. Knepley           }
30084961fc4SMatthew G. Knepley         }
30184961fc4SMatthew G. Knepley         ierr = MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
30284961fc4SMatthew G. Knepley         ierr = MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
30384961fc4SMatthew G. Knepley 
30484961fc4SMatthew G. Knepley         ierr = PetscBTCreate(numProcs, &seenProcs);CHKERRQ(ierr);
30584961fc4SMatthew G. Knepley         ierr = PetscBTMemzero(numProcs, seenProcs);CHKERRQ(ierr);
30684961fc4SMatthew G. Knepley         ierr = PetscBTCreate(numProcs, &flippedProcs);CHKERRQ(ierr);
30784961fc4SMatthew G. Knepley         ierr = PetscBTMemzero(numProcs, flippedProcs);CHKERRQ(ierr);
30884961fc4SMatthew G. Knepley         ierr = PetscMalloc1(numProcs,&procFIFO);CHKERRQ(ierr);
30984961fc4SMatthew G. Knepley         pTop = pBottom = 0;
31084961fc4SMatthew G. Knepley         for (p = 0; p < numProcs; ++p) {
31184961fc4SMatthew G. Knepley           if (PetscBTLookup(seenProcs, p)) continue;
31284961fc4SMatthew G. Knepley           /* Initialize FIFO with next proc */
31384961fc4SMatthew G. Knepley           procFIFO[pBottom++] = p;
31484961fc4SMatthew G. Knepley           ierr = PetscBTSet(seenProcs, p);CHKERRQ(ierr);
31584961fc4SMatthew G. Knepley           /* Consider each proc in FIFO */
31684961fc4SMatthew G. Knepley           while (pTop < pBottom) {
31784961fc4SMatthew G. Knepley             const PetscScalar *ornt;
31884961fc4SMatthew G. Knepley             const PetscInt    *neighbors;
31984961fc4SMatthew G. Knepley             PetscInt           proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors;
32084961fc4SMatthew G. Knepley 
32184961fc4SMatthew G. Knepley             proc     = procFIFO[pTop++];
32284961fc4SMatthew G. Knepley             flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0;
32384961fc4SMatthew G. Knepley             ierr = MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt);CHKERRQ(ierr);
32484961fc4SMatthew G. Knepley             /* Loop over neighboring procs */
32584961fc4SMatthew G. Knepley             for (n = 0; n < numNeighbors; ++n) {
32684961fc4SMatthew G. Knepley               nproc    = neighbors[n];
32784961fc4SMatthew G. Knepley               mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1;
32884961fc4SMatthew G. Knepley               seen     = PetscBTLookup(seenProcs, nproc);
32984961fc4SMatthew G. Knepley               flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0;
33084961fc4SMatthew G. Knepley 
33184961fc4SMatthew G. Knepley               if (mismatch ^ (flippedA ^ flippedB)) {
33284961fc4SMatthew 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);
33384961fc4SMatthew G. Knepley                 if (!flippedB) {
33484961fc4SMatthew G. Knepley                   ierr = PetscBTSet(flippedProcs, nproc);CHKERRQ(ierr);
33584961fc4SMatthew G. Knepley               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
33684961fc4SMatthew 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");
33784961fc4SMatthew G. Knepley               if (!seen) {
33884961fc4SMatthew G. Knepley                 procFIFO[pBottom++] = nproc;
33984961fc4SMatthew G. Knepley                 ierr = PetscBTSet(seenProcs, nproc);CHKERRQ(ierr);
34084961fc4SMatthew G. Knepley               }
34184961fc4SMatthew G. Knepley             }
34284961fc4SMatthew G. Knepley           }
34384961fc4SMatthew G. Knepley         }
34484961fc4SMatthew G. Knepley         ierr = PetscFree(procFIFO);CHKERRQ(ierr);
34584961fc4SMatthew G. Knepley         ierr = MatDestroy(&G);CHKERRQ(ierr);
34684961fc4SMatthew G. Knepley 
34784961fc4SMatthew G. Knepley         ierr = PetscFree2(recvcounts,displs);CHKERRQ(ierr);
34884961fc4SMatthew G. Knepley         ierr = PetscFree2(adj,val);CHKERRQ(ierr);
34984961fc4SMatthew G. Knepley         {
35084961fc4SMatthew G. Knepley           PetscBool *flips;
35184961fc4SMatthew G. Knepley 
35284961fc4SMatthew G. Knepley           ierr = PetscMalloc1(numProcs,&flips);CHKERRQ(ierr);
35384961fc4SMatthew G. Knepley           for (p = 0; p < numProcs; ++p) {
35484961fc4SMatthew G. Knepley             flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE;
35584961fc4SMatthew G. Knepley             if (debug && flips[p]) {ierr = PetscPrintf(comm, "Flipping Proc %d:\n", p);}
35684961fc4SMatthew G. Knepley           }
35784961fc4SMatthew G. Knepley           ierr = MPI_Scatter(flips, 1, MPIU_BOOL, &flipped, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
35884961fc4SMatthew G. Knepley           ierr = PetscFree(flips);CHKERRQ(ierr);
35984961fc4SMatthew G. Knepley         }
36084961fc4SMatthew G. Knepley         ierr = PetscBTDestroy(&seenProcs);CHKERRQ(ierr);
36184961fc4SMatthew G. Knepley         ierr = PetscBTDestroy(&flippedProcs);CHKERRQ(ierr);
36284961fc4SMatthew G. Knepley       }
36384961fc4SMatthew G. Knepley       ierr = PetscFree4(match,nranks,rornt,lornt);CHKERRQ(ierr);
36484961fc4SMatthew G. Knepley       ierr = PetscFree(neighbors);CHKERRQ(ierr);
36584961fc4SMatthew G. Knepley       if (flipped) {for (c = cStart; c < cEnd; ++c) {ierr = PetscBTNegate(flippedCells, c-cStart);CHKERRQ(ierr);}}
36684961fc4SMatthew G. Knepley     }
36784961fc4SMatthew G. Knepley   }
36884961fc4SMatthew G. Knepley   if (flg) {
36984961fc4SMatthew G. Knepley     PetscViewer v;
37084961fc4SMatthew G. Knepley     PetscMPIInt rank;
37184961fc4SMatthew G. Knepley 
37284961fc4SMatthew G. Knepley     ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
373*7e09831bSMatthew G. Knepley     ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr);
37484961fc4SMatthew G. Knepley     ierr = PetscViewerASCIISynchronizedAllow(v, PETSC_TRUE);CHKERRQ(ierr);
37584961fc4SMatthew G. Knepley     ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank);CHKERRQ(ierr);
37684961fc4SMatthew G. Knepley     ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr);
37784961fc4SMatthew G. Knepley   }
37884961fc4SMatthew G. Knepley   /* Reverse flipped cells in the mesh */
37984961fc4SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
38084961fc4SMatthew G. Knepley     if (PetscBTLookup(flippedCells, c-cStart)) {ierr = DMPlexReverseCell(dm, c);CHKERRQ(ierr);}
38184961fc4SMatthew G. Knepley   }
38284961fc4SMatthew G. Knepley   ierr = PetscBTDestroy(&seenCells);CHKERRQ(ierr);
38384961fc4SMatthew G. Knepley   ierr = PetscBTDestroy(&flippedCells);CHKERRQ(ierr);
38484961fc4SMatthew G. Knepley   ierr = PetscBTDestroy(&seenFaces);CHKERRQ(ierr);
38584961fc4SMatthew G. Knepley   ierr = PetscFree(faceFIFO);CHKERRQ(ierr);
38684961fc4SMatthew G. Knepley   PetscFunctionReturn(0);
38784961fc4SMatthew G. Knepley }
388