1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 #include <petscsf.h> 3 4 /*@ 5 DMPlexOrientPoint - Act with the given orientation on the cone points of this mesh point, and update its use in the mesh. 6 7 Not Collective 8 9 Input Parameters: 10 + dm - The `DM` 11 . p - The mesh point 12 - o - The orientation 13 14 Level: intermediate 15 16 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexOrient()`, `DMPlexGetCone()`, `DMPlexGetConeOrientation()`, `DMPlexInterpolate()`, `DMPlexGetChart()` 17 @*/ 18 PetscErrorCode DMPlexOrientPoint(DM dm, PetscInt p, PetscInt o) 19 { 20 DMPolytopeType ct; 21 const PetscInt *arr, *cone, *ornt, *support; 22 PetscInt *newcone, *newornt; 23 PetscInt coneSize, c, supportSize, s; 24 25 PetscFunctionBegin; 26 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 27 PetscCall(DMPlexGetCellType(dm, p, &ct)); 28 arr = DMPolytopeTypeGetArrangement(ct, o); 29 if (!arr) PetscFunctionReturn(PETSC_SUCCESS); 30 PetscCall(DMPlexGetConeSize(dm, p, &coneSize)); 31 PetscCall(DMPlexGetCone(dm, p, &cone)); 32 PetscCall(DMPlexGetConeOrientation(dm, p, &ornt)); 33 PetscCall(DMGetWorkArray(dm, coneSize, MPIU_INT, &newcone)); 34 PetscCall(DMGetWorkArray(dm, coneSize, MPIU_INT, &newornt)); 35 for (c = 0; c < coneSize; ++c) { 36 DMPolytopeType ft; 37 PetscInt nO; 38 39 PetscCall(DMPlexGetCellType(dm, cone[c], &ft)); 40 nO = DMPolytopeTypeGetNumArrangements(ft) / 2; 41 newcone[c] = cone[arr[c * 2 + 0]]; 42 newornt[c] = DMPolytopeTypeComposeOrientation(ft, arr[c * 2 + 1], ornt[arr[c * 2 + 0]]); 43 PetscCheck(!newornt[c] || !(newornt[c] >= nO || newornt[c] < -nO), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid orientation %" PetscInt_FMT " not in [%" PetscInt_FMT ",%" PetscInt_FMT ") for %s %" PetscInt_FMT, newornt[c], -nO, nO, DMPolytopeTypes[ft], cone[c]); 44 } 45 PetscCall(DMPlexSetCone(dm, p, newcone)); 46 PetscCall(DMPlexSetConeOrientation(dm, p, newornt)); 47 PetscCall(DMRestoreWorkArray(dm, coneSize, MPIU_INT, &newcone)); 48 PetscCall(DMRestoreWorkArray(dm, coneSize, MPIU_INT, &newornt)); 49 /* Update orientation of this point in the support points */ 50 PetscCall(DMPlexGetSupportSize(dm, p, &supportSize)); 51 PetscCall(DMPlexGetSupport(dm, p, &support)); 52 for (s = 0; s < supportSize; ++s) { 53 PetscCall(DMPlexGetConeSize(dm, support[s], &coneSize)); 54 PetscCall(DMPlexGetCone(dm, support[s], &cone)); 55 PetscCall(DMPlexGetConeOrientation(dm, support[s], &ornt)); 56 for (c = 0; c < coneSize; ++c) { 57 PetscInt po; 58 59 if (cone[c] != p) continue; 60 /* ornt[c] * 0 = target = po * o so that po = ornt[c] * o^{-1} */ 61 po = DMPolytopeTypeComposeOrientationInv(ct, ornt[c], o); 62 PetscCall(DMPlexInsertConeOrientation(dm, support[s], c, po)); 63 } 64 } 65 PetscFunctionReturn(PETSC_SUCCESS); 66 } 67 68 static PetscInt GetPointIndex(PetscInt point, PetscInt pStart, PetscInt pEnd, const PetscInt points[]) 69 { 70 if (points) { 71 PetscInt loc; 72 73 PetscCall(PetscFindInt(point, pEnd - pStart, points, &loc)); 74 if (loc >= 0) return loc; 75 } else { 76 if (point >= pStart && point < pEnd) return point - pStart; 77 } 78 return -1; 79 } 80 81 /* 82 - Checks face match 83 - Flips non-matching 84 - Inserts faces of support cells in FIFO 85 */ 86 static PetscErrorCode DMPlexCheckFace_Internal(DM dm, PetscInt *faceFIFO, PetscInt *fTop, PetscInt *fBottom, IS cellIS, IS faceIS, PetscBT seenCells, PetscBT flippedCells, PetscBT seenFaces) 87 { 88 const PetscInt *supp, *coneA, *coneB, *coneOA, *coneOB; 89 PetscInt suppSize, Ns = 0, coneSizeA, coneSizeB, posA = -1, posB = -1; 90 PetscInt face, dim, indC[3], indS[3], seenA, flippedA, seenB, flippedB, mismatch; 91 const PetscInt *cells, *faces; 92 PetscInt cStart, cEnd, fStart, fEnd; 93 94 PetscFunctionBegin; 95 face = faceFIFO[(*fTop)++]; 96 PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells)); 97 PetscCall(ISGetPointRange(faceIS, &fStart, &fEnd, &faces)); 98 PetscCall(DMPlexGetPointDepth(dm, cells ? cells[cStart] : cStart, &dim)); 99 PetscCall(DMPlexGetSupportSize(dm, face, &suppSize)); 100 PetscCall(DMPlexGetSupport(dm, face, &supp)); 101 // Filter the support 102 for (PetscInt s = 0; s < suppSize; ++s) { 103 // Filter support 104 indC[Ns] = GetPointIndex(supp[s], cStart, cEnd, cells); 105 indS[Ns] = s; 106 if (indC[Ns] >= 0) ++Ns; 107 } 108 if (Ns < 2) PetscFunctionReturn(PETSC_SUCCESS); 109 PetscCheck(Ns == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %" PetscInt_FMT, Ns); 110 PetscCheck(indC[0] >= 0 && indC[1] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Support cells %" PetscInt_FMT " (%" PetscInt_FMT ") and %" PetscInt_FMT " (%" PetscInt_FMT ") are not both valid", supp[0], indC[0], supp[1], indC[1]); 111 seenA = PetscBTLookup(seenCells, indC[0]); 112 flippedA = PetscBTLookup(flippedCells, indC[0]) ? 1 : 0; 113 seenB = PetscBTLookup(seenCells, indC[1]); 114 flippedB = PetscBTLookup(flippedCells, indC[1]) ? 1 : 0; 115 116 PetscCall(DMPlexGetConeSize(dm, supp[indS[0]], &coneSizeA)); 117 PetscCall(DMPlexGetConeSize(dm, supp[indS[1]], &coneSizeB)); 118 PetscCall(DMPlexGetCone(dm, supp[indS[0]], &coneA)); 119 PetscCall(DMPlexGetCone(dm, supp[indS[1]], &coneB)); 120 PetscCall(DMPlexGetConeOrientation(dm, supp[indS[0]], &coneOA)); 121 PetscCall(DMPlexGetConeOrientation(dm, supp[indS[1]], &coneOB)); 122 for (PetscInt c = 0; c < coneSizeA; ++c) { 123 const PetscInt indF = GetPointIndex(coneA[c], fStart, fEnd, faces); 124 125 // Filter cone 126 if (indF < 0) continue; 127 if (!PetscBTLookup(seenFaces, indF)) { 128 faceFIFO[(*fBottom)++] = coneA[c]; 129 PetscCall(PetscBTSet(seenFaces, indF)); 130 } 131 if (coneA[c] == face) posA = c; 132 PetscCheck(*fBottom <= fEnd - fStart, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %" PetscInt_FMT " was pushed exceeding capacity %" PetscInt_FMT " > %" PetscInt_FMT, coneA[c], *fBottom, fEnd - fStart); 133 } 134 PetscCheck(posA >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " could not be located in cell %" PetscInt_FMT, face, supp[indS[0]]); 135 for (PetscInt c = 0; c < coneSizeB; ++c) { 136 const PetscInt indF = GetPointIndex(coneB[c], fStart, fEnd, faces); 137 138 // Filter cone 139 if (indF < 0) continue; 140 if (!PetscBTLookup(seenFaces, indF)) { 141 faceFIFO[(*fBottom)++] = coneB[c]; 142 PetscCall(PetscBTSet(seenFaces, indF)); 143 } 144 if (coneB[c] == face) posB = c; 145 PetscCheck(*fBottom <= fEnd - fStart, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %" PetscInt_FMT " was pushed exceeding capacity %" PetscInt_FMT " > %" PetscInt_FMT, coneA[c], *fBottom, fEnd - fStart); 146 } 147 PetscCheck(posB >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " could not be located in cell %" PetscInt_FMT, face, supp[indS[1]]); 148 149 if (dim == 1) { 150 mismatch = posA == posB; 151 } else { 152 mismatch = coneOA[posA] == coneOB[posB]; 153 } 154 155 if (mismatch ^ (flippedA ^ flippedB)) { 156 PetscCheck(!seenA || !seenB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Previously seen cells %" PetscInt_FMT " and %" PetscInt_FMT " do not match: Fault mesh is non-orientable", supp[indS[0]], supp[indS[1]]); 157 if (!seenA && !flippedA) { 158 PetscCall(PetscBTSet(flippedCells, indC[0])); 159 } else if (!seenB && !flippedB) { 160 PetscCall(PetscBTSet(flippedCells, indC[1])); 161 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable"); 162 } else PetscCheck(!mismatch || !flippedA || !flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable"); 163 PetscCall(PetscBTSet(seenCells, indC[0])); 164 PetscCall(PetscBTSet(seenCells, indC[1])); 165 PetscFunctionReturn(PETSC_SUCCESS); 166 } 167 168 static PetscErrorCode DMPlexCheckFace_Old_Internal(DM dm, PetscInt *faceFIFO, PetscInt *fTop, PetscInt *fBottom, PetscInt cStart, PetscInt fStart, PetscInt fEnd, PetscBT seenCells, PetscBT flippedCells, PetscBT seenFaces) 169 { 170 const PetscInt *support, *coneA, *coneB, *coneOA, *coneOB; 171 PetscInt supportSize, coneSizeA, coneSizeB, posA = -1, posB = -1; 172 PetscInt face, dim, seenA, flippedA, seenB, flippedB, mismatch, c; 173 174 PetscFunctionBegin; 175 face = faceFIFO[(*fTop)++]; 176 PetscCall(DMGetDimension(dm, &dim)); 177 PetscCall(DMPlexGetSupportSize(dm, face, &supportSize)); 178 PetscCall(DMPlexGetSupport(dm, face, &support)); 179 if (supportSize < 2) PetscFunctionReturn(PETSC_SUCCESS); 180 PetscCheck(supportSize == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %" PetscInt_FMT, supportSize); 181 seenA = PetscBTLookup(seenCells, support[0] - cStart); 182 flippedA = PetscBTLookup(flippedCells, support[0] - cStart) ? 1 : 0; 183 seenB = PetscBTLookup(seenCells, support[1] - cStart); 184 flippedB = PetscBTLookup(flippedCells, support[1] - cStart) ? 1 : 0; 185 186 PetscCall(DMPlexGetConeSize(dm, support[0], &coneSizeA)); 187 PetscCall(DMPlexGetConeSize(dm, support[1], &coneSizeB)); 188 PetscCall(DMPlexGetCone(dm, support[0], &coneA)); 189 PetscCall(DMPlexGetCone(dm, support[1], &coneB)); 190 PetscCall(DMPlexGetConeOrientation(dm, support[0], &coneOA)); 191 PetscCall(DMPlexGetConeOrientation(dm, support[1], &coneOB)); 192 for (c = 0; c < coneSizeA; ++c) { 193 if (!PetscBTLookup(seenFaces, coneA[c] - fStart)) { 194 faceFIFO[(*fBottom)++] = coneA[c]; 195 PetscCall(PetscBTSet(seenFaces, coneA[c] - fStart)); 196 } 197 if (coneA[c] == face) posA = c; 198 PetscCheck(*fBottom <= fEnd - fStart, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %" PetscInt_FMT " was pushed exceeding capacity %" PetscInt_FMT " > %" PetscInt_FMT, coneA[c], *fBottom, fEnd - fStart); 199 } 200 PetscCheck(posA >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " could not be located in cell %" PetscInt_FMT, face, support[0]); 201 for (c = 0; c < coneSizeB; ++c) { 202 if (!PetscBTLookup(seenFaces, coneB[c] - fStart)) { 203 faceFIFO[(*fBottom)++] = coneB[c]; 204 PetscCall(PetscBTSet(seenFaces, coneB[c] - fStart)); 205 } 206 if (coneB[c] == face) posB = c; 207 PetscCheck(*fBottom <= fEnd - fStart, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %" PetscInt_FMT " was pushed exceeding capacity %" PetscInt_FMT " > %" PetscInt_FMT, coneA[c], *fBottom, fEnd - fStart); 208 } 209 PetscCheck(posB >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " could not be located in cell %" PetscInt_FMT, face, support[1]); 210 211 if (dim == 1) { 212 mismatch = posA == posB; 213 } else { 214 mismatch = coneOA[posA] == coneOB[posB]; 215 } 216 217 if (mismatch ^ (flippedA ^ flippedB)) { 218 PetscCheck(!seenA || !seenB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Previously seen cells %" PetscInt_FMT " and %" PetscInt_FMT " do not match: Fault mesh is non-orientable", support[0], support[1]); 219 if (!seenA && !flippedA) { 220 PetscCall(PetscBTSet(flippedCells, support[0] - cStart)); 221 } else if (!seenB && !flippedB) { 222 PetscCall(PetscBTSet(flippedCells, support[1] - cStart)); 223 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable"); 224 } else PetscCheck(!mismatch || !flippedA || !flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable"); 225 PetscCall(PetscBTSet(seenCells, support[0] - cStart)); 226 PetscCall(PetscBTSet(seenCells, support[1] - cStart)); 227 PetscFunctionReturn(PETSC_SUCCESS); 228 } 229 230 /* 231 DMPlexOrient_Serial - Compute valid orientation for local connected components 232 233 Not collective 234 235 Input Parameters: 236 + dm - The `DM` 237 - cellHeight - The height of k-cells to be oriented 238 239 Output Parameters: 240 + Ncomp - The number of connected component 241 . cellComp - The connected component for each local cell 242 . faceComp - The connected component for each local face 243 - flippedCells - Marked cells should be inverted 244 245 Level: developer 246 247 .seealso: `DMPlexOrient()`` 248 */ 249 static PetscErrorCode DMPlexOrient_Serial(DM dm, IS cellIS, IS faceIS, PetscInt *Ncomp, PetscInt cellComp[], PetscInt faceComp[], PetscBT flippedCells) 250 { 251 PetscBT seenCells, seenFaces; 252 PetscInt *faceFIFO; 253 const PetscInt *cells = NULL, *faces = NULL; 254 PetscInt cStart = 0, cEnd = 0, fStart = 0, fEnd = 0; 255 256 PetscFunctionBegin; 257 if (cellIS) PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells)); 258 if (faceIS) PetscCall(ISGetPointRange(faceIS, &fStart, &fEnd, &faces)); 259 PetscCall(PetscBTCreate(cEnd - cStart, &seenCells)); 260 PetscCall(PetscBTMemzero(cEnd - cStart, seenCells)); 261 PetscCall(PetscBTCreate(fEnd - fStart, &seenFaces)); 262 PetscCall(PetscBTMemzero(fEnd - fStart, seenFaces)); 263 PetscCall(PetscMalloc1(fEnd - fStart, &faceFIFO)); 264 *Ncomp = 0; 265 for (PetscInt c = 0; c < cEnd - cStart; ++c) cellComp[c] = -1; 266 do { 267 PetscInt cc, fTop, fBottom; 268 269 // Look for first unmarked cell 270 for (cc = cStart; cc < cEnd; ++cc) 271 if (cellComp[cc - cStart] < 0) break; 272 if (cc >= cEnd) break; 273 // Initialize FIFO with first cell in component 274 { 275 const PetscInt cell = cells ? cells[cc] : cc; 276 const PetscInt *cone; 277 PetscInt coneSize; 278 279 fTop = fBottom = 0; 280 PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 281 PetscCall(DMPlexGetCone(dm, cell, &cone)); 282 for (PetscInt c = 0; c < coneSize; ++c) { 283 // Cell faces are guaranteed to be in the face set 284 faceFIFO[fBottom++] = cone[c]; 285 PetscCall(PetscBTSet(seenFaces, GetPointIndex(cone[c], fStart, fEnd, faces))); 286 } 287 PetscCall(PetscBTSet(seenCells, cc - cStart)); 288 } 289 // Consider each face in FIFO 290 while (fTop < fBottom) PetscCall(DMPlexCheckFace_Internal(dm, faceFIFO, &fTop, &fBottom, cellIS, faceIS, seenCells, flippedCells, seenFaces)); 291 // Set component for cells and faces 292 for (PetscInt c = 0; c < cEnd - cStart; ++c) { 293 if (PetscBTLookup(seenCells, c)) cellComp[c] = *Ncomp; 294 } 295 for (PetscInt f = 0; f < fEnd - fStart; ++f) { 296 if (PetscBTLookup(seenFaces, f)) faceComp[f] = *Ncomp; 297 } 298 // Wipe seenCells and seenFaces for next component 299 PetscCall(PetscBTMemzero(fEnd - fStart, seenFaces)); 300 PetscCall(PetscBTMemzero(cEnd - cStart, seenCells)); 301 ++(*Ncomp); 302 } while (1); 303 PetscCall(PetscBTDestroy(&seenCells)); 304 PetscCall(PetscBTDestroy(&seenFaces)); 305 PetscCall(PetscFree(faceFIFO)); 306 PetscFunctionReturn(PETSC_SUCCESS); 307 } 308 309 /*@ 310 DMPlexOrient - Give a consistent orientation to the input mesh 311 312 Input Parameter: 313 . dm - The `DM` 314 315 Note: 316 The orientation data for the `DM` are change in-place. 317 318 This routine will fail for non-orientable surfaces, such as the Moebius strip. 319 320 Level: advanced 321 322 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMCreate()` 323 @*/ 324 PetscErrorCode DMPlexOrient(DM dm) 325 { 326 #if 0 327 IS cellIS, faceIS; 328 329 PetscFunctionBegin; 330 PetscCall(DMPlexGetAllCells_Internal(dm, &cellIS)); 331 PetscCall(DMPlexGetAllFaces_Internal(dm, &faceIS)); 332 PetscCall(DMPlexOrientCells_Internal(dm, cellIS, faceIS)); 333 PetscCall(ISDestroy(&cellIS)); 334 PetscCall(ISDestroy(&faceIS)); 335 PetscFunctionReturn(PETSC_SUCCESS); 336 #else 337 MPI_Comm comm; 338 PetscSF sf; 339 const PetscInt *lpoints; 340 const PetscSFNode *rpoints; 341 PetscSFNode *rorntComp = NULL, *lorntComp = NULL; 342 PetscInt *numNeighbors, **neighbors, *locSupport = NULL; 343 PetscSFNode *nrankComp; 344 PetscBool *match, *flipped; 345 PetscBT seenCells, flippedCells, seenFaces; 346 PetscInt *faceFIFO, fTop, fBottom, *cellComp, *faceComp; 347 PetscInt numLeaves, numRoots, dim, h, cStart, cEnd, c, cell, fStart, fEnd, face, off, totNeighbors = 0; 348 PetscMPIInt rank, size, numComponents, comp = 0; 349 PetscBool flg, flg2; 350 PetscViewer viewer = NULL, selfviewer = NULL; 351 352 PetscFunctionBegin; 353 PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 354 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 355 PetscCallMPI(MPI_Comm_size(comm, &size)); 356 PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-orientation_view", &flg)); 357 PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-orientation_view_synchronized", &flg2)); 358 PetscCall(DMGetPointSF(dm, &sf)); 359 PetscCall(PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints)); 360 /* Truth Table 361 mismatch flips do action mismatch flipA ^ flipB action 362 F 0 flips no F F F 363 F 1 flip yes F T T 364 F 2 flips no T F T 365 T 0 flips yes T T F 366 T 1 flip no 367 T 2 flips yes 368 */ 369 PetscCall(DMGetDimension(dm, &dim)); 370 PetscCall(DMPlexGetVTKCellHeight(dm, &h)); 371 PetscCall(DMPlexGetHeightStratum(dm, h, &cStart, &cEnd)); 372 PetscCall(DMPlexGetHeightStratum(dm, h + 1, &fStart, &fEnd)); 373 PetscCall(PetscBTCreate(cEnd - cStart, &seenCells)); 374 PetscCall(PetscBTMemzero(cEnd - cStart, seenCells)); 375 PetscCall(PetscBTCreate(cEnd - cStart, &flippedCells)); 376 PetscCall(PetscBTMemzero(cEnd - cStart, flippedCells)); 377 PetscCall(PetscBTCreate(fEnd - fStart, &seenFaces)); 378 PetscCall(PetscBTMemzero(fEnd - fStart, seenFaces)); 379 PetscCall(PetscCalloc3(fEnd - fStart, &faceFIFO, cEnd - cStart, &cellComp, fEnd - fStart, &faceComp)); 380 /* 381 OLD STYLE 382 - Add an integer array over cells and faces (component) for connected component number 383 Foreach component 384 - Mark the initial cell as seen 385 - Process component as usual 386 - Set component for all seenCells 387 - Wipe seenCells and seenFaces (flippedCells can stay) 388 - Generate parallel adjacency for component using SF and seenFaces 389 - Collect numComponents adj data from each proc to 0 390 - Build same serial graph 391 - Use same solver 392 - Use Scatterv to send back flipped flags for each component 393 - Negate flippedCells by component 394 395 NEW STYLE 396 - Create the adj on each process 397 - Bootstrap to complete graph on proc 0 398 */ 399 /* Loop over components */ 400 for (cell = cStart; cell < cEnd; ++cell) cellComp[cell - cStart] = -1; 401 do { 402 /* Look for first unmarked cell */ 403 for (cell = cStart; cell < cEnd; ++cell) 404 if (cellComp[cell - cStart] < 0) break; 405 if (cell >= cEnd) break; 406 /* Initialize FIFO with first cell in component */ 407 { 408 const PetscInt *cone; 409 PetscInt coneSize; 410 411 fTop = fBottom = 0; 412 PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 413 PetscCall(DMPlexGetCone(dm, cell, &cone)); 414 for (c = 0; c < coneSize; ++c) { 415 faceFIFO[fBottom++] = cone[c]; 416 PetscCall(PetscBTSet(seenFaces, cone[c] - fStart)); 417 } 418 PetscCall(PetscBTSet(seenCells, cell - cStart)); 419 } 420 /* Consider each face in FIFO */ 421 while (fTop < fBottom) PetscCall(DMPlexCheckFace_Old_Internal(dm, faceFIFO, &fTop, &fBottom, cStart, fStart, fEnd, seenCells, flippedCells, seenFaces)); 422 /* Set component for cells and faces */ 423 for (cell = 0; cell < cEnd - cStart; ++cell) { 424 if (PetscBTLookup(seenCells, cell)) cellComp[cell] = comp; 425 } 426 for (face = 0; face < fEnd - fStart; ++face) { 427 if (PetscBTLookup(seenFaces, face)) faceComp[face] = comp; 428 } 429 /* Wipe seenCells and seenFaces for next component */ 430 PetscCall(PetscBTMemzero(fEnd - fStart, seenFaces)); 431 PetscCall(PetscBTMemzero(cEnd - cStart, seenCells)); 432 ++comp; 433 } while (1); 434 numComponents = comp; 435 if (flg) { 436 PetscViewer v; 437 438 PetscCall(PetscViewerASCIIGetStdout(comm, &v)); 439 PetscCall(PetscViewerASCIIPushSynchronized(v)); 440 PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank)); 441 PetscCall(PetscBTView(cEnd - cStart, flippedCells, v)); 442 PetscCall(PetscViewerFlush(v)); 443 PetscCall(PetscViewerASCIIPopSynchronized(v)); 444 } 445 /* Now all subdomains are oriented, but we need a consistent parallel orientation */ 446 if (numLeaves >= 0) { 447 PetscInt maxSupportSize, neighbor; 448 449 /* Store orientations of boundary faces*/ 450 PetscCall(DMPlexGetMaxSizes(dm, NULL, &maxSupportSize)); 451 PetscCall(PetscCalloc3(numRoots, &rorntComp, numRoots, &lorntComp, maxSupportSize, &locSupport)); 452 for (face = fStart; face < fEnd; ++face) { 453 const PetscInt *cone, *support, *ornt; 454 PetscInt coneSize, supportSize, Ns = 0, s, l; 455 456 PetscCall(DMPlexGetSupportSize(dm, face, &supportSize)); 457 /* Ignore overlapping cells */ 458 PetscCall(DMPlexGetSupport(dm, face, &support)); 459 for (s = 0; s < supportSize; ++s) { 460 PetscCall(PetscFindInt(support[s], numLeaves, lpoints, &l)); 461 if (l >= 0) continue; 462 locSupport[Ns++] = support[s]; 463 } 464 if (Ns != 1) continue; 465 neighbor = locSupport[0]; 466 PetscCall(DMPlexGetCone(dm, neighbor, &cone)); 467 PetscCall(DMPlexGetConeSize(dm, neighbor, &coneSize)); 468 PetscCall(DMPlexGetConeOrientation(dm, neighbor, &ornt)); 469 for (c = 0; c < coneSize; ++c) 470 if (cone[c] == face) break; 471 if (dim == 1) { 472 /* Use cone position instead, shifted to -1 or 1 */ 473 if (PetscBTLookup(flippedCells, neighbor - cStart)) rorntComp[face].rank = 1 - c * 2; 474 else rorntComp[face].rank = c * 2 - 1; 475 } else { 476 if (PetscBTLookup(flippedCells, neighbor - cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 : 1; 477 else rorntComp[face].rank = ornt[c] < 0 ? 1 : -1; 478 } 479 rorntComp[face].index = faceComp[face - fStart]; 480 } 481 /* Communicate boundary edge orientations */ 482 PetscCall(PetscSFBcastBegin(sf, MPIU_2INT, rorntComp, lorntComp, MPI_REPLACE)); 483 PetscCall(PetscSFBcastEnd(sf, MPIU_2INT, rorntComp, lorntComp, MPI_REPLACE)); 484 } 485 /* Get process adjacency */ 486 PetscCall(PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors)); 487 viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm)); 488 if (flg2) PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 489 PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &selfviewer)); 490 for (comp = 0; comp < numComponents; ++comp) { 491 PetscInt l, n; 492 493 numNeighbors[comp] = 0; 494 PetscCall(PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp])); 495 /* I know this is p^2 time in general, but for bounded degree its alright */ 496 for (l = 0; l < numLeaves; ++l) { 497 const PetscInt face = lpoints[l]; 498 499 /* Find a representative face (edge) separating pairs of procs */ 500 if ((face >= fStart) && (face < fEnd) && (faceComp[face - fStart] == comp) && rorntComp[face].rank) { 501 const PetscInt rrank = rpoints[l].rank; 502 const PetscInt rcomp = lorntComp[face].index; 503 504 for (n = 0; n < numNeighbors[comp]; ++n) 505 if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break; 506 if (n >= numNeighbors[comp]) { 507 PetscInt supportSize; 508 509 PetscCall(DMPlexGetSupportSize(dm, face, &supportSize)); 510 PetscCheck(supportSize == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %" PetscInt_FMT, supportSize); 511 if (flg) 512 PetscCall(PetscViewerASCIIPrintf(selfviewer, "[%d]: component %d, Found representative leaf %" PetscInt_FMT " (face %" PetscInt_FMT ") connecting to face %" PetscInt_FMT " on (%" PetscInt_FMT ", %" PetscInt_FMT ") with orientation %" PetscInt_FMT "\n", rank, comp, l, face, 513 rpoints[l].index, rrank, rcomp, lorntComp[face].rank)); 514 neighbors[comp][numNeighbors[comp]++] = l; 515 } 516 } 517 } 518 totNeighbors += numNeighbors[comp]; 519 } 520 PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &selfviewer)); 521 if (flg2) PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 522 PetscCall(PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match)); 523 for (comp = 0, off = 0; comp < numComponents; ++comp) { 524 PetscInt n; 525 526 for (n = 0; n < numNeighbors[comp]; ++n, ++off) { 527 const PetscInt face = lpoints[neighbors[comp][n]]; 528 const PetscInt o = rorntComp[face].rank * lorntComp[face].rank; 529 530 if (o < 0) match[off] = PETSC_TRUE; 531 else if (o > 0) match[off] = PETSC_FALSE; 532 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid face %" PetscInt_FMT " (%" PetscInt_FMT ", %" PetscInt_FMT ") neighbor: %" PetscInt_FMT " comp: %d", face, rorntComp[face].rank, lorntComp[face].rank, neighbors[comp][n], comp); 533 nrankComp[off].rank = rpoints[neighbors[comp][n]].rank; 534 nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index; 535 } 536 PetscCall(PetscFree(neighbors[comp])); 537 } 538 /* Collect the graph on 0 */ 539 if (numLeaves >= 0) { 540 Mat G; 541 PetscBT seenProcs, flippedProcs; 542 PetscInt *procFIFO, pTop, pBottom; 543 PetscInt *N = NULL, *Noff; 544 PetscSFNode *adj = NULL; 545 PetscBool *val = NULL; 546 PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc, p, o; 547 PetscMPIInt size = 0; 548 549 PetscCall(PetscCalloc1(numComponents, &flipped)); 550 if (rank == 0) PetscCallMPI(MPI_Comm_size(comm, &size)); 551 PetscCall(PetscCalloc4(size, &recvcounts, size + 1, &displs, size, &Nc, size + 1, &Noff)); 552 PetscCallMPI(MPI_Gather(&numComponents, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm)); 553 for (p = 0; p < size; ++p) displs[p + 1] = displs[p] + Nc[p]; 554 if (rank == 0) PetscCall(PetscMalloc1(displs[size], &N)); 555 PetscCallMPI(MPI_Gatherv(numNeighbors, numComponents, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm)); 556 for (p = 0, o = 0; p < size; ++p) { 557 recvcounts[p] = 0; 558 for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o]; 559 displs[p + 1] = displs[p] + recvcounts[p]; 560 } 561 if (rank == 0) PetscCall(PetscMalloc2(displs[size], &adj, displs[size], &val)); 562 PetscCallMPI(MPI_Gatherv(nrankComp, totNeighbors, MPIU_2INT, adj, recvcounts, displs, MPIU_2INT, 0, comm)); 563 PetscCallMPI(MPI_Gatherv(match, totNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm)); 564 PetscCall(PetscFree2(numNeighbors, neighbors)); 565 if (rank == 0) { 566 for (p = 1; p <= size; ++p) Noff[p] = Noff[p - 1] + Nc[p - 1]; 567 if (flg) { 568 PetscInt n; 569 570 for (p = 0, off = 0; p < size; ++p) { 571 for (c = 0; c < Nc[p]; ++c) { 572 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %" PetscInt_FMT ":\n", p, c)); 573 for (n = 0; n < N[Noff[p] + c]; ++n, ++off) PetscCall(PetscPrintf(PETSC_COMM_SELF, " edge (%" PetscInt_FMT ", %" PetscInt_FMT ") (%s):\n", adj[off].rank, adj[off].index, PetscBools[val[off]])); 574 } 575 } 576 } 577 /* Symmetrize the graph */ 578 PetscCall(MatCreate(PETSC_COMM_SELF, &G)); 579 PetscCall(MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size])); 580 PetscCall(MatSetUp(G)); 581 for (p = 0, off = 0; p < size; ++p) { 582 for (c = 0; c < Nc[p]; ++c) { 583 const PetscInt r = Noff[p] + c; 584 PetscInt n; 585 586 for (n = 0; n < N[r]; ++n, ++off) { 587 const PetscInt q = Noff[adj[off].rank] + adj[off].index; 588 const PetscScalar o = val[off] ? 1.0 : 0.0; 589 590 PetscCall(MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES)); 591 PetscCall(MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES)); 592 } 593 } 594 } 595 PetscCall(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY)); 596 PetscCall(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY)); 597 598 PetscCall(PetscBTCreate(Noff[size], &seenProcs)); 599 PetscCall(PetscBTMemzero(Noff[size], seenProcs)); 600 PetscCall(PetscBTCreate(Noff[size], &flippedProcs)); 601 PetscCall(PetscBTMemzero(Noff[size], flippedProcs)); 602 PetscCall(PetscMalloc1(Noff[size], &procFIFO)); 603 pTop = pBottom = 0; 604 for (p = 0; p < Noff[size]; ++p) { 605 if (PetscBTLookup(seenProcs, p)) continue; 606 /* Initialize FIFO with next proc */ 607 procFIFO[pBottom++] = p; 608 PetscCall(PetscBTSet(seenProcs, p)); 609 /* Consider each proc in FIFO */ 610 while (pTop < pBottom) { 611 const PetscScalar *ornt; 612 const PetscInt *neighbors; 613 PetscInt proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n; 614 615 proc = procFIFO[pTop++]; 616 flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0; 617 PetscCall(MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt)); 618 /* Loop over neighboring procs */ 619 for (n = 0; n < numNeighbors; ++n) { 620 nproc = neighbors[n]; 621 mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1; 622 seen = PetscBTLookup(seenProcs, nproc); 623 flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0; 624 625 if (mismatch ^ (flippedA ^ flippedB)) { 626 PetscCheck(!seen, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Previously seen procs %" PetscInt_FMT " and %" PetscInt_FMT " do not match: Fault mesh is non-orientable", proc, nproc); 627 if (!flippedB) { 628 PetscCall(PetscBTSet(flippedProcs, nproc)); 629 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable"); 630 } else PetscCheck(!mismatch || !flippedA || !flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable"); 631 if (!seen) { 632 procFIFO[pBottom++] = nproc; 633 PetscCall(PetscBTSet(seenProcs, nproc)); 634 } 635 } 636 } 637 } 638 PetscCall(PetscFree(procFIFO)); 639 PetscCall(MatDestroy(&G)); 640 PetscCall(PetscFree2(adj, val)); 641 PetscCall(PetscBTDestroy(&seenProcs)); 642 } 643 /* Scatter flip flags */ 644 { 645 PetscBool *flips = NULL; 646 647 if (rank == 0) { 648 PetscCall(PetscMalloc1(Noff[size], &flips)); 649 for (p = 0; p < Noff[size]; ++p) { 650 flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE; 651 if (flg && flips[p]) PetscCall(PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p)); 652 } 653 for (p = 0; p < size; ++p) displs[p + 1] = displs[p] + Nc[p]; 654 } 655 PetscCallMPI(MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, numComponents, MPIU_BOOL, 0, comm)); 656 PetscCall(PetscFree(flips)); 657 } 658 if (rank == 0) PetscCall(PetscBTDestroy(&flippedProcs)); 659 PetscCall(PetscFree(N)); 660 PetscCall(PetscFree4(recvcounts, displs, Nc, Noff)); 661 PetscCall(PetscFree2(nrankComp, match)); 662 663 /* Decide whether to flip cells in each component */ 664 for (c = 0; c < cEnd - cStart; ++c) { 665 if (flipped[cellComp[c]]) PetscCall(PetscBTNegate(flippedCells, c)); 666 } 667 PetscCall(PetscFree(flipped)); 668 } 669 if (flg) { 670 PetscViewer v; 671 672 PetscCall(PetscViewerASCIIGetStdout(comm, &v)); 673 PetscCall(PetscViewerASCIIPushSynchronized(v)); 674 PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank)); 675 PetscCall(PetscBTView(cEnd - cStart, flippedCells, v)); 676 PetscCall(PetscViewerFlush(v)); 677 PetscCall(PetscViewerASCIIPopSynchronized(v)); 678 } 679 /* Reverse flipped cells in the mesh */ 680 for (c = cStart; c < cEnd; ++c) { 681 if (PetscBTLookup(flippedCells, c - cStart)) PetscCall(DMPlexOrientPoint(dm, c, -1)); 682 } 683 PetscCall(PetscBTDestroy(&seenCells)); 684 PetscCall(PetscBTDestroy(&flippedCells)); 685 PetscCall(PetscBTDestroy(&seenFaces)); 686 PetscCall(PetscFree2(numNeighbors, neighbors)); 687 PetscCall(PetscFree3(rorntComp, lorntComp, locSupport)); 688 PetscCall(PetscFree3(faceFIFO, cellComp, faceComp)); 689 PetscFunctionReturn(PETSC_SUCCESS); 690 #endif 691 } 692 693 static PetscErrorCode CreateCellAndFaceIS_Private(DM dm, DMLabel label, IS *cellIS, IS *faceIS) 694 { 695 IS valueIS; 696 const PetscInt *values; 697 PetscInt Nv, depth = 0; 698 699 PetscFunctionBegin; 700 PetscCall(DMLabelGetValueIS(label, &valueIS)); 701 PetscCall(ISGetLocalSize(valueIS, &Nv)); 702 PetscCall(ISGetIndices(valueIS, &values)); 703 for (PetscInt v = 0; v < Nv; ++v) { 704 const PetscInt val = values[v] < 0 || values[v] >= 100 ? 0 : values[v]; 705 706 depth = PetscMax(val, depth); 707 } 708 PetscCall(ISDestroy(&valueIS)); 709 PetscCheck(depth >= 1 || !Nv, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Depth for interface must be at least 1, not %" PetscInt_FMT, depth); 710 PetscCall(DMLabelGetStratumIS(label, depth, cellIS)); 711 PetscCall(DMLabelGetStratumIS(label, depth - 1, faceIS)); 712 if (!(*cellIS)) PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, cellIS)); 713 if (!(*faceIS)) PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, faceIS)); 714 PetscFunctionReturn(PETSC_SUCCESS); 715 } 716 717 PetscErrorCode DMPlexOrientLabel(DM dm, DMLabel label) 718 { 719 IS cellIS, faceIS; 720 721 PetscFunctionBegin; 722 PetscCall(CreateCellAndFaceIS_Private(dm, label, &cellIS, &faceIS)); 723 PetscCall(DMPlexOrientCells_Internal(dm, cellIS, faceIS)); 724 PetscCall(ISDestroy(&cellIS)); 725 PetscCall(ISDestroy(&faceIS)); 726 PetscFunctionReturn(PETSC_SUCCESS); 727 } 728 729 PetscErrorCode DMPlexOrientCells_Internal(DM dm, IS cellIS, IS faceIS) 730 { 731 MPI_Comm comm; 732 PetscSF sf; 733 const PetscInt *lpoints; 734 const PetscSFNode *rpoints; 735 PetscSFNode *rorntComp = NULL, *lorntComp = NULL; 736 PetscInt *numNeighbors, **neighbors, *locSupp = NULL; 737 PetscSFNode *nrankComp; 738 PetscBool *match, *flipped; 739 PetscBT flippedCells; 740 PetscInt *cellComp, *faceComp; 741 const PetscInt *cells = NULL, *faces = NULL; 742 PetscInt cStart = 0, cEnd = 0, fStart = 0, fEnd = 0; 743 PetscInt numLeaves, numRoots, dim, Ncomp, totNeighbors = 0; 744 PetscMPIInt rank, size; 745 PetscBool view, viewSync; 746 PetscViewer viewer = NULL, selfviewer = NULL; 747 748 PetscFunctionBegin; 749 PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 750 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 751 PetscCallMPI(MPI_Comm_size(comm, &size)); 752 PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-orientation_view", &view)); 753 PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-orientation_view_synchronized", &viewSync)); 754 755 if (cellIS) PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells)); 756 if (faceIS) PetscCall(ISGetPointRange(faceIS, &fStart, &fEnd, &faces)); 757 PetscCall(DMGetPointSF(dm, &sf)); 758 PetscCall(PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints)); 759 /* Truth Table 760 mismatch flips do action mismatch flipA ^ flipB action 761 F 0 flips no F F F 762 F 1 flip yes F T T 763 F 2 flips no T F T 764 T 0 flips yes T T F 765 T 1 flip no 766 T 2 flips yes 767 */ 768 PetscCall(DMGetDimension(dm, &dim)); 769 PetscCall(PetscBTCreate(cEnd - cStart, &flippedCells)); 770 PetscCall(PetscBTMemzero(cEnd - cStart, flippedCells)); 771 PetscCall(PetscCalloc2(cEnd - cStart, &cellComp, fEnd - fStart, &faceComp)); 772 /* 773 OLD STYLE 774 - Add an integer array over cells and faces (component) for connected component number 775 Foreach component 776 - Mark the initial cell as seen 777 - Process component as usual 778 - Set component for all seenCells 779 - Wipe seenCells and seenFaces (flippedCells can stay) 780 - Generate parallel adjacency for component using SF and seenFaces 781 - Collect Ncomp adj data from each proc to 0 782 - Build same serial graph 783 - Use same solver 784 - Use Scatterv to send back flipped flags for each component 785 - Negate flippedCells by component 786 787 NEW STYLE 788 - Create the adj on each process 789 - Bootstrap to complete graph on proc 0 790 */ 791 PetscCall(DMPlexOrient_Serial(dm, cellIS, faceIS, &Ncomp, cellComp, faceComp, flippedCells)); 792 if (view) { 793 PetscViewer v; 794 795 PetscCall(PetscViewerASCIIGetStdout(comm, &v)); 796 PetscCall(PetscViewerASCIIPushSynchronized(v)); 797 PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank)); 798 PetscCall(PetscBTView(cEnd - cStart, flippedCells, v)); 799 PetscCall(PetscViewerFlush(v)); 800 PetscCall(PetscViewerASCIIPopSynchronized(v)); 801 } 802 /* Now all subdomains are oriented, but we need a consistent parallel orientation */ 803 // TODO: This all has to be rewritten to filter cones/supports to the ISes 804 if (numLeaves >= 0) { 805 PetscInt maxSuppSize, neighbor; 806 807 // Store orientations of boundary faces 808 PetscCall(DMPlexGetMaxSizes(dm, NULL, &maxSuppSize)); 809 PetscCall(PetscCalloc3(numRoots, &rorntComp, numRoots, &lorntComp, maxSuppSize, &locSupp)); 810 for (PetscInt f = fStart; f < fEnd; ++f) { 811 const PetscInt face = faces ? faces[f] : f; 812 const PetscInt *cone, *supp, *ornt; 813 PetscInt coneSize, suppSize, nind, c, Ns = 0; 814 815 PetscCall(DMPlexGetSupportSize(dm, face, &suppSize)); 816 PetscCall(DMPlexGetSupport(dm, face, &supp)); 817 for (PetscInt s = 0; s < suppSize; ++s) { 818 PetscInt ind, l; 819 820 // Filter support 821 ind = GetPointIndex(supp[s], cStart, cEnd, cells); 822 if (ind < 0) continue; 823 // Ignore overlapping cells 824 PetscCall(PetscFindInt(supp[s], numLeaves, lpoints, &l)); 825 if (l >= 0) continue; 826 locSupp[Ns++] = supp[s]; 827 } 828 PetscCheck(Ns < maxSuppSize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Index %" PetscInt_FMT " exceeds array size %" PetscInt_FMT, Ns, maxSuppSize); 829 if (Ns != 1) continue; 830 neighbor = locSupp[0]; 831 nind = GetPointIndex(neighbor, cStart, cEnd, cells); 832 PetscCall(DMPlexGetCone(dm, neighbor, &cone)); 833 PetscCall(DMPlexGetConeSize(dm, neighbor, &coneSize)); 834 PetscCall(DMPlexGetConeOrientation(dm, neighbor, &ornt)); 835 for (c = 0; c < coneSize; ++c) 836 if (cone[c] == face) break; 837 if (dim == 1) { 838 /* Use cone position instead, shifted to -1 or 1 */ 839 if (PetscBTLookup(flippedCells, nind)) rorntComp[face].rank = 1 - c * 2; 840 else rorntComp[face].rank = c * 2 - 1; 841 } else { 842 if (PetscBTLookup(flippedCells, nind)) rorntComp[face].rank = ornt[c] < 0 ? -1 : 1; 843 else rorntComp[face].rank = ornt[c] < 0 ? 1 : -1; 844 } 845 rorntComp[face].index = faceComp[GetPointIndex(face, fStart, fEnd, faces)]; 846 } 847 // Communicate boundary edge orientations 848 PetscCall(PetscSFBcastBegin(sf, MPIU_2INT, rorntComp, lorntComp, MPI_REPLACE)); 849 PetscCall(PetscSFBcastEnd(sf, MPIU_2INT, rorntComp, lorntComp, MPI_REPLACE)); 850 } 851 /* Get process adjacency */ 852 PetscCall(PetscMalloc2(Ncomp, &numNeighbors, Ncomp, &neighbors)); 853 viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm)); 854 if (viewSync) PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 855 PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &selfviewer)); 856 for (PetscInt comp = 0; comp < Ncomp; ++comp) { 857 PetscInt n; 858 859 numNeighbors[comp] = 0; 860 PetscCall(PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp])); 861 /* I know this is p^2 time in general, but for bounded degree its alright */ 862 for (PetscInt l = 0; l < numLeaves; ++l) { 863 const PetscInt face = lpoints[l]; 864 PetscInt find; 865 866 /* Find a representative face (edge) separating pairs of procs */ 867 find = GetPointIndex(face, fStart, fEnd, faces); 868 if ((find >= 0) && (faceComp[find] == comp) && rorntComp[face].rank) { 869 const PetscInt rrank = rpoints[l].rank; 870 const PetscInt rcomp = lorntComp[face].index; 871 872 for (n = 0; n < numNeighbors[comp]; ++n) 873 if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break; 874 if (n >= numNeighbors[comp]) { 875 const PetscInt *supp; 876 PetscInt suppSize, Ns = 0; 877 878 PetscCall(DMPlexGetSupport(dm, face, &supp)); 879 PetscCall(DMPlexGetSupportSize(dm, face, &suppSize)); 880 for (PetscInt s = 0; s < suppSize; ++s) { 881 // Filter support 882 if (GetPointIndex(supp[s], cStart, cEnd, cells) >= 0) ++Ns; 883 } 884 PetscCheck(Ns == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary face %" PetscInt_FMT " should see one cell, not %" PetscInt_FMT, face, Ns); 885 if (view) 886 PetscCall(PetscViewerASCIIPrintf(selfviewer, "[%d]: component %" PetscInt_FMT ", Found representative leaf %" PetscInt_FMT " (face %" PetscInt_FMT ") connecting to face %" PetscInt_FMT " on (%" PetscInt_FMT ", %" PetscInt_FMT ") with orientation %" PetscInt_FMT "\n", rank, comp, l, face, 887 rpoints[l].index, rrank, rcomp, lorntComp[face].rank)); 888 neighbors[comp][numNeighbors[comp]++] = l; 889 } 890 } 891 } 892 totNeighbors += numNeighbors[comp]; 893 } 894 PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &selfviewer)); 895 if (viewSync) PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 896 PetscCall(PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match)); 897 for (PetscInt comp = 0, off = 0; comp < Ncomp; ++comp) { 898 for (PetscInt n = 0; n < numNeighbors[comp]; ++n, ++off) { 899 const PetscInt face = lpoints[neighbors[comp][n]]; 900 const PetscInt o = rorntComp[face].rank * lorntComp[face].rank; 901 902 if (o < 0) match[off] = PETSC_TRUE; 903 else if (o > 0) match[off] = PETSC_FALSE; 904 else 905 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid face %" PetscInt_FMT " (%" PetscInt_FMT ", %" PetscInt_FMT ") neighbor: %" PetscInt_FMT " comp: %" PetscInt_FMT, face, rorntComp[face].rank, lorntComp[face].rank, neighbors[comp][n], comp); 906 nrankComp[off].rank = rpoints[neighbors[comp][n]].rank; 907 nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index; 908 } 909 PetscCall(PetscFree(neighbors[comp])); 910 } 911 /* Collect the graph on 0 */ 912 if (numLeaves >= 0) { 913 Mat G; 914 PetscBT seenProcs, flippedProcs; 915 PetscInt *procFIFO, pTop, pBottom; 916 PetscInt *N = NULL, *Noff; 917 PetscSFNode *adj = NULL; 918 PetscBool *val = NULL; 919 PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc; 920 PetscMPIInt size = 0; 921 922 PetscCall(PetscCalloc1(Ncomp, &flipped)); 923 if (rank == 0) PetscCallMPI(MPI_Comm_size(comm, &size)); 924 PetscCall(PetscCalloc4(size, &recvcounts, size + 1, &displs, size, &Nc, size + 1, &Noff)); 925 PetscCallMPI(MPI_Gather(&Ncomp, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm)); 926 for (PetscInt p = 0; p < size; ++p) displs[p + 1] = displs[p] + Nc[p]; 927 if (rank == 0) PetscCall(PetscMalloc1(displs[size], &N)); 928 PetscCallMPI(MPI_Gatherv(numNeighbors, Ncomp, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm)); 929 for (PetscInt p = 0, o = 0; p < size; ++p) { 930 recvcounts[p] = 0; 931 for (PetscInt c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o]; 932 displs[p + 1] = displs[p] + recvcounts[p]; 933 } 934 if (rank == 0) PetscCall(PetscMalloc2(displs[size], &adj, displs[size], &val)); 935 PetscCallMPI(MPI_Gatherv(nrankComp, totNeighbors, MPIU_2INT, adj, recvcounts, displs, MPIU_2INT, 0, comm)); 936 PetscCallMPI(MPI_Gatherv(match, totNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm)); 937 PetscCall(PetscFree2(numNeighbors, neighbors)); 938 if (rank == 0) { 939 for (PetscInt p = 1; p <= size; ++p) Noff[p] = Noff[p - 1] + Nc[p - 1]; 940 if (view) { 941 for (PetscInt p = 0, off = 0; p < size; ++p) { 942 for (PetscInt c = 0; c < Nc[p]; ++c) { 943 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Proc %" PetscInt_FMT " Comp %" PetscInt_FMT ":\n", p, c)); 944 for (PetscInt n = 0; n < N[Noff[p] + c]; ++n, ++off) PetscCall(PetscPrintf(PETSC_COMM_SELF, " edge (%" PetscInt_FMT ", %" PetscInt_FMT ") (%s):\n", adj[off].rank, adj[off].index, PetscBools[val[off]])); 945 } 946 } 947 } 948 /* Symmetrize the graph */ 949 PetscCall(MatCreate(PETSC_COMM_SELF, &G)); 950 PetscCall(MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size])); 951 PetscCall(MatSetUp(G)); 952 for (PetscInt p = 0, off = 0; p < size; ++p) { 953 for (PetscInt c = 0; c < Nc[p]; ++c) { 954 const PetscInt r = Noff[p] + c; 955 956 for (PetscInt n = 0; n < N[r]; ++n, ++off) { 957 const PetscInt q = Noff[adj[off].rank] + adj[off].index; 958 const PetscScalar o = val[off] ? 1.0 : 0.0; 959 960 PetscCall(MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES)); 961 PetscCall(MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES)); 962 } 963 } 964 } 965 PetscCall(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY)); 966 PetscCall(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY)); 967 968 PetscCall(PetscBTCreate(Noff[size], &seenProcs)); 969 PetscCall(PetscBTMemzero(Noff[size], seenProcs)); 970 PetscCall(PetscBTCreate(Noff[size], &flippedProcs)); 971 PetscCall(PetscBTMemzero(Noff[size], flippedProcs)); 972 PetscCall(PetscMalloc1(Noff[size], &procFIFO)); 973 pTop = pBottom = 0; 974 for (PetscInt p = 0; p < Noff[size]; ++p) { 975 if (PetscBTLookup(seenProcs, p)) continue; 976 /* Initialize FIFO with next proc */ 977 procFIFO[pBottom++] = p; 978 PetscCall(PetscBTSet(seenProcs, p)); 979 /* Consider each proc in FIFO */ 980 while (pTop < pBottom) { 981 const PetscScalar *ornt; 982 const PetscInt *neighbors; 983 PetscInt proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors; 984 985 proc = procFIFO[pTop++]; 986 flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0; 987 PetscCall(MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt)); 988 /* Loop over neighboring procs */ 989 for (PetscInt n = 0; n < numNeighbors; ++n) { 990 nproc = neighbors[n]; 991 mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1; 992 seen = PetscBTLookup(seenProcs, nproc); 993 flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0; 994 995 if (mismatch ^ (flippedA ^ flippedB)) { 996 PetscCheck(!seen, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Previously seen procs %" PetscInt_FMT " and %" PetscInt_FMT " do not match: Fault mesh is non-orientable", proc, nproc); 997 if (!flippedB) { 998 PetscCall(PetscBTSet(flippedProcs, nproc)); 999 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable"); 1000 } else PetscCheck(!mismatch || !flippedA || !flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable"); 1001 if (!seen) { 1002 procFIFO[pBottom++] = nproc; 1003 PetscCall(PetscBTSet(seenProcs, nproc)); 1004 } 1005 } 1006 } 1007 } 1008 PetscCall(PetscFree(procFIFO)); 1009 PetscCall(MatDestroy(&G)); 1010 PetscCall(PetscFree2(adj, val)); 1011 PetscCall(PetscBTDestroy(&seenProcs)); 1012 } 1013 /* Scatter flip flags */ 1014 { 1015 PetscBool *flips = NULL; 1016 1017 if (rank == 0) { 1018 PetscCall(PetscMalloc1(Noff[size], &flips)); 1019 for (PetscInt p = 0; p < Noff[size]; ++p) { 1020 flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE; 1021 if (view && flips[p]) PetscCall(PetscPrintf(comm, "Flipping Proc+Comp %" PetscInt_FMT ":\n", p)); 1022 } 1023 for (PetscInt p = 0; p < size; ++p) displs[p + 1] = displs[p] + Nc[p]; 1024 } 1025 PetscCallMPI(MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, Ncomp, MPIU_BOOL, 0, comm)); 1026 PetscCall(PetscFree(flips)); 1027 } 1028 if (rank == 0) PetscCall(PetscBTDestroy(&flippedProcs)); 1029 PetscCall(PetscFree(N)); 1030 PetscCall(PetscFree4(recvcounts, displs, Nc, Noff)); 1031 PetscCall(PetscFree2(nrankComp, match)); 1032 1033 /* Decide whether to flip cells in each component */ 1034 for (PetscInt c = 0; c < cEnd - cStart; ++c) { 1035 if (flipped[cellComp[c]]) PetscCall(PetscBTNegate(flippedCells, c)); 1036 } 1037 PetscCall(PetscFree(flipped)); 1038 } 1039 if (view) { 1040 PetscViewer v; 1041 1042 PetscCall(PetscViewerASCIIGetStdout(comm, &v)); 1043 PetscCall(PetscViewerASCIIPushSynchronized(v)); 1044 PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank)); 1045 PetscCall(PetscBTView(cEnd - cStart, flippedCells, v)); 1046 PetscCall(PetscViewerFlush(v)); 1047 PetscCall(PetscViewerASCIIPopSynchronized(v)); 1048 } 1049 // Reverse flipped cells in the mesh 1050 PetscViewer v; 1051 const PetscInt *degree = NULL; 1052 PetscInt *points; 1053 PetscInt pStart, pEnd; 1054 1055 if (view) { 1056 PetscCall(PetscViewerASCIIGetStdout(comm, &v)); 1057 PetscCall(PetscViewerASCIIPushSynchronized(v)); 1058 } 1059 PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 1060 if (numRoots >= 0) { 1061 PetscCall(PetscSFComputeDegreeBegin(sf, °ree)); 1062 PetscCall(PetscSFComputeDegreeEnd(sf, °ree)); 1063 } 1064 PetscCall(PetscCalloc1(pEnd - pStart, &points)); 1065 for (PetscInt c = cStart; c < cEnd; ++c) { 1066 if (PetscBTLookup(flippedCells, c - cStart)) { 1067 const PetscInt cell = cells ? cells[c] : c; 1068 1069 PetscCall(DMPlexOrientPoint(dm, cell, -1)); 1070 if (degree && degree[cell]) points[cell] = 1; 1071 if (view) PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]Flipping cell %" PetscInt_FMT "%s\n", rank, cell, degree && degree[cell] ? " and sending to overlap" : "")); 1072 } 1073 } 1074 // Must propagate flips for cells in the overlap 1075 if (numRoots >= 0) { 1076 PetscCall(PetscSFBcastBegin(sf, MPIU_INT, points, points, MPI_SUM)); 1077 PetscCall(PetscSFBcastEnd(sf, MPIU_INT, points, points, MPI_SUM)); 1078 } 1079 for (PetscInt c = cStart; c < cEnd; ++c) { 1080 const PetscInt cell = cells ? cells[c] : c; 1081 1082 if (points[cell] && !PetscBTLookup(flippedCells, c - cStart)) { 1083 PetscCall(DMPlexOrientPoint(dm, cell, -1)); 1084 if (view) PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]Flipping cell %" PetscInt_FMT " through overlap\n", rank, cell)); 1085 } 1086 } 1087 if (view) { 1088 PetscCall(PetscViewerFlush(v)); 1089 PetscCall(PetscViewerASCIIPopSynchronized(v)); 1090 } 1091 PetscCall(PetscFree(points)); 1092 PetscCall(PetscBTDestroy(&flippedCells)); 1093 PetscCall(PetscFree2(numNeighbors, neighbors)); 1094 PetscCall(PetscFree3(rorntComp, lorntComp, locSupp)); 1095 PetscCall(PetscFree2(cellComp, faceComp)); 1096 PetscFunctionReturn(PETSC_SUCCESS); 1097 } 1098