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: `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 = DMPolytopeTypeGetArrangment(ct, o); 29 PetscCall(DMPlexGetConeSize(dm, p, &coneSize)); 30 PetscCall(DMPlexGetCone(dm, p, &cone)); 31 PetscCall(DMPlexGetConeOrientation(dm, p, &ornt)); 32 PetscCall(DMGetWorkArray(dm, coneSize, MPIU_INT, &newcone)); 33 PetscCall(DMGetWorkArray(dm, coneSize, MPIU_INT, &newornt)); 34 for (c = 0; c < coneSize; ++c) { 35 DMPolytopeType ft; 36 PetscInt nO; 37 38 PetscCall(DMPlexGetCellType(dm, cone[c], &ft)); 39 nO = DMPolytopeTypeGetNumArrangments(ft)/2; 40 newcone[c] = cone[arr[c*2+0]]; 41 newornt[c] = DMPolytopeTypeComposeOrientation(ft, arr[c*2+1], ornt[arr[c*2+0]]); 42 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]); 43 } 44 PetscCall(DMPlexSetCone(dm, p, newcone)); 45 PetscCall(DMPlexSetConeOrientation(dm, p, newornt)); 46 PetscCall(DMRestoreWorkArray(dm, coneSize, MPIU_INT, &newcone)); 47 PetscCall(DMRestoreWorkArray(dm, coneSize, MPIU_INT, &newornt)); 48 /* Update orientation of this point in the support points */ 49 PetscCall(DMPlexGetSupportSize(dm, p, &supportSize)); 50 PetscCall(DMPlexGetSupport(dm, p, &support)); 51 for (s = 0; s < supportSize; ++s) { 52 PetscCall(DMPlexGetConeSize(dm, support[s], &coneSize)); 53 PetscCall(DMPlexGetCone(dm, support[s], &cone)); 54 PetscCall(DMPlexGetConeOrientation(dm, support[s], &ornt)); 55 for (c = 0; c < coneSize; ++c) { 56 PetscInt po; 57 58 if (cone[c] != p) continue; 59 /* ornt[c] * 0 = target = po * o so that po = ornt[c] * o^{-1} */ 60 po = DMPolytopeTypeComposeOrientationInv(ct, ornt[c], o); 61 PetscCall(DMPlexInsertConeOrientation(dm, support[s], c, po)); 62 } 63 } 64 PetscFunctionReturn(0); 65 } 66 67 /* 68 - Checks face match 69 - Flips non-matching 70 - Inserts faces of support cells in FIFO 71 */ 72 static PetscErrorCode DMPlexCheckFace_Internal(DM dm, PetscInt *faceFIFO, PetscInt *fTop, PetscInt *fBottom, PetscInt cStart, PetscInt fStart, PetscInt fEnd, PetscBT seenCells, PetscBT flippedCells, PetscBT seenFaces) 73 { 74 const PetscInt *support, *coneA, *coneB, *coneOA, *coneOB; 75 PetscInt supportSize, coneSizeA, coneSizeB, posA = -1, posB = -1; 76 PetscInt face, dim, seenA, flippedA, seenB, flippedB, mismatch, c; 77 78 PetscFunctionBegin; 79 face = faceFIFO[(*fTop)++]; 80 PetscCall(DMGetDimension(dm, &dim)); 81 PetscCall(DMPlexGetSupportSize(dm, face, &supportSize)); 82 PetscCall(DMPlexGetSupport(dm, face, &support)); 83 if (supportSize < 2) PetscFunctionReturn(0); 84 PetscCheck(supportSize == 2,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %" PetscInt_FMT, supportSize); 85 seenA = PetscBTLookup(seenCells, support[0]-cStart); 86 flippedA = PetscBTLookup(flippedCells, support[0]-cStart) ? 1 : 0; 87 seenB = PetscBTLookup(seenCells, support[1]-cStart); 88 flippedB = PetscBTLookup(flippedCells, support[1]-cStart) ? 1 : 0; 89 90 PetscCall(DMPlexGetConeSize(dm, support[0], &coneSizeA)); 91 PetscCall(DMPlexGetConeSize(dm, support[1], &coneSizeB)); 92 PetscCall(DMPlexGetCone(dm, support[0], &coneA)); 93 PetscCall(DMPlexGetCone(dm, support[1], &coneB)); 94 PetscCall(DMPlexGetConeOrientation(dm, support[0], &coneOA)); 95 PetscCall(DMPlexGetConeOrientation(dm, support[1], &coneOB)); 96 for (c = 0; c < coneSizeA; ++c) { 97 if (!PetscBTLookup(seenFaces, coneA[c]-fStart)) { 98 faceFIFO[(*fBottom)++] = coneA[c]; 99 PetscCall(PetscBTSet(seenFaces, coneA[c]-fStart)); 100 } 101 if (coneA[c] == face) posA = c; 102 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); 103 } 104 PetscCheck(posA >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " could not be located in cell %" PetscInt_FMT, face, support[0]); 105 for (c = 0; c < coneSizeB; ++c) { 106 if (!PetscBTLookup(seenFaces, coneB[c]-fStart)) { 107 faceFIFO[(*fBottom)++] = coneB[c]; 108 PetscCall(PetscBTSet(seenFaces, coneB[c]-fStart)); 109 } 110 if (coneB[c] == face) posB = c; 111 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); 112 } 113 PetscCheck(posB >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " could not be located in cell %" PetscInt_FMT, face, support[1]); 114 115 if (dim == 1) { 116 mismatch = posA == posB; 117 } else { 118 mismatch = coneOA[posA] == coneOB[posB]; 119 } 120 121 if (mismatch ^ (flippedA ^ flippedB)) { 122 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]); 123 if (!seenA && !flippedA) { 124 PetscCall(PetscBTSet(flippedCells, support[0]-cStart)); 125 } else if (!seenB && !flippedB) { 126 PetscCall(PetscBTSet(flippedCells, support[1]-cStart)); 127 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable"); 128 } else PetscCheck(!mismatch || !flippedA || !flippedB,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable"); 129 PetscCall(PetscBTSet(seenCells, support[0]-cStart)); 130 PetscCall(PetscBTSet(seenCells, support[1]-cStart)); 131 PetscFunctionReturn(0); 132 } 133 134 /*@ 135 DMPlexOrient - Give a consistent orientation to the input mesh 136 137 Input Parameters: 138 . dm - The DM 139 140 Note: The orientation data for the DM are change in-place. 141 $ This routine will fail for non-orientable surfaces, such as the Moebius strip. 142 143 Level: advanced 144 145 .seealso: `DMCreate()`, `DMPLEX` 146 @*/ 147 PetscErrorCode DMPlexOrient(DM dm) 148 { 149 MPI_Comm comm; 150 PetscSF sf; 151 const PetscInt *lpoints; 152 const PetscSFNode *rpoints; 153 PetscSFNode *rorntComp = NULL, *lorntComp = NULL; 154 PetscInt *numNeighbors, **neighbors, *locSupport = NULL; 155 PetscSFNode *nrankComp; 156 PetscBool *match, *flipped; 157 PetscBT seenCells, flippedCells, seenFaces; 158 PetscInt *faceFIFO, fTop, fBottom, *cellComp, *faceComp; 159 PetscInt numLeaves, numRoots, dim, h, cStart, cEnd, c, cell, fStart, fEnd, face, off, totNeighbors = 0; 160 PetscMPIInt rank, size, numComponents, comp = 0; 161 PetscBool flg, flg2; 162 PetscViewer viewer = NULL, selfviewer = NULL; 163 164 PetscFunctionBegin; 165 PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 166 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 167 PetscCallMPI(MPI_Comm_size(comm, &size)); 168 PetscCall(PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view", &flg)); 169 PetscCall(PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view_synchronized", &flg2)); 170 PetscCall(DMGetPointSF(dm, &sf)); 171 PetscCall(PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints)); 172 /* Truth Table 173 mismatch flips do action mismatch flipA ^ flipB action 174 F 0 flips no F F F 175 F 1 flip yes F T T 176 F 2 flips no T F T 177 T 0 flips yes T T F 178 T 1 flip no 179 T 2 flips yes 180 */ 181 PetscCall(DMGetDimension(dm, &dim)); 182 PetscCall(DMPlexGetVTKCellHeight(dm, &h)); 183 PetscCall(DMPlexGetHeightStratum(dm, h, &cStart, &cEnd)); 184 PetscCall(DMPlexGetHeightStratum(dm, h+1, &fStart, &fEnd)); 185 PetscCall(PetscBTCreate(cEnd - cStart, &seenCells)); 186 PetscCall(PetscBTMemzero(cEnd - cStart, seenCells)); 187 PetscCall(PetscBTCreate(cEnd - cStart, &flippedCells)); 188 PetscCall(PetscBTMemzero(cEnd - cStart, flippedCells)); 189 PetscCall(PetscBTCreate(fEnd - fStart, &seenFaces)); 190 PetscCall(PetscBTMemzero(fEnd - fStart, seenFaces)); 191 PetscCall(PetscCalloc3(fEnd - fStart, &faceFIFO, cEnd-cStart, &cellComp, fEnd-fStart, &faceComp)); 192 /* 193 OLD STYLE 194 - Add an integer array over cells and faces (component) for connected component number 195 Foreach component 196 - Mark the initial cell as seen 197 - Process component as usual 198 - Set component for all seenCells 199 - Wipe seenCells and seenFaces (flippedCells can stay) 200 - Generate parallel adjacency for component using SF and seenFaces 201 - Collect numComponents adj data from each proc to 0 202 - Build same serial graph 203 - Use same solver 204 - Use Scatterv to to send back flipped flags for each component 205 - Negate flippedCells by component 206 207 NEW STYLE 208 - Create the adj on each process 209 - Bootstrap to complete graph on proc 0 210 */ 211 /* Loop over components */ 212 for (cell = cStart; cell < cEnd; ++cell) cellComp[cell-cStart] = -1; 213 do { 214 /* Look for first unmarked cell */ 215 for (cell = cStart; cell < cEnd; ++cell) if (cellComp[cell-cStart] < 0) break; 216 if (cell >= cEnd) break; 217 /* Initialize FIFO with first cell in component */ 218 { 219 const PetscInt *cone; 220 PetscInt coneSize; 221 222 fTop = fBottom = 0; 223 PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 224 PetscCall(DMPlexGetCone(dm, cell, &cone)); 225 for (c = 0; c < coneSize; ++c) { 226 faceFIFO[fBottom++] = cone[c]; 227 PetscCall(PetscBTSet(seenFaces, cone[c]-fStart)); 228 } 229 PetscCall(PetscBTSet(seenCells, cell-cStart)); 230 } 231 /* Consider each face in FIFO */ 232 while (fTop < fBottom) { 233 PetscCall(DMPlexCheckFace_Internal(dm, faceFIFO, &fTop, &fBottom, cStart, fStart, fEnd, seenCells, flippedCells, seenFaces)); 234 } 235 /* Set component for cells and faces */ 236 for (cell = 0; cell < cEnd-cStart; ++cell) { 237 if (PetscBTLookup(seenCells, cell)) cellComp[cell] = comp; 238 } 239 for (face = 0; face < fEnd-fStart; ++face) { 240 if (PetscBTLookup(seenFaces, face)) faceComp[face] = comp; 241 } 242 /* Wipe seenCells and seenFaces for next component */ 243 PetscCall(PetscBTMemzero(fEnd - fStart, seenFaces)); 244 PetscCall(PetscBTMemzero(cEnd - cStart, seenCells)); 245 ++comp; 246 } while (1); 247 numComponents = comp; 248 if (flg) { 249 PetscViewer v; 250 251 PetscCall(PetscViewerASCIIGetStdout(comm, &v)); 252 PetscCall(PetscViewerASCIIPushSynchronized(v)); 253 PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank)); 254 PetscCall(PetscBTView(cEnd-cStart, flippedCells, v)); 255 PetscCall(PetscViewerFlush(v)); 256 PetscCall(PetscViewerASCIIPopSynchronized(v)); 257 } 258 /* Now all subdomains are oriented, but we need a consistent parallel orientation */ 259 if (numLeaves >= 0) { 260 PetscInt maxSupportSize, neighbor; 261 262 /* Store orientations of boundary faces*/ 263 PetscCall(DMPlexGetMaxSizes(dm, NULL, &maxSupportSize)); 264 PetscCall(PetscCalloc3(numRoots, &rorntComp, numRoots, &lorntComp, maxSupportSize, &locSupport)); 265 for (face = fStart; face < fEnd; ++face) { 266 const PetscInt *cone, *support, *ornt; 267 PetscInt coneSize, supportSize, Ns = 0, s, l; 268 269 PetscCall(DMPlexGetSupportSize(dm, face, &supportSize)); 270 /* Ignore overlapping cells */ 271 PetscCall(DMPlexGetSupport(dm, face, &support)); 272 for (s = 0; s < supportSize; ++s) { 273 PetscCall(PetscFindInt(support[s], numLeaves, lpoints, &l)); 274 if (l >= 0) continue; 275 locSupport[Ns++] = support[s]; 276 } 277 if (Ns != 1) continue; 278 neighbor = locSupport[0]; 279 PetscCall(DMPlexGetCone(dm, neighbor, &cone)); 280 PetscCall(DMPlexGetConeSize(dm, neighbor, &coneSize)); 281 PetscCall(DMPlexGetConeOrientation(dm, neighbor, &ornt)); 282 for (c = 0; c < coneSize; ++c) if (cone[c] == face) break; 283 if (dim == 1) { 284 /* Use cone position instead, shifted to -1 or 1 */ 285 if (PetscBTLookup(flippedCells, neighbor-cStart)) rorntComp[face].rank = 1-c*2; 286 else rorntComp[face].rank = c*2-1; 287 } else { 288 if (PetscBTLookup(flippedCells, neighbor-cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 : 1; 289 else rorntComp[face].rank = ornt[c] < 0 ? 1 : -1; 290 } 291 rorntComp[face].index = faceComp[face-fStart]; 292 } 293 /* Communicate boundary edge orientations */ 294 PetscCall(PetscSFBcastBegin(sf, MPIU_2INT, rorntComp, lorntComp,MPI_REPLACE)); 295 PetscCall(PetscSFBcastEnd(sf, MPIU_2INT, rorntComp, lorntComp,MPI_REPLACE)); 296 } 297 /* Get process adjacency */ 298 PetscCall(PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors)); 299 viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm)); 300 if (flg2) PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 301 PetscCall(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&selfviewer)); 302 for (comp = 0; comp < numComponents; ++comp) { 303 PetscInt l, n; 304 305 numNeighbors[comp] = 0; 306 PetscCall(PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp])); 307 /* I know this is p^2 time in general, but for bounded degree its alright */ 308 for (l = 0; l < numLeaves; ++l) { 309 const PetscInt face = lpoints[l]; 310 311 /* Find a representative face (edge) separating pairs of procs */ 312 if ((face >= fStart) && (face < fEnd) && (faceComp[face-fStart] == comp) && rorntComp[face].rank) { 313 const PetscInt rrank = rpoints[l].rank; 314 const PetscInt rcomp = lorntComp[face].index; 315 316 for (n = 0; n < numNeighbors[comp]; ++n) if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break; 317 if (n >= numNeighbors[comp]) { 318 PetscInt supportSize; 319 320 PetscCall(DMPlexGetSupportSize(dm, face, &supportSize)); 321 PetscCheck(supportSize == 1,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %" PetscInt_FMT, supportSize); 322 if (flg) 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, rpoints[l].index, rrank, rcomp, lorntComp[face].rank)); 323 neighbors[comp][numNeighbors[comp]++] = l; 324 } 325 } 326 } 327 totNeighbors += numNeighbors[comp]; 328 } 329 PetscCall(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&selfviewer)); 330 PetscCall(PetscViewerFlush(viewer)); 331 if (flg2) PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 332 PetscCall(PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match)); 333 for (comp = 0, off = 0; comp < numComponents; ++comp) { 334 PetscInt n; 335 336 for (n = 0; n < numNeighbors[comp]; ++n, ++off) { 337 const PetscInt face = lpoints[neighbors[comp][n]]; 338 const PetscInt o = rorntComp[face].rank*lorntComp[face].rank; 339 340 if (o < 0) match[off] = PETSC_TRUE; 341 else if (o > 0) match[off] = PETSC_FALSE; 342 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); 343 nrankComp[off].rank = rpoints[neighbors[comp][n]].rank; 344 nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index; 345 } 346 PetscCall(PetscFree(neighbors[comp])); 347 } 348 /* Collect the graph on 0 */ 349 if (numLeaves >= 0) { 350 Mat G; 351 PetscBT seenProcs, flippedProcs; 352 PetscInt *procFIFO, pTop, pBottom; 353 PetscInt *N = NULL, *Noff; 354 PetscSFNode *adj = NULL; 355 PetscBool *val = NULL; 356 PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc, p, o; 357 PetscMPIInt size = 0; 358 359 PetscCall(PetscCalloc1(numComponents, &flipped)); 360 if (rank == 0) PetscCallMPI(MPI_Comm_size(comm, &size)); 361 PetscCall(PetscCalloc4(size, &recvcounts, size+1, &displs, size, &Nc, size+1, &Noff)); 362 PetscCallMPI(MPI_Gather(&numComponents, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm)); 363 for (p = 0; p < size; ++p) { 364 displs[p+1] = displs[p] + Nc[p]; 365 } 366 if (rank == 0) PetscCall(PetscMalloc1(displs[size],&N)); 367 PetscCallMPI(MPI_Gatherv(numNeighbors, numComponents, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm)); 368 for (p = 0, o = 0; p < size; ++p) { 369 recvcounts[p] = 0; 370 for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o]; 371 displs[p+1] = displs[p] + recvcounts[p]; 372 } 373 if (rank == 0) PetscCall(PetscMalloc2(displs[size], &adj, displs[size], &val)); 374 PetscCallMPI(MPI_Gatherv(nrankComp, totNeighbors, MPIU_2INT, adj, recvcounts, displs, MPIU_2INT, 0, comm)); 375 PetscCallMPI(MPI_Gatherv(match, totNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm)); 376 PetscCall(PetscFree2(numNeighbors, neighbors)); 377 if (rank == 0) { 378 for (p = 1; p <= size; ++p) {Noff[p] = Noff[p-1] + Nc[p-1];} 379 if (flg) { 380 PetscInt n; 381 382 for (p = 0, off = 0; p < size; ++p) { 383 for (c = 0; c < Nc[p]; ++c) { 384 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %" PetscInt_FMT ":\n", p, c)); 385 for (n = 0; n < N[Noff[p]+c]; ++n, ++off) { 386 PetscCall(PetscPrintf(PETSC_COMM_SELF, " edge (%" PetscInt_FMT ", %" PetscInt_FMT ") (%s):\n", adj[off].rank, adj[off].index, PetscBools[val[off]])); 387 } 388 } 389 } 390 } 391 /* Symmetrize the graph */ 392 PetscCall(MatCreate(PETSC_COMM_SELF, &G)); 393 PetscCall(MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size])); 394 PetscCall(MatSetUp(G)); 395 for (p = 0, off = 0; p < size; ++p) { 396 for (c = 0; c < Nc[p]; ++c) { 397 const PetscInt r = Noff[p]+c; 398 PetscInt n; 399 400 for (n = 0; n < N[r]; ++n, ++off) { 401 const PetscInt q = Noff[adj[off].rank] + adj[off].index; 402 const PetscScalar o = val[off] ? 1.0 : 0.0; 403 404 PetscCall(MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES)); 405 PetscCall(MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES)); 406 } 407 } 408 } 409 PetscCall(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY)); 410 PetscCall(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY)); 411 412 PetscCall(PetscBTCreate(Noff[size], &seenProcs)); 413 PetscCall(PetscBTMemzero(Noff[size], seenProcs)); 414 PetscCall(PetscBTCreate(Noff[size], &flippedProcs)); 415 PetscCall(PetscBTMemzero(Noff[size], flippedProcs)); 416 PetscCall(PetscMalloc1(Noff[size], &procFIFO)); 417 pTop = pBottom = 0; 418 for (p = 0; p < Noff[size]; ++p) { 419 if (PetscBTLookup(seenProcs, p)) continue; 420 /* Initialize FIFO with next proc */ 421 procFIFO[pBottom++] = p; 422 PetscCall(PetscBTSet(seenProcs, p)); 423 /* Consider each proc in FIFO */ 424 while (pTop < pBottom) { 425 const PetscScalar *ornt; 426 const PetscInt *neighbors; 427 PetscInt proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n; 428 429 proc = procFIFO[pTop++]; 430 flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0; 431 PetscCall(MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt)); 432 /* Loop over neighboring procs */ 433 for (n = 0; n < numNeighbors; ++n) { 434 nproc = neighbors[n]; 435 mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1; 436 seen = PetscBTLookup(seenProcs, nproc); 437 flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0; 438 439 if (mismatch ^ (flippedA ^ flippedB)) { 440 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); 441 if (!flippedB) { 442 PetscCall(PetscBTSet(flippedProcs, nproc)); 443 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable"); 444 } else PetscCheck(!mismatch || !flippedA || !flippedB,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable"); 445 if (!seen) { 446 procFIFO[pBottom++] = nproc; 447 PetscCall(PetscBTSet(seenProcs, nproc)); 448 } 449 } 450 } 451 } 452 PetscCall(PetscFree(procFIFO)); 453 PetscCall(MatDestroy(&G)); 454 PetscCall(PetscFree2(adj, val)); 455 PetscCall(PetscBTDestroy(&seenProcs)); 456 } 457 /* Scatter flip flags */ 458 { 459 PetscBool *flips = NULL; 460 461 if (rank == 0) { 462 PetscCall(PetscMalloc1(Noff[size], &flips)); 463 for (p = 0; p < Noff[size]; ++p) { 464 flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE; 465 if (flg && flips[p]) PetscCall(PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p)); 466 } 467 for (p = 0; p < size; ++p) { 468 displs[p+1] = displs[p] + Nc[p]; 469 } 470 } 471 PetscCallMPI(MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, numComponents, MPIU_BOOL, 0, comm)); 472 PetscCall(PetscFree(flips)); 473 } 474 if (rank == 0) PetscCall(PetscBTDestroy(&flippedProcs)); 475 PetscCall(PetscFree(N)); 476 PetscCall(PetscFree4(recvcounts, displs, Nc, Noff)); 477 PetscCall(PetscFree2(nrankComp, match)); 478 479 /* Decide whether to flip cells in each component */ 480 for (c = 0; c < cEnd-cStart; ++c) {if (flipped[cellComp[c]]) PetscCall(PetscBTNegate(flippedCells, c));} 481 PetscCall(PetscFree(flipped)); 482 } 483 if (flg) { 484 PetscViewer v; 485 486 PetscCall(PetscViewerASCIIGetStdout(comm, &v)); 487 PetscCall(PetscViewerASCIIPushSynchronized(v)); 488 PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank)); 489 PetscCall(PetscBTView(cEnd-cStart, flippedCells, v)); 490 PetscCall(PetscViewerFlush(v)); 491 PetscCall(PetscViewerASCIIPopSynchronized(v)); 492 } 493 /* Reverse flipped cells in the mesh */ 494 for (c = cStart; c < cEnd; ++c) { 495 if (PetscBTLookup(flippedCells, c-cStart)) { 496 PetscCall(DMPlexOrientPoint(dm, c, -1)); 497 } 498 } 499 PetscCall(PetscBTDestroy(&seenCells)); 500 PetscCall(PetscBTDestroy(&flippedCells)); 501 PetscCall(PetscBTDestroy(&seenFaces)); 502 PetscCall(PetscFree2(numNeighbors, neighbors)); 503 PetscCall(PetscFree3(rorntComp, lorntComp, locSupport)); 504 PetscCall(PetscFree3(faceFIFO, cellComp, faceComp)); 505 PetscFunctionReturn(0); 506 } 507