1*84961fc4SMatthew G. Knepley #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2*84961fc4SMatthew G. Knepley #include <petscsf.h> 3*84961fc4SMatthew G. Knepley 4*84961fc4SMatthew G. Knepley #undef __FUNCT__ 5*84961fc4SMatthew G. Knepley #define __FUNCT__ "DMPlexReverseCell" 6*84961fc4SMatthew G. Knepley /*@ 7*84961fc4SMatthew G. Knepley DMPlexReverseCell - Give a mesh cell the opposite orientation 8*84961fc4SMatthew G. Knepley 9*84961fc4SMatthew G. Knepley Input Parameters: 10*84961fc4SMatthew G. Knepley + dm - The DM 11*84961fc4SMatthew G. Knepley - cell - The cell number 12*84961fc4SMatthew G. Knepley 13*84961fc4SMatthew G. Knepley Note: The modification of the DM is done in-place. 14*84961fc4SMatthew G. Knepley 15*84961fc4SMatthew G. Knepley Level: advanced 16*84961fc4SMatthew G. Knepley 17*84961fc4SMatthew G. Knepley .seealso: DMPlexOrient(), DMCreate(), DMPLEX 18*84961fc4SMatthew G. Knepley @*/ 19*84961fc4SMatthew G. Knepley PetscErrorCode DMPlexReverseCell(DM dm, PetscInt cell) 20*84961fc4SMatthew G. Knepley { 21*84961fc4SMatthew G. Knepley /* Note that the reverse orientation ro of a face with orientation o is: 22*84961fc4SMatthew G. Knepley 23*84961fc4SMatthew G. Knepley ro = o >= 0 ? -(faceSize - o) : faceSize + o 24*84961fc4SMatthew G. Knepley 25*84961fc4SMatthew G. Knepley where faceSize is the size of the cone for the face. 26*84961fc4SMatthew G. Knepley */ 27*84961fc4SMatthew G. Knepley const PetscInt *cone, *coneO, *support; 28*84961fc4SMatthew G. Knepley PetscInt *revcone, *revconeO; 29*84961fc4SMatthew G. Knepley PetscInt maxConeSize, coneSize, supportSize, faceSize, cp, sp; 30*84961fc4SMatthew G. Knepley PetscErrorCode ierr; 31*84961fc4SMatthew G. Knepley 32*84961fc4SMatthew G. Knepley PetscFunctionBegin; 33*84961fc4SMatthew G. Knepley ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr); 34*84961fc4SMatthew G. Knepley ierr = DMGetWorkArray(dm, maxConeSize, PETSC_INT, &revcone);CHKERRQ(ierr); 35*84961fc4SMatthew G. Knepley ierr = DMGetWorkArray(dm, maxConeSize, PETSC_INT, &revconeO);CHKERRQ(ierr); 36*84961fc4SMatthew G. Knepley /* Reverse cone, and reverse orientations of faces */ 37*84961fc4SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 38*84961fc4SMatthew G. Knepley ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr); 39*84961fc4SMatthew G. Knepley ierr = DMPlexGetConeOrientation(dm, cell, &coneO);CHKERRQ(ierr); 40*84961fc4SMatthew G. Knepley for (cp = 0; cp < coneSize; ++cp) { 41*84961fc4SMatthew G. Knepley const PetscInt rcp = coneSize-cp-1; 42*84961fc4SMatthew G. Knepley 43*84961fc4SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cone[rcp], &faceSize);CHKERRQ(ierr); 44*84961fc4SMatthew G. Knepley revcone[cp] = cone[rcp]; 45*84961fc4SMatthew G. Knepley revconeO[cp] = coneO[rcp] >= 0 ? -(faceSize-coneO[rcp]) : faceSize+coneO[rcp]; 46*84961fc4SMatthew G. Knepley } 47*84961fc4SMatthew G. Knepley ierr = DMPlexSetCone(dm, cell, revcone);CHKERRQ(ierr); 48*84961fc4SMatthew G. Knepley ierr = DMPlexSetConeOrientation(dm, cell, revconeO);CHKERRQ(ierr); 49*84961fc4SMatthew G. Knepley /* Reverse orientation of this cell in the support hypercells */ 50*84961fc4SMatthew G. Knepley faceSize = coneSize; 51*84961fc4SMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, cell, &supportSize);CHKERRQ(ierr); 52*84961fc4SMatthew G. Knepley ierr = DMPlexGetSupport(dm, cell, &support);CHKERRQ(ierr); 53*84961fc4SMatthew G. Knepley for (sp = 0; sp < supportSize; ++sp) { 54*84961fc4SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, support[sp], &coneSize);CHKERRQ(ierr); 55*84961fc4SMatthew G. Knepley ierr = DMPlexGetCone(dm, support[sp], &cone);CHKERRQ(ierr); 56*84961fc4SMatthew G. Knepley ierr = DMPlexGetConeOrientation(dm, support[sp], &coneO);CHKERRQ(ierr); 57*84961fc4SMatthew G. Knepley for (cp = 0; cp < coneSize; ++cp) { 58*84961fc4SMatthew G. Knepley if (cone[cp] != cell) continue; 59*84961fc4SMatthew G. Knepley ierr = DMPlexInsertConeOrientation(dm, support[sp], cp, coneO[cp] >= 0 ? -(faceSize-coneO[cp]) : faceSize+coneO[cp]);CHKERRQ(ierr); 60*84961fc4SMatthew G. Knepley } 61*84961fc4SMatthew G. Knepley } 62*84961fc4SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, maxConeSize, PETSC_INT, &revcone);CHKERRQ(ierr); 63*84961fc4SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, maxConeSize, PETSC_INT, &revconeO);CHKERRQ(ierr); 64*84961fc4SMatthew G. Knepley PetscFunctionReturn(0); 65*84961fc4SMatthew G. Knepley } 66*84961fc4SMatthew G. Knepley 67*84961fc4SMatthew G. Knepley #undef __FUNCT__ 68*84961fc4SMatthew G. Knepley #define __FUNCT__ "DMPlexOrient" 69*84961fc4SMatthew G. Knepley /*@ 70*84961fc4SMatthew G. Knepley DMPlexOrient - Give a consistent orientation to the input mesh 71*84961fc4SMatthew G. Knepley 72*84961fc4SMatthew G. Knepley Input Parameters: 73*84961fc4SMatthew G. Knepley . dm - The DM 74*84961fc4SMatthew G. Knepley 75*84961fc4SMatthew G. Knepley Note: The orientation data for the DM are change in-place. 76*84961fc4SMatthew G. Knepley $ This routine will fail for non-orientable surfaces, such as the Moebius strip. 77*84961fc4SMatthew G. Knepley 78*84961fc4SMatthew G. Knepley Level: advanced 79*84961fc4SMatthew G. Knepley 80*84961fc4SMatthew G. Knepley .seealso: DMCreate(), DMPLEX 81*84961fc4SMatthew G. Knepley @*/ 82*84961fc4SMatthew G. Knepley PetscErrorCode DMPlexOrient(DM dm) 83*84961fc4SMatthew G. Knepley { 84*84961fc4SMatthew G. Knepley MPI_Comm comm; 85*84961fc4SMatthew G. Knepley PetscBT seenCells, flippedCells, seenFaces; 86*84961fc4SMatthew G. Knepley PetscInt *faceFIFO, fTop, fBottom; 87*84961fc4SMatthew G. Knepley PetscInt dim, h, cStart, cEnd, c, fStart, fEnd, face; 88*84961fc4SMatthew G. Knepley PetscBool flg; 89*84961fc4SMatthew G. Knepley PetscErrorCode ierr; 90*84961fc4SMatthew G. Knepley 91*84961fc4SMatthew G. Knepley PetscFunctionBegin; 92*84961fc4SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 93*84961fc4SMatthew G. Knepley ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-orientation_view", &flg);CHKERRQ(ierr); 94*84961fc4SMatthew G. Knepley /* Truth Table 95*84961fc4SMatthew G. Knepley mismatch flips do action mismatch flipA ^ flipB action 96*84961fc4SMatthew G. Knepley F 0 flips no F F F 97*84961fc4SMatthew G. Knepley F 1 flip yes F T T 98*84961fc4SMatthew G. Knepley F 2 flips no T F T 99*84961fc4SMatthew G. Knepley T 0 flips yes T T F 100*84961fc4SMatthew G. Knepley T 1 flip no 101*84961fc4SMatthew G. Knepley T 2 flips yes 102*84961fc4SMatthew G. Knepley */ 103*84961fc4SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 104*84961fc4SMatthew G. Knepley ierr = DMPlexGetVTKCellHeight(dm, &h);CHKERRQ(ierr); 105*84961fc4SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, h, &cStart, &cEnd);CHKERRQ(ierr); 106*84961fc4SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, h+1, &fStart, &fEnd);CHKERRQ(ierr); 107*84961fc4SMatthew G. Knepley ierr = PetscBTCreate(cEnd - cStart, &seenCells);CHKERRQ(ierr); 108*84961fc4SMatthew G. Knepley ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr); 109*84961fc4SMatthew G. Knepley ierr = PetscBTCreate(cEnd - cStart, &flippedCells);CHKERRQ(ierr); 110*84961fc4SMatthew G. Knepley ierr = PetscBTMemzero(cEnd - cStart, flippedCells);CHKERRQ(ierr); 111*84961fc4SMatthew G. Knepley ierr = PetscBTCreate(fEnd - fStart, &seenFaces);CHKERRQ(ierr); 112*84961fc4SMatthew G. Knepley ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr); 113*84961fc4SMatthew G. Knepley ierr = PetscMalloc1((fEnd - fStart), &faceFIFO);CHKERRQ(ierr); 114*84961fc4SMatthew G. Knepley fTop = fBottom = 0; 115*84961fc4SMatthew G. Knepley /* Initialize FIFO with first cell */ 116*84961fc4SMatthew G. Knepley if (cEnd > cStart) { 117*84961fc4SMatthew G. Knepley const PetscInt *cone; 118*84961fc4SMatthew G. Knepley PetscInt coneSize; 119*84961fc4SMatthew G. Knepley 120*84961fc4SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cStart, &coneSize);CHKERRQ(ierr); 121*84961fc4SMatthew G. Knepley ierr = DMPlexGetCone(dm, cStart, &cone);CHKERRQ(ierr); 122*84961fc4SMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 123*84961fc4SMatthew G. Knepley faceFIFO[fBottom++] = cone[c]; 124*84961fc4SMatthew G. Knepley ierr = PetscBTSet(seenFaces, cone[c]-fStart);CHKERRQ(ierr); 125*84961fc4SMatthew G. Knepley } 126*84961fc4SMatthew G. Knepley } 127*84961fc4SMatthew G. Knepley /* Consider each face in FIFO */ 128*84961fc4SMatthew G. Knepley while (fTop < fBottom) { 129*84961fc4SMatthew G. Knepley const PetscInt *support, *coneA, *coneB, *coneOA, *coneOB; 130*84961fc4SMatthew G. Knepley PetscInt supportSize, coneSizeA, coneSizeB, posA = -1, posB = -1; 131*84961fc4SMatthew G. Knepley PetscInt seenA, flippedA, seenB, flippedB, mismatch; 132*84961fc4SMatthew G. Knepley 133*84961fc4SMatthew G. Knepley face = faceFIFO[fTop++]; 134*84961fc4SMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); 135*84961fc4SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr); 136*84961fc4SMatthew G. Knepley if (supportSize < 2) continue; 137*84961fc4SMatthew G. Knepley if (supportSize != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %d", supportSize); 138*84961fc4SMatthew G. Knepley seenA = PetscBTLookup(seenCells, support[0]-cStart); 139*84961fc4SMatthew G. Knepley flippedA = PetscBTLookup(flippedCells, support[0]-cStart) ? 1 : 0; 140*84961fc4SMatthew G. Knepley seenB = PetscBTLookup(seenCells, support[1]-cStart); 141*84961fc4SMatthew G. Knepley flippedB = PetscBTLookup(flippedCells, support[1]-cStart) ? 1 : 0; 142*84961fc4SMatthew G. Knepley 143*84961fc4SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, support[0], &coneSizeA);CHKERRQ(ierr); 144*84961fc4SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, support[1], &coneSizeB);CHKERRQ(ierr); 145*84961fc4SMatthew G. Knepley ierr = DMPlexGetCone(dm, support[0], &coneA);CHKERRQ(ierr); 146*84961fc4SMatthew G. Knepley ierr = DMPlexGetCone(dm, support[1], &coneB);CHKERRQ(ierr); 147*84961fc4SMatthew G. Knepley ierr = DMPlexGetConeOrientation(dm, support[0], &coneOA);CHKERRQ(ierr); 148*84961fc4SMatthew G. Knepley ierr = DMPlexGetConeOrientation(dm, support[1], &coneOB);CHKERRQ(ierr); 149*84961fc4SMatthew G. Knepley for (c = 0; c < coneSizeA; ++c) { 150*84961fc4SMatthew G. Knepley if (!PetscBTLookup(seenFaces, coneA[c]-fStart)) { 151*84961fc4SMatthew G. Knepley faceFIFO[fBottom++] = coneA[c]; 152*84961fc4SMatthew G. Knepley ierr = PetscBTSet(seenFaces, coneA[c]-fStart);CHKERRQ(ierr); 153*84961fc4SMatthew G. Knepley } 154*84961fc4SMatthew G. Knepley if (coneA[c] == face) posA = c; 155*84961fc4SMatthew 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); 156*84961fc4SMatthew G. Knepley } 157*84961fc4SMatthew 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]); 158*84961fc4SMatthew G. Knepley for (c = 0; c < coneSizeB; ++c) { 159*84961fc4SMatthew G. Knepley if (!PetscBTLookup(seenFaces, coneB[c]-fStart)) { 160*84961fc4SMatthew G. Knepley faceFIFO[fBottom++] = coneB[c]; 161*84961fc4SMatthew G. Knepley ierr = PetscBTSet(seenFaces, coneB[c]-fStart);CHKERRQ(ierr); 162*84961fc4SMatthew G. Knepley } 163*84961fc4SMatthew G. Knepley if (coneB[c] == face) posB = c; 164*84961fc4SMatthew 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); 165*84961fc4SMatthew G. Knepley } 166*84961fc4SMatthew 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]); 167*84961fc4SMatthew G. Knepley 168*84961fc4SMatthew G. Knepley if (dim == 1) { 169*84961fc4SMatthew G. Knepley mismatch = posA == posB; 170*84961fc4SMatthew G. Knepley } else { 171*84961fc4SMatthew G. Knepley mismatch = coneOA[posA] == coneOB[posB]; 172*84961fc4SMatthew G. Knepley } 173*84961fc4SMatthew G. Knepley 174*84961fc4SMatthew G. Knepley if (mismatch ^ (flippedA ^ flippedB)) { 175*84961fc4SMatthew 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]); 176*84961fc4SMatthew G. Knepley if (!seenA && !flippedA) { 177*84961fc4SMatthew G. Knepley ierr = PetscBTSet(flippedCells, support[0]-cStart);CHKERRQ(ierr); 178*84961fc4SMatthew G. Knepley } else if (!seenB && !flippedB) { 179*84961fc4SMatthew G. Knepley ierr = PetscBTSet(flippedCells, support[1]-cStart);CHKERRQ(ierr); 180*84961fc4SMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable"); 181*84961fc4SMatthew 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"); 182*84961fc4SMatthew G. Knepley ierr = PetscBTSet(seenCells, support[0]-cStart);CHKERRQ(ierr); 183*84961fc4SMatthew G. Knepley ierr = PetscBTSet(seenCells, support[1]-cStart);CHKERRQ(ierr); 184*84961fc4SMatthew G. Knepley } 185*84961fc4SMatthew G. Knepley if (flg) { 186*84961fc4SMatthew G. Knepley PetscViewer v; 187*84961fc4SMatthew G. Knepley PetscMPIInt rank; 188*84961fc4SMatthew G. Knepley 189*84961fc4SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 190*84961fc4SMatthew G. Knepley ierr = PetscViewerASCIIGetStdout(PETSC_COMM_SELF, &v);CHKERRQ(ierr); 191*84961fc4SMatthew G. Knepley ierr = PetscViewerASCIISynchronizedAllow(v, PETSC_TRUE);CHKERRQ(ierr); 192*84961fc4SMatthew G. Knepley ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank);CHKERRQ(ierr); 193*84961fc4SMatthew G. Knepley ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr); 194*84961fc4SMatthew G. Knepley } 195*84961fc4SMatthew G. Knepley /* Now all subdomains are oriented, but we need a consistent parallel orientation */ 196*84961fc4SMatthew G. Knepley { 197*84961fc4SMatthew G. Knepley /* Find a representative face (edge) separating pairs of procs */ 198*84961fc4SMatthew G. Knepley PetscSF sf; 199*84961fc4SMatthew G. Knepley const PetscInt *lpoints; 200*84961fc4SMatthew G. Knepley const PetscSFNode *rpoints; 201*84961fc4SMatthew G. Knepley PetscInt *neighbors, *nranks; 202*84961fc4SMatthew G. Knepley PetscInt numLeaves, numRoots, numNeighbors = 0, l, n; 203*84961fc4SMatthew G. Knepley 204*84961fc4SMatthew G. Knepley ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 205*84961fc4SMatthew G. Knepley ierr = PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints);CHKERRQ(ierr); 206*84961fc4SMatthew G. Knepley if (numLeaves >= 0) { 207*84961fc4SMatthew G. Knepley const PetscInt *cone, *ornt, *support; 208*84961fc4SMatthew G. Knepley PetscInt coneSize, supportSize; 209*84961fc4SMatthew G. Knepley int *rornt, *lornt; /* PetscSF cannot handle smaller than int */ 210*84961fc4SMatthew G. Knepley PetscBool *match, flipped = PETSC_FALSE; 211*84961fc4SMatthew G. Knepley 212*84961fc4SMatthew G. Knepley ierr = PetscMalloc1(numLeaves,&neighbors);CHKERRQ(ierr); 213*84961fc4SMatthew G. Knepley /* I know this is p^2 time in general, but for bounded degree its alright */ 214*84961fc4SMatthew G. Knepley for (l = 0; l < numLeaves; ++l) { 215*84961fc4SMatthew G. Knepley const PetscInt face = lpoints[l]; 216*84961fc4SMatthew G. Knepley if ((face >= fStart) && (face < fEnd)) { 217*84961fc4SMatthew G. Knepley const PetscInt rank = rpoints[l].rank; 218*84961fc4SMatthew G. Knepley for (n = 0; n < numNeighbors; ++n) if (rank == rpoints[neighbors[n]].rank) break; 219*84961fc4SMatthew G. Knepley if (n >= numNeighbors) { 220*84961fc4SMatthew G. Knepley PetscInt supportSize; 221*84961fc4SMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); 222*84961fc4SMatthew G. Knepley if (supportSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %d", supportSize); 223*84961fc4SMatthew G. Knepley neighbors[numNeighbors++] = l; 224*84961fc4SMatthew G. Knepley } 225*84961fc4SMatthew G. Knepley } 226*84961fc4SMatthew G. Knepley } 227*84961fc4SMatthew G. Knepley ierr = PetscCalloc4(numNeighbors,&match,numNeighbors,&nranks,numRoots,&rornt,numRoots,&lornt);CHKERRQ(ierr); 228*84961fc4SMatthew G. Knepley for (face = fStart; face < fEnd; ++face) { 229*84961fc4SMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); 230*84961fc4SMatthew G. Knepley if (supportSize != 1) continue; 231*84961fc4SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr); 232*84961fc4SMatthew G. Knepley 233*84961fc4SMatthew G. Knepley ierr = DMPlexGetCone(dm, support[0], &cone);CHKERRQ(ierr); 234*84961fc4SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, support[0], &coneSize);CHKERRQ(ierr); 235*84961fc4SMatthew G. Knepley ierr = DMPlexGetConeOrientation(dm, support[0], &ornt);CHKERRQ(ierr); 236*84961fc4SMatthew G. Knepley for (c = 0; c < coneSize; ++c) if (cone[c] == face) break; 237*84961fc4SMatthew G. Knepley if (dim == 1) { 238*84961fc4SMatthew G. Knepley /* Use cone position instead, shifted to -1 or 1 */ 239*84961fc4SMatthew G. Knepley rornt[face] = c*2-1; 240*84961fc4SMatthew G. Knepley } else { 241*84961fc4SMatthew G. Knepley if (PetscBTLookup(flippedCells, support[0]-cStart)) rornt[face] = ornt[c] < 0 ? -1 : 1; 242*84961fc4SMatthew G. Knepley else rornt[face] = ornt[c] < 0 ? 1 : -1; 243*84961fc4SMatthew G. Knepley } 244*84961fc4SMatthew G. Knepley } 245*84961fc4SMatthew G. Knepley /* Mark each edge with match or nomatch */ 246*84961fc4SMatthew G. Knepley ierr = PetscSFBcastBegin(sf, MPI_INT, rornt, lornt);CHKERRQ(ierr); 247*84961fc4SMatthew G. Knepley ierr = PetscSFBcastEnd(sf, MPI_INT, rornt, lornt);CHKERRQ(ierr); 248*84961fc4SMatthew G. Knepley for (n = 0; n < numNeighbors; ++n) { 249*84961fc4SMatthew G. Knepley const PetscInt face = lpoints[neighbors[n]]; 250*84961fc4SMatthew G. Knepley 251*84961fc4SMatthew G. Knepley if (rornt[face]*lornt[face] < 0) match[n] = PETSC_TRUE; 252*84961fc4SMatthew G. Knepley else match[n] = PETSC_FALSE; 253*84961fc4SMatthew G. Knepley nranks[n] = rpoints[neighbors[n]].rank; 254*84961fc4SMatthew G. Knepley } 255*84961fc4SMatthew G. Knepley /* Collect the graph on 0 */ 256*84961fc4SMatthew G. Knepley { 257*84961fc4SMatthew G. Knepley Mat G; 258*84961fc4SMatthew G. Knepley PetscBT seenProcs, flippedProcs; 259*84961fc4SMatthew G. Knepley PetscInt *procFIFO, pTop, pBottom; 260*84961fc4SMatthew G. Knepley PetscInt *adj = NULL; 261*84961fc4SMatthew G. Knepley PetscBool *val = NULL; 262*84961fc4SMatthew G. Knepley PetscMPIInt *recvcounts = NULL, *displs = NULL, p; 263*84961fc4SMatthew G. Knepley PetscMPIInt N = numNeighbors, numProcs = 0, rank; 264*84961fc4SMatthew G. Knepley PetscInt debug = 0; 265*84961fc4SMatthew G. Knepley 266*84961fc4SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 267*84961fc4SMatthew G. Knepley if (!rank) {ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);} 268*84961fc4SMatthew G. Knepley ierr = PetscCalloc2(numProcs,&recvcounts,numProcs+1,&displs);CHKERRQ(ierr); 269*84961fc4SMatthew G. Knepley ierr = MPI_Gather(&N, 1, MPI_INT, recvcounts, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 270*84961fc4SMatthew G. Knepley for (p = 0; p < numProcs; ++p) { 271*84961fc4SMatthew G. Knepley displs[p+1] = displs[p] + recvcounts[p]; 272*84961fc4SMatthew G. Knepley } 273*84961fc4SMatthew G. Knepley if (!rank) {ierr = PetscMalloc2(displs[numProcs],&adj,displs[numProcs],&val);CHKERRQ(ierr);} 274*84961fc4SMatthew G. Knepley ierr = MPI_Gatherv(nranks, numNeighbors, MPIU_INT, adj, recvcounts, displs, MPIU_INT, 0, comm);CHKERRQ(ierr); 275*84961fc4SMatthew G. Knepley ierr = MPI_Gatherv(match, numNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm);CHKERRQ(ierr); 276*84961fc4SMatthew G. Knepley if (debug) { 277*84961fc4SMatthew G. Knepley for (p = 0; p < numProcs; ++p) { 278*84961fc4SMatthew G. Knepley ierr = PetscPrintf(comm, "Proc %d:\n", p); 279*84961fc4SMatthew G. Knepley for (n = 0; n < recvcounts[p]; ++n) { 280*84961fc4SMatthew G. Knepley ierr = PetscPrintf(comm, " edge %d (%d):\n", adj[displs[p]+n], val[displs[p]+n]); 281*84961fc4SMatthew G. Knepley } 282*84961fc4SMatthew G. Knepley } 283*84961fc4SMatthew G. Knepley } 284*84961fc4SMatthew G. Knepley /* Symmetrize the graph */ 285*84961fc4SMatthew G. Knepley ierr = MatCreate(PETSC_COMM_SELF, &G);CHKERRQ(ierr); 286*84961fc4SMatthew G. Knepley ierr = MatSetSizes(G, numProcs, numProcs, numProcs, numProcs);CHKERRQ(ierr); 287*84961fc4SMatthew G. Knepley ierr = MatSetUp(G);CHKERRQ(ierr); 288*84961fc4SMatthew G. Knepley for (p = 0; p < numProcs; ++p) { 289*84961fc4SMatthew G. Knepley for (n = 0; n < recvcounts[p]; ++n) { 290*84961fc4SMatthew G. Knepley const PetscInt r = p; 291*84961fc4SMatthew G. Knepley const PetscInt q = adj[displs[p]+n]; 292*84961fc4SMatthew G. Knepley const PetscScalar o = val[displs[p]+n] ? 1.0 : 0.0; 293*84961fc4SMatthew G. Knepley 294*84961fc4SMatthew G. Knepley ierr = MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES);CHKERRQ(ierr); 295*84961fc4SMatthew G. Knepley ierr = MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES);CHKERRQ(ierr); 296*84961fc4SMatthew G. Knepley } 297*84961fc4SMatthew G. Knepley } 298*84961fc4SMatthew G. Knepley ierr = MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 299*84961fc4SMatthew G. Knepley ierr = MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 300*84961fc4SMatthew G. Knepley 301*84961fc4SMatthew G. Knepley ierr = PetscBTCreate(numProcs, &seenProcs);CHKERRQ(ierr); 302*84961fc4SMatthew G. Knepley ierr = PetscBTMemzero(numProcs, seenProcs);CHKERRQ(ierr); 303*84961fc4SMatthew G. Knepley ierr = PetscBTCreate(numProcs, &flippedProcs);CHKERRQ(ierr); 304*84961fc4SMatthew G. Knepley ierr = PetscBTMemzero(numProcs, flippedProcs);CHKERRQ(ierr); 305*84961fc4SMatthew G. Knepley ierr = PetscMalloc1(numProcs,&procFIFO);CHKERRQ(ierr); 306*84961fc4SMatthew G. Knepley pTop = pBottom = 0; 307*84961fc4SMatthew G. Knepley for (p = 0; p < numProcs; ++p) { 308*84961fc4SMatthew G. Knepley if (PetscBTLookup(seenProcs, p)) continue; 309*84961fc4SMatthew G. Knepley /* Initialize FIFO with next proc */ 310*84961fc4SMatthew G. Knepley procFIFO[pBottom++] = p; 311*84961fc4SMatthew G. Knepley ierr = PetscBTSet(seenProcs, p);CHKERRQ(ierr); 312*84961fc4SMatthew G. Knepley /* Consider each proc in FIFO */ 313*84961fc4SMatthew G. Knepley while (pTop < pBottom) { 314*84961fc4SMatthew G. Knepley const PetscScalar *ornt; 315*84961fc4SMatthew G. Knepley const PetscInt *neighbors; 316*84961fc4SMatthew G. Knepley PetscInt proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors; 317*84961fc4SMatthew G. Knepley 318*84961fc4SMatthew G. Knepley proc = procFIFO[pTop++]; 319*84961fc4SMatthew G. Knepley flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0; 320*84961fc4SMatthew G. Knepley ierr = MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt);CHKERRQ(ierr); 321*84961fc4SMatthew G. Knepley /* Loop over neighboring procs */ 322*84961fc4SMatthew G. Knepley for (n = 0; n < numNeighbors; ++n) { 323*84961fc4SMatthew G. Knepley nproc = neighbors[n]; 324*84961fc4SMatthew G. Knepley mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1; 325*84961fc4SMatthew G. Knepley seen = PetscBTLookup(seenProcs, nproc); 326*84961fc4SMatthew G. Knepley flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0; 327*84961fc4SMatthew G. Knepley 328*84961fc4SMatthew G. Knepley if (mismatch ^ (flippedA ^ flippedB)) { 329*84961fc4SMatthew 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); 330*84961fc4SMatthew G. Knepley if (!flippedB) { 331*84961fc4SMatthew G. Knepley ierr = PetscBTSet(flippedProcs, nproc);CHKERRQ(ierr); 332*84961fc4SMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable"); 333*84961fc4SMatthew 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"); 334*84961fc4SMatthew G. Knepley if (!seen) { 335*84961fc4SMatthew G. Knepley procFIFO[pBottom++] = nproc; 336*84961fc4SMatthew G. Knepley ierr = PetscBTSet(seenProcs, nproc);CHKERRQ(ierr); 337*84961fc4SMatthew G. Knepley } 338*84961fc4SMatthew G. Knepley } 339*84961fc4SMatthew G. Knepley } 340*84961fc4SMatthew G. Knepley } 341*84961fc4SMatthew G. Knepley ierr = PetscFree(procFIFO);CHKERRQ(ierr); 342*84961fc4SMatthew G. Knepley ierr = MatDestroy(&G);CHKERRQ(ierr); 343*84961fc4SMatthew G. Knepley 344*84961fc4SMatthew G. Knepley ierr = PetscFree2(recvcounts,displs);CHKERRQ(ierr); 345*84961fc4SMatthew G. Knepley ierr = PetscFree2(adj,val);CHKERRQ(ierr); 346*84961fc4SMatthew G. Knepley { 347*84961fc4SMatthew G. Knepley PetscBool *flips; 348*84961fc4SMatthew G. Knepley 349*84961fc4SMatthew G. Knepley ierr = PetscMalloc1(numProcs,&flips);CHKERRQ(ierr); 350*84961fc4SMatthew G. Knepley for (p = 0; p < numProcs; ++p) { 351*84961fc4SMatthew G. Knepley flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE; 352*84961fc4SMatthew G. Knepley if (debug && flips[p]) {ierr = PetscPrintf(comm, "Flipping Proc %d:\n", p);} 353*84961fc4SMatthew G. Knepley } 354*84961fc4SMatthew G. Knepley ierr = MPI_Scatter(flips, 1, MPIU_BOOL, &flipped, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr); 355*84961fc4SMatthew G. Knepley ierr = PetscFree(flips);CHKERRQ(ierr); 356*84961fc4SMatthew G. Knepley } 357*84961fc4SMatthew G. Knepley ierr = PetscBTDestroy(&seenProcs);CHKERRQ(ierr); 358*84961fc4SMatthew G. Knepley ierr = PetscBTDestroy(&flippedProcs);CHKERRQ(ierr); 359*84961fc4SMatthew G. Knepley } 360*84961fc4SMatthew G. Knepley ierr = PetscFree4(match,nranks,rornt,lornt);CHKERRQ(ierr); 361*84961fc4SMatthew G. Knepley ierr = PetscFree(neighbors);CHKERRQ(ierr); 362*84961fc4SMatthew G. Knepley if (flipped) {for (c = cStart; c < cEnd; ++c) {ierr = PetscBTNegate(flippedCells, c-cStart);CHKERRQ(ierr);}} 363*84961fc4SMatthew G. Knepley } 364*84961fc4SMatthew G. Knepley } 365*84961fc4SMatthew G. Knepley if (flg) { 366*84961fc4SMatthew G. Knepley PetscViewer v; 367*84961fc4SMatthew G. Knepley PetscMPIInt rank; 368*84961fc4SMatthew G. Knepley 369*84961fc4SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 370*84961fc4SMatthew G. Knepley ierr = PetscViewerASCIIGetStdout(PETSC_COMM_SELF, &v);CHKERRQ(ierr); 371*84961fc4SMatthew G. Knepley ierr = PetscViewerASCIISynchronizedAllow(v, PETSC_TRUE);CHKERRQ(ierr); 372*84961fc4SMatthew G. Knepley ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank);CHKERRQ(ierr); 373*84961fc4SMatthew G. Knepley ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr); 374*84961fc4SMatthew G. Knepley } 375*84961fc4SMatthew G. Knepley /* Reverse flipped cells in the mesh */ 376*84961fc4SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 377*84961fc4SMatthew G. Knepley if (PetscBTLookup(flippedCells, c-cStart)) {ierr = DMPlexReverseCell(dm, c);CHKERRQ(ierr);} 378*84961fc4SMatthew G. Knepley } 379*84961fc4SMatthew G. Knepley ierr = PetscBTDestroy(&seenCells);CHKERRQ(ierr); 380*84961fc4SMatthew G. Knepley ierr = PetscBTDestroy(&flippedCells);CHKERRQ(ierr); 381*84961fc4SMatthew G. Knepley ierr = PetscBTDestroy(&seenFaces);CHKERRQ(ierr); 382*84961fc4SMatthew G. Knepley ierr = PetscFree(faceFIFO);CHKERRQ(ierr); 383*84961fc4SMatthew G. Knepley PetscFunctionReturn(0); 384*84961fc4SMatthew G. Knepley } 385