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__ 68*7b310ce2SMatthew G. Knepley #define __FUNCT__ "DMPlexCheckFace_Internal" 69*7b310ce2SMatthew G. Knepley /* 70*7b310ce2SMatthew G. Knepley - Checks face match 71*7b310ce2SMatthew G. Knepley - Flips non-matching 72*7b310ce2SMatthew G. Knepley - Inserts faces of support cells in FIFO 73*7b310ce2SMatthew G. Knepley */ 74*7b310ce2SMatthew G. Knepley static PetscErrorCode DMPlexCheckFace_Internal(DM dm, PetscInt face, PetscInt *faceFIFO, PetscInt *fTop, PetscInt *fBottom, PetscInt cStart, PetscInt fStart, PetscInt fEnd, PetscBT seenCells, PetscBT flippedCells, PetscBT seenFaces) 75*7b310ce2SMatthew G. Knepley { 76*7b310ce2SMatthew G. Knepley const PetscInt *support, *coneA, *coneB, *coneOA, *coneOB; 77*7b310ce2SMatthew G. Knepley PetscInt supportSize, coneSizeA, coneSizeB, posA = -1, posB = -1; 78*7b310ce2SMatthew G. Knepley PetscInt dim, seenA, flippedA, seenB, flippedB, mismatch, c; 79*7b310ce2SMatthew G. Knepley PetscErrorCode ierr; 80*7b310ce2SMatthew G. Knepley 81*7b310ce2SMatthew G. Knepley PetscFunctionBegin; 82*7b310ce2SMatthew G. Knepley face = faceFIFO[(*fTop)++]; 83*7b310ce2SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 84*7b310ce2SMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); 85*7b310ce2SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr); 86*7b310ce2SMatthew G. Knepley if (supportSize < 2) PetscFunctionReturn(0); 87*7b310ce2SMatthew G. Knepley if (supportSize != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %d", supportSize); 88*7b310ce2SMatthew G. Knepley seenA = PetscBTLookup(seenCells, support[0]-cStart); 89*7b310ce2SMatthew G. Knepley flippedA = PetscBTLookup(flippedCells, support[0]-cStart) ? 1 : 0; 90*7b310ce2SMatthew G. Knepley seenB = PetscBTLookup(seenCells, support[1]-cStart); 91*7b310ce2SMatthew G. Knepley flippedB = PetscBTLookup(flippedCells, support[1]-cStart) ? 1 : 0; 92*7b310ce2SMatthew G. Knepley 93*7b310ce2SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, support[0], &coneSizeA);CHKERRQ(ierr); 94*7b310ce2SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, support[1], &coneSizeB);CHKERRQ(ierr); 95*7b310ce2SMatthew G. Knepley ierr = DMPlexGetCone(dm, support[0], &coneA);CHKERRQ(ierr); 96*7b310ce2SMatthew G. Knepley ierr = DMPlexGetCone(dm, support[1], &coneB);CHKERRQ(ierr); 97*7b310ce2SMatthew G. Knepley ierr = DMPlexGetConeOrientation(dm, support[0], &coneOA);CHKERRQ(ierr); 98*7b310ce2SMatthew G. Knepley ierr = DMPlexGetConeOrientation(dm, support[1], &coneOB);CHKERRQ(ierr); 99*7b310ce2SMatthew G. Knepley for (c = 0; c < coneSizeA; ++c) { 100*7b310ce2SMatthew G. Knepley if (!PetscBTLookup(seenFaces, coneA[c]-fStart)) { 101*7b310ce2SMatthew G. Knepley faceFIFO[(*fBottom)++] = coneA[c]; 102*7b310ce2SMatthew G. Knepley ierr = PetscBTSet(seenFaces, coneA[c]-fStart);CHKERRQ(ierr); 103*7b310ce2SMatthew G. Knepley } 104*7b310ce2SMatthew G. Knepley if (coneA[c] == face) posA = c; 105*7b310ce2SMatthew 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); 106*7b310ce2SMatthew G. Knepley } 107*7b310ce2SMatthew 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]); 108*7b310ce2SMatthew G. Knepley for (c = 0; c < coneSizeB; ++c) { 109*7b310ce2SMatthew G. Knepley if (!PetscBTLookup(seenFaces, coneB[c]-fStart)) { 110*7b310ce2SMatthew G. Knepley faceFIFO[(*fBottom)++] = coneB[c]; 111*7b310ce2SMatthew G. Knepley ierr = PetscBTSet(seenFaces, coneB[c]-fStart);CHKERRQ(ierr); 112*7b310ce2SMatthew G. Knepley } 113*7b310ce2SMatthew G. Knepley if (coneB[c] == face) posB = c; 114*7b310ce2SMatthew 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); 115*7b310ce2SMatthew G. Knepley } 116*7b310ce2SMatthew 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]); 117*7b310ce2SMatthew G. Knepley 118*7b310ce2SMatthew G. Knepley if (dim == 1) { 119*7b310ce2SMatthew G. Knepley mismatch = posA == posB; 120*7b310ce2SMatthew G. Knepley } else { 121*7b310ce2SMatthew G. Knepley mismatch = coneOA[posA] == coneOB[posB]; 122*7b310ce2SMatthew G. Knepley } 123*7b310ce2SMatthew G. Knepley 124*7b310ce2SMatthew G. Knepley if (mismatch ^ (flippedA ^ flippedB)) { 125*7b310ce2SMatthew 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]); 126*7b310ce2SMatthew G. Knepley if (!seenA && !flippedA) { 127*7b310ce2SMatthew G. Knepley ierr = PetscBTSet(flippedCells, support[0]-cStart);CHKERRQ(ierr); 128*7b310ce2SMatthew G. Knepley } else if (!seenB && !flippedB) { 129*7b310ce2SMatthew G. Knepley ierr = PetscBTSet(flippedCells, support[1]-cStart);CHKERRQ(ierr); 130*7b310ce2SMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable"); 131*7b310ce2SMatthew 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"); 132*7b310ce2SMatthew G. Knepley ierr = PetscBTSet(seenCells, support[0]-cStart);CHKERRQ(ierr); 133*7b310ce2SMatthew G. Knepley ierr = PetscBTSet(seenCells, support[1]-cStart);CHKERRQ(ierr); 134*7b310ce2SMatthew G. Knepley PetscFunctionReturn(0); 135*7b310ce2SMatthew G. Knepley } 136*7b310ce2SMatthew G. Knepley 137*7b310ce2SMatthew G. Knepley #undef __FUNCT__ 13884961fc4SMatthew G. Knepley #define __FUNCT__ "DMPlexOrient" 13984961fc4SMatthew G. Knepley /*@ 14084961fc4SMatthew G. Knepley DMPlexOrient - Give a consistent orientation to the input mesh 14184961fc4SMatthew G. Knepley 14284961fc4SMatthew G. Knepley Input Parameters: 14384961fc4SMatthew G. Knepley . dm - The DM 14484961fc4SMatthew G. Knepley 14584961fc4SMatthew G. Knepley Note: The orientation data for the DM are change in-place. 14684961fc4SMatthew G. Knepley $ This routine will fail for non-orientable surfaces, such as the Moebius strip. 14784961fc4SMatthew G. Knepley 14884961fc4SMatthew G. Knepley Level: advanced 14984961fc4SMatthew G. Knepley 15084961fc4SMatthew G. Knepley .seealso: DMCreate(), DMPLEX 15184961fc4SMatthew G. Knepley @*/ 15284961fc4SMatthew G. Knepley PetscErrorCode DMPlexOrient(DM dm) 15384961fc4SMatthew G. Knepley { 15484961fc4SMatthew G. Knepley MPI_Comm comm; 15584961fc4SMatthew G. Knepley PetscBT seenCells, flippedCells, seenFaces; 15684961fc4SMatthew G. Knepley PetscInt *faceFIFO, fTop, fBottom; 15784961fc4SMatthew G. Knepley PetscInt dim, h, cStart, cEnd, c, fStart, fEnd, face; 15884961fc4SMatthew G. Knepley PetscBool flg; 15984961fc4SMatthew G. Knepley PetscErrorCode ierr; 16084961fc4SMatthew G. Knepley 16184961fc4SMatthew G. Knepley PetscFunctionBegin; 16284961fc4SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 16384961fc4SMatthew G. Knepley ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-orientation_view", &flg);CHKERRQ(ierr); 16484961fc4SMatthew G. Knepley /* Truth Table 16584961fc4SMatthew G. Knepley mismatch flips do action mismatch flipA ^ flipB action 16684961fc4SMatthew G. Knepley F 0 flips no F F F 16784961fc4SMatthew G. Knepley F 1 flip yes F T T 16884961fc4SMatthew G. Knepley F 2 flips no T F T 16984961fc4SMatthew G. Knepley T 0 flips yes T T F 17084961fc4SMatthew G. Knepley T 1 flip no 17184961fc4SMatthew G. Knepley T 2 flips yes 17284961fc4SMatthew G. Knepley */ 17384961fc4SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 17484961fc4SMatthew G. Knepley ierr = DMPlexGetVTKCellHeight(dm, &h);CHKERRQ(ierr); 17584961fc4SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, h, &cStart, &cEnd);CHKERRQ(ierr); 17684961fc4SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, h+1, &fStart, &fEnd);CHKERRQ(ierr); 17784961fc4SMatthew G. Knepley ierr = PetscBTCreate(cEnd - cStart, &seenCells);CHKERRQ(ierr); 17884961fc4SMatthew G. Knepley ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr); 17984961fc4SMatthew G. Knepley ierr = PetscBTCreate(cEnd - cStart, &flippedCells);CHKERRQ(ierr); 18084961fc4SMatthew G. Knepley ierr = PetscBTMemzero(cEnd - cStart, flippedCells);CHKERRQ(ierr); 18184961fc4SMatthew G. Knepley ierr = PetscBTCreate(fEnd - fStart, &seenFaces);CHKERRQ(ierr); 18284961fc4SMatthew G. Knepley ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr); 18384961fc4SMatthew G. Knepley ierr = PetscMalloc1((fEnd - fStart), &faceFIFO);CHKERRQ(ierr); 18484961fc4SMatthew G. Knepley fTop = fBottom = 0; 18584961fc4SMatthew G. Knepley /* Initialize FIFO with first cell */ 18684961fc4SMatthew G. Knepley if (cEnd > cStart) { 18784961fc4SMatthew G. Knepley const PetscInt *cone; 18884961fc4SMatthew G. Knepley PetscInt coneSize; 18984961fc4SMatthew G. Knepley 19084961fc4SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cStart, &coneSize);CHKERRQ(ierr); 19184961fc4SMatthew G. Knepley ierr = DMPlexGetCone(dm, cStart, &cone);CHKERRQ(ierr); 19284961fc4SMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 19384961fc4SMatthew G. Knepley faceFIFO[fBottom++] = cone[c]; 19484961fc4SMatthew G. Knepley ierr = PetscBTSet(seenFaces, cone[c]-fStart);CHKERRQ(ierr); 19584961fc4SMatthew G. Knepley } 19684961fc4SMatthew G. Knepley } 19784961fc4SMatthew G. Knepley /* Consider each face in FIFO */ 19884961fc4SMatthew G. Knepley while (fTop < fBottom) { 199*7b310ce2SMatthew G. Knepley ierr = DMPlexCheckFace_Internal(dm, face, faceFIFO, &fTop, &fBottom, cStart, fStart, fEnd, seenCells, flippedCells, seenFaces);CHKERRQ(ierr); 20084961fc4SMatthew G. Knepley } 20184961fc4SMatthew G. Knepley if (flg) { 20284961fc4SMatthew G. Knepley PetscViewer v; 20384961fc4SMatthew G. Knepley PetscMPIInt rank; 20484961fc4SMatthew G. Knepley 20584961fc4SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 2067e09831bSMatthew G. Knepley ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr); 2077e09831bSMatthew G. Knepley ierr = PetscViewerASCIISynchronizedAllow(v, PETSC_TRUE);CHKERRQ(ierr); 2087e09831bSMatthew G. Knepley ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial seen faces:\n", rank);CHKERRQ(ierr); 2097e09831bSMatthew G. Knepley ierr = PetscBTView(fEnd-fStart, seenFaces, v);CHKERRQ(ierr); 21084961fc4SMatthew G. Knepley ierr = PetscViewerASCIISynchronizedAllow(v, PETSC_TRUE);CHKERRQ(ierr); 21184961fc4SMatthew G. Knepley ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank);CHKERRQ(ierr); 21284961fc4SMatthew G. Knepley ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr); 21384961fc4SMatthew G. Knepley } 21484961fc4SMatthew G. Knepley /* Now all subdomains are oriented, but we need a consistent parallel orientation */ 21584961fc4SMatthew G. Knepley { 21684961fc4SMatthew G. Knepley /* Find a representative face (edge) separating pairs of procs */ 21784961fc4SMatthew G. Knepley PetscSF sf; 21884961fc4SMatthew G. Knepley const PetscInt *lpoints; 21984961fc4SMatthew G. Knepley const PetscSFNode *rpoints; 22084961fc4SMatthew G. Knepley PetscInt *neighbors, *nranks; 22184961fc4SMatthew G. Knepley PetscInt numLeaves, numRoots, numNeighbors = 0, l, n; 22284961fc4SMatthew G. Knepley 22384961fc4SMatthew G. Knepley ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 22484961fc4SMatthew G. Knepley ierr = PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints);CHKERRQ(ierr); 22584961fc4SMatthew G. Knepley if (numLeaves >= 0) { 22684961fc4SMatthew G. Knepley const PetscInt *cone, *ornt, *support; 22784961fc4SMatthew G. Knepley PetscInt coneSize, supportSize; 22884961fc4SMatthew G. Knepley int *rornt, *lornt; /* PetscSF cannot handle smaller than int */ 22984961fc4SMatthew G. Knepley PetscBool *match, flipped = PETSC_FALSE; 23084961fc4SMatthew G. Knepley 23184961fc4SMatthew G. Knepley ierr = PetscMalloc1(numLeaves,&neighbors);CHKERRQ(ierr); 23284961fc4SMatthew G. Knepley /* I know this is p^2 time in general, but for bounded degree its alright */ 23384961fc4SMatthew G. Knepley for (l = 0; l < numLeaves; ++l) { 23484961fc4SMatthew G. Knepley const PetscInt face = lpoints[l]; 23584961fc4SMatthew G. Knepley if ((face >= fStart) && (face < fEnd)) { 23684961fc4SMatthew G. Knepley const PetscInt rank = rpoints[l].rank; 23784961fc4SMatthew G. Knepley for (n = 0; n < numNeighbors; ++n) if (rank == rpoints[neighbors[n]].rank) break; 23884961fc4SMatthew G. Knepley if (n >= numNeighbors) { 23984961fc4SMatthew G. Knepley PetscInt supportSize; 24084961fc4SMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); 24184961fc4SMatthew G. Knepley if (supportSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %d", supportSize); 24284961fc4SMatthew G. Knepley neighbors[numNeighbors++] = l; 24384961fc4SMatthew G. Knepley } 24484961fc4SMatthew G. Knepley } 24584961fc4SMatthew G. Knepley } 24684961fc4SMatthew G. Knepley ierr = PetscCalloc4(numNeighbors,&match,numNeighbors,&nranks,numRoots,&rornt,numRoots,&lornt);CHKERRQ(ierr); 24784961fc4SMatthew G. Knepley for (face = fStart; face < fEnd; ++face) { 24884961fc4SMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); 24984961fc4SMatthew G. Knepley if (supportSize != 1) continue; 25084961fc4SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr); 25184961fc4SMatthew G. Knepley 25284961fc4SMatthew G. Knepley ierr = DMPlexGetCone(dm, support[0], &cone);CHKERRQ(ierr); 25384961fc4SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, support[0], &coneSize);CHKERRQ(ierr); 25484961fc4SMatthew G. Knepley ierr = DMPlexGetConeOrientation(dm, support[0], &ornt);CHKERRQ(ierr); 25584961fc4SMatthew G. Knepley for (c = 0; c < coneSize; ++c) if (cone[c] == face) break; 25684961fc4SMatthew G. Knepley if (dim == 1) { 25784961fc4SMatthew G. Knepley /* Use cone position instead, shifted to -1 or 1 */ 25884961fc4SMatthew G. Knepley rornt[face] = c*2-1; 25984961fc4SMatthew G. Knepley } else { 26084961fc4SMatthew G. Knepley if (PetscBTLookup(flippedCells, support[0]-cStart)) rornt[face] = ornt[c] < 0 ? -1 : 1; 26184961fc4SMatthew G. Knepley else rornt[face] = ornt[c] < 0 ? 1 : -1; 26284961fc4SMatthew G. Knepley } 26384961fc4SMatthew G. Knepley } 26484961fc4SMatthew G. Knepley /* Mark each edge with match or nomatch */ 26584961fc4SMatthew G. Knepley ierr = PetscSFBcastBegin(sf, MPI_INT, rornt, lornt);CHKERRQ(ierr); 26684961fc4SMatthew G. Knepley ierr = PetscSFBcastEnd(sf, MPI_INT, rornt, lornt);CHKERRQ(ierr); 26784961fc4SMatthew G. Knepley for (n = 0; n < numNeighbors; ++n) { 26884961fc4SMatthew G. Knepley const PetscInt face = lpoints[neighbors[n]]; 26984961fc4SMatthew G. Knepley 27084961fc4SMatthew G. Knepley if (rornt[face]*lornt[face] < 0) match[n] = PETSC_TRUE; 27184961fc4SMatthew G. Knepley else match[n] = PETSC_FALSE; 27284961fc4SMatthew G. Knepley nranks[n] = rpoints[neighbors[n]].rank; 27384961fc4SMatthew G. Knepley } 27484961fc4SMatthew G. Knepley /* Collect the graph on 0 */ 27584961fc4SMatthew G. Knepley { 27684961fc4SMatthew G. Knepley Mat G; 27784961fc4SMatthew G. Knepley PetscBT seenProcs, flippedProcs; 27884961fc4SMatthew G. Knepley PetscInt *procFIFO, pTop, pBottom; 27984961fc4SMatthew G. Knepley PetscInt *adj = NULL; 28084961fc4SMatthew G. Knepley PetscBool *val = NULL; 28184961fc4SMatthew G. Knepley PetscMPIInt *recvcounts = NULL, *displs = NULL, p; 28284961fc4SMatthew G. Knepley PetscMPIInt N = numNeighbors, numProcs = 0, rank; 28384961fc4SMatthew G. Knepley PetscInt debug = 0; 28484961fc4SMatthew G. Knepley 28584961fc4SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 28684961fc4SMatthew G. Knepley if (!rank) {ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);} 28784961fc4SMatthew G. Knepley ierr = PetscCalloc2(numProcs,&recvcounts,numProcs+1,&displs);CHKERRQ(ierr); 28884961fc4SMatthew G. Knepley ierr = MPI_Gather(&N, 1, MPI_INT, recvcounts, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 28984961fc4SMatthew G. Knepley for (p = 0; p < numProcs; ++p) { 29084961fc4SMatthew G. Knepley displs[p+1] = displs[p] + recvcounts[p]; 29184961fc4SMatthew G. Knepley } 29284961fc4SMatthew G. Knepley if (!rank) {ierr = PetscMalloc2(displs[numProcs],&adj,displs[numProcs],&val);CHKERRQ(ierr);} 29384961fc4SMatthew G. Knepley ierr = MPI_Gatherv(nranks, numNeighbors, MPIU_INT, adj, recvcounts, displs, MPIU_INT, 0, comm);CHKERRQ(ierr); 29484961fc4SMatthew G. Knepley ierr = MPI_Gatherv(match, numNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm);CHKERRQ(ierr); 29584961fc4SMatthew G. Knepley if (debug) { 29684961fc4SMatthew G. Knepley for (p = 0; p < numProcs; ++p) { 29784961fc4SMatthew G. Knepley ierr = PetscPrintf(comm, "Proc %d:\n", p); 29884961fc4SMatthew G. Knepley for (n = 0; n < recvcounts[p]; ++n) { 29984961fc4SMatthew G. Knepley ierr = PetscPrintf(comm, " edge %d (%d):\n", adj[displs[p]+n], val[displs[p]+n]); 30084961fc4SMatthew G. Knepley } 30184961fc4SMatthew G. Knepley } 30284961fc4SMatthew G. Knepley } 30384961fc4SMatthew G. Knepley /* Symmetrize the graph */ 30484961fc4SMatthew G. Knepley ierr = MatCreate(PETSC_COMM_SELF, &G);CHKERRQ(ierr); 30584961fc4SMatthew G. Knepley ierr = MatSetSizes(G, numProcs, numProcs, numProcs, numProcs);CHKERRQ(ierr); 30684961fc4SMatthew G. Knepley ierr = MatSetUp(G);CHKERRQ(ierr); 30784961fc4SMatthew G. Knepley for (p = 0; p < numProcs; ++p) { 30884961fc4SMatthew G. Knepley for (n = 0; n < recvcounts[p]; ++n) { 30984961fc4SMatthew G. Knepley const PetscInt r = p; 31084961fc4SMatthew G. Knepley const PetscInt q = adj[displs[p]+n]; 31184961fc4SMatthew G. Knepley const PetscScalar o = val[displs[p]+n] ? 1.0 : 0.0; 31284961fc4SMatthew G. Knepley 31384961fc4SMatthew G. Knepley ierr = MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES);CHKERRQ(ierr); 31484961fc4SMatthew G. Knepley ierr = MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES);CHKERRQ(ierr); 31584961fc4SMatthew G. Knepley } 31684961fc4SMatthew G. Knepley } 31784961fc4SMatthew G. Knepley ierr = MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 31884961fc4SMatthew G. Knepley ierr = MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 31984961fc4SMatthew G. Knepley 32084961fc4SMatthew G. Knepley ierr = PetscBTCreate(numProcs, &seenProcs);CHKERRQ(ierr); 32184961fc4SMatthew G. Knepley ierr = PetscBTMemzero(numProcs, seenProcs);CHKERRQ(ierr); 32284961fc4SMatthew G. Knepley ierr = PetscBTCreate(numProcs, &flippedProcs);CHKERRQ(ierr); 32384961fc4SMatthew G. Knepley ierr = PetscBTMemzero(numProcs, flippedProcs);CHKERRQ(ierr); 32484961fc4SMatthew G. Knepley ierr = PetscMalloc1(numProcs,&procFIFO);CHKERRQ(ierr); 32584961fc4SMatthew G. Knepley pTop = pBottom = 0; 32684961fc4SMatthew G. Knepley for (p = 0; p < numProcs; ++p) { 32784961fc4SMatthew G. Knepley if (PetscBTLookup(seenProcs, p)) continue; 32884961fc4SMatthew G. Knepley /* Initialize FIFO with next proc */ 32984961fc4SMatthew G. Knepley procFIFO[pBottom++] = p; 33084961fc4SMatthew G. Knepley ierr = PetscBTSet(seenProcs, p);CHKERRQ(ierr); 33184961fc4SMatthew G. Knepley /* Consider each proc in FIFO */ 33284961fc4SMatthew G. Knepley while (pTop < pBottom) { 33384961fc4SMatthew G. Knepley const PetscScalar *ornt; 33484961fc4SMatthew G. Knepley const PetscInt *neighbors; 33584961fc4SMatthew G. Knepley PetscInt proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors; 33684961fc4SMatthew G. Knepley 33784961fc4SMatthew G. Knepley proc = procFIFO[pTop++]; 33884961fc4SMatthew G. Knepley flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0; 33984961fc4SMatthew G. Knepley ierr = MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt);CHKERRQ(ierr); 34084961fc4SMatthew G. Knepley /* Loop over neighboring procs */ 34184961fc4SMatthew G. Knepley for (n = 0; n < numNeighbors; ++n) { 34284961fc4SMatthew G. Knepley nproc = neighbors[n]; 34384961fc4SMatthew G. Knepley mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1; 34484961fc4SMatthew G. Knepley seen = PetscBTLookup(seenProcs, nproc); 34584961fc4SMatthew G. Knepley flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0; 34684961fc4SMatthew G. Knepley 34784961fc4SMatthew G. Knepley if (mismatch ^ (flippedA ^ flippedB)) { 34884961fc4SMatthew 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); 34984961fc4SMatthew G. Knepley if (!flippedB) { 35084961fc4SMatthew G. Knepley ierr = PetscBTSet(flippedProcs, nproc);CHKERRQ(ierr); 35184961fc4SMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable"); 35284961fc4SMatthew 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"); 35384961fc4SMatthew G. Knepley if (!seen) { 35484961fc4SMatthew G. Knepley procFIFO[pBottom++] = nproc; 35584961fc4SMatthew G. Knepley ierr = PetscBTSet(seenProcs, nproc);CHKERRQ(ierr); 35684961fc4SMatthew G. Knepley } 35784961fc4SMatthew G. Knepley } 35884961fc4SMatthew G. Knepley } 35984961fc4SMatthew G. Knepley } 36084961fc4SMatthew G. Knepley ierr = PetscFree(procFIFO);CHKERRQ(ierr); 36184961fc4SMatthew G. Knepley ierr = MatDestroy(&G);CHKERRQ(ierr); 36284961fc4SMatthew G. Knepley 36384961fc4SMatthew G. Knepley ierr = PetscFree2(recvcounts,displs);CHKERRQ(ierr); 36484961fc4SMatthew G. Knepley ierr = PetscFree2(adj,val);CHKERRQ(ierr); 36584961fc4SMatthew G. Knepley { 36684961fc4SMatthew G. Knepley PetscBool *flips; 36784961fc4SMatthew G. Knepley 36884961fc4SMatthew G. Knepley ierr = PetscMalloc1(numProcs,&flips);CHKERRQ(ierr); 36984961fc4SMatthew G. Knepley for (p = 0; p < numProcs; ++p) { 37084961fc4SMatthew G. Knepley flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE; 37184961fc4SMatthew G. Knepley if (debug && flips[p]) {ierr = PetscPrintf(comm, "Flipping Proc %d:\n", p);} 37284961fc4SMatthew G. Knepley } 37384961fc4SMatthew G. Knepley ierr = MPI_Scatter(flips, 1, MPIU_BOOL, &flipped, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr); 37484961fc4SMatthew G. Knepley ierr = PetscFree(flips);CHKERRQ(ierr); 37584961fc4SMatthew G. Knepley } 37684961fc4SMatthew G. Knepley ierr = PetscBTDestroy(&seenProcs);CHKERRQ(ierr); 37784961fc4SMatthew G. Knepley ierr = PetscBTDestroy(&flippedProcs);CHKERRQ(ierr); 37884961fc4SMatthew G. Knepley } 37984961fc4SMatthew G. Knepley ierr = PetscFree4(match,nranks,rornt,lornt);CHKERRQ(ierr); 38084961fc4SMatthew G. Knepley ierr = PetscFree(neighbors);CHKERRQ(ierr); 38184961fc4SMatthew G. Knepley if (flipped) {for (c = cStart; c < cEnd; ++c) {ierr = PetscBTNegate(flippedCells, c-cStart);CHKERRQ(ierr);}} 38284961fc4SMatthew G. Knepley } 38384961fc4SMatthew G. Knepley } 38484961fc4SMatthew G. Knepley if (flg) { 38584961fc4SMatthew G. Knepley PetscViewer v; 38684961fc4SMatthew G. Knepley PetscMPIInt rank; 38784961fc4SMatthew G. Knepley 38884961fc4SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 3897e09831bSMatthew G. Knepley ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr); 39084961fc4SMatthew G. Knepley ierr = PetscViewerASCIISynchronizedAllow(v, PETSC_TRUE);CHKERRQ(ierr); 39184961fc4SMatthew G. Knepley ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank);CHKERRQ(ierr); 39284961fc4SMatthew G. Knepley ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr); 39384961fc4SMatthew G. Knepley } 39484961fc4SMatthew G. Knepley /* Reverse flipped cells in the mesh */ 39584961fc4SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 39684961fc4SMatthew G. Knepley if (PetscBTLookup(flippedCells, c-cStart)) {ierr = DMPlexReverseCell(dm, c);CHKERRQ(ierr);} 39784961fc4SMatthew G. Knepley } 39884961fc4SMatthew G. Knepley ierr = PetscBTDestroy(&seenCells);CHKERRQ(ierr); 39984961fc4SMatthew G. Knepley ierr = PetscBTDestroy(&flippedCells);CHKERRQ(ierr); 40084961fc4SMatthew G. Knepley ierr = PetscBTDestroy(&seenFaces);CHKERRQ(ierr); 40184961fc4SMatthew G. Knepley ierr = PetscFree(faceFIFO);CHKERRQ(ierr); 40284961fc4SMatthew G. Knepley PetscFunctionReturn(0); 40384961fc4SMatthew G. Knepley } 404