1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 #include <petscsf.h> 3 4 /*@ 5 DMPlexCompareOrientations - Compare the cone of the given DAG point (cell) with the given reference cone (with the same cone points modulo order), and return relative orientation. 6 7 Not Collective 8 9 Input Parameters: 10 + dm - The DM (DMPLEX) 11 . p - The DAG point whose cone is compared 12 . masterConeSize - Number of the reference cone points passed (at least 2 and at most size of the cone of p) 13 - masterCone - The reference cone points 14 15 Output Parameters: 16 + start - The new starting point within the cone of p to make it conforming with the reference cone 17 - reverse - The flag whether the order of the cone points should be reversed 18 19 Level: advanced 20 21 .seealso: DMPlexOrient(), DMPlexOrientCell() 22 @*/ 23 PetscErrorCode DMPlexCompareOrientations(DM dm, PetscInt p, PetscInt masterConeSize, const PetscInt masterCone[], PetscInt *start, PetscBool *reverse) 24 { 25 PetscInt coneSize; 26 const PetscInt *cone; 27 PetscInt i, start_; 28 PetscBool reverse_; 29 PetscErrorCode ierr; 30 31 PetscFunctionBegin; 32 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 33 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 34 if (coneSize < 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %D has no cone", p); 35 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 36 if (masterConeSize < 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %D: masterConeSize must be at least 2", p); 37 if (masterConeSize > coneSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %D: masterConeSize must be at most coneSize", p); 38 start_ = 0; 39 for (i=0; i<coneSize; i++) { 40 if (cone[i] == masterCone[0]) { 41 start_ = i; 42 break; 43 } 44 } 45 if (PetscUnlikely(i==coneSize)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Point %D: starting point of reference cone not found in slave cone", p); 46 reverse_ = PETSC_FALSE; 47 for (i=0; i<masterConeSize; i++) {if (cone[(start_+i)%coneSize] != masterCone[i]) break;} 48 if (i != masterConeSize) { 49 reverse_ = PETSC_TRUE; 50 for (i=0; i<masterConeSize; i++) {if (cone[(coneSize+start_-i)%coneSize] != masterCone[i]) break;} 51 if (i < masterConeSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Point %D: cone has non-conforming order of points with respect to reference cone", p); 52 } 53 if (start) *start = start_; 54 if (reverse) *reverse = reverse_; 55 if (PetscUnlikely(cone[start_] != masterCone[0])) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D: cone[%d] = %d != %d = masterCone[0]", p, start_, cone[start_], masterCone[0]); 56 PetscFunctionReturn(0); 57 } 58 59 /*@ 60 DMPlexOrientCell - Set the desired order of cone points of this DAG point, and fix orientations accordingly. 61 62 Not Collective 63 64 Input Parameters: 65 + dm - The DM 66 . p - The DAG point (from interval given by DMPlexGetChart()) 67 . masterConeSize - Number of specified cone points (at least 2) 68 - masterCone - Specified cone points, i.e. ordered subset of current cone in DAG numbering (not cone-local numbering) 69 70 Level: intermediate 71 72 .seealso: DMPlexOrient(), DMPlexGetCone(), DMPlexGetConeOrientation(), DMPlexInterpolate(), DMPlexGetChart() 73 @*/ 74 PetscErrorCode DMPlexOrientCell(DM dm, PetscInt p, PetscInt masterConeSize, const PetscInt masterCone[]) 75 { 76 PetscInt coneSize; 77 PetscInt start1=0; 78 PetscBool reverse1=PETSC_FALSE; 79 PetscErrorCode ierr; 80 81 PetscFunctionBegin; 82 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 83 if (masterConeSize) PetscValidIntPointer(masterCone,4); 84 if (masterConeSize == 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "masterConeSize cannot be 1"); 85 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 86 if (!coneSize) PetscFunctionReturn(0); /* do nothing for points with no cone */ 87 ierr = DMPlexCompareOrientations(dm, p, masterConeSize, masterCone, &start1, &reverse1);CHKERRQ(ierr); 88 ierr = DMPlexOrientCell_Internal(dm, p, start1, reverse1);CHKERRQ(ierr); 89 #if defined(PETSC_USE_DEBUG) 90 { 91 PetscInt c; 92 const PetscInt *cone; 93 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 94 for (c = 0; c < masterConeSize; c++) { 95 if (PetscUnlikely(cone[c] != masterCone[c])) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_PLIB, "The algorithm above is wrong as cone[%d] = %d != %d = masterCone[%d]", c, cone[c], masterCone[c], c); 96 } 97 } 98 #endif 99 PetscFunctionReturn(0); 100 } 101 102 PetscErrorCode DMPlexOrientCell_Internal(DM dm, PetscInt p, PetscInt start1, PetscBool reverse1) 103 { 104 PetscInt i, j, k, maxConeSize, coneSize, coneConeSize, supportSize, supportConeSize; 105 PetscInt start0, start; 106 PetscBool reverse0, reverse; 107 PetscInt newornt; 108 const PetscInt *cone=NULL, *support=NULL, *supportCone=NULL, *ornts=NULL; 109 PetscInt *newcone=NULL, *newornts=NULL; 110 PetscErrorCode ierr; 111 112 PetscFunctionBegin; 113 if (!start1 && !reverse1) PetscFunctionReturn(0); 114 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 115 if (!coneSize) PetscFunctionReturn(0); /* do nothing for points with no cone */ 116 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 117 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr); 118 /* permute p's cone and orientations */ 119 ierr = DMPlexGetConeOrientation(dm, p, &ornts);CHKERRQ(ierr); 120 ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &newcone);CHKERRQ(ierr); 121 ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &newornts);CHKERRQ(ierr); 122 ierr = DMPlexFixFaceOrientations_Permute_Private(coneSize, cone, start1, reverse1, newcone);CHKERRQ(ierr); 123 ierr = DMPlexFixFaceOrientations_Permute_Private(coneSize, ornts, start1, reverse1, newornts);CHKERRQ(ierr); 124 /* if direction of p (face) is flipped, flip also p's cone points (edges) */ 125 if (reverse1) { 126 for (i=0; i<coneSize; i++) { 127 ierr = DMPlexGetConeSize(dm, cone[i], &coneConeSize);CHKERRQ(ierr); 128 ierr = DMPlexFixFaceOrientations_Translate_Private(newornts[i], &start0, &reverse0);CHKERRQ(ierr); 129 ierr = DMPlexFixFaceOrientations_Combine_Private(coneConeSize, start0, reverse0, 1, PETSC_FALSE, &start, &reverse);CHKERRQ(ierr); 130 ierr = DMPlexFixFaceOrientations_TranslateBack_Private(coneConeSize, start, reverse, &newornts[i]);CHKERRQ(ierr); 131 } 132 } 133 ierr = DMPlexSetConeOrientation(dm, p, newornts);CHKERRQ(ierr); 134 /* fix oriention of p within cones of p's support points */ 135 ierr = DMPlexGetSupport(dm, p, &support);CHKERRQ(ierr); 136 ierr = DMPlexGetSupportSize(dm, p, &supportSize);CHKERRQ(ierr); 137 for (j=0; j<supportSize; j++) { 138 ierr = DMPlexGetCone(dm, support[j], &supportCone);CHKERRQ(ierr); 139 ierr = DMPlexGetConeSize(dm, support[j], &supportConeSize);CHKERRQ(ierr); 140 ierr = DMPlexGetConeOrientation(dm, support[j], &ornts);CHKERRQ(ierr); 141 for (k=0; k<supportConeSize; k++) { 142 if (supportCone[k] != p) continue; 143 ierr = DMPlexFixFaceOrientations_Translate_Private(ornts[k], &start0, &reverse0);CHKERRQ(ierr); 144 ierr = DMPlexFixFaceOrientations_Combine_Private(coneSize, start0, reverse0, start1, reverse1, &start, &reverse);CHKERRQ(ierr); 145 ierr = DMPlexFixFaceOrientations_TranslateBack_Private(coneSize, start, reverse, &newornt);CHKERRQ(ierr); 146 ierr = DMPlexInsertConeOrientation(dm, support[j], k, newornt);CHKERRQ(ierr); 147 } 148 } 149 /* rewrite cone */ 150 ierr = DMPlexSetCone(dm, p, newcone);CHKERRQ(ierr); 151 ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &newcone);CHKERRQ(ierr); 152 ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &newornts);CHKERRQ(ierr); 153 PetscFunctionReturn(0); 154 } 155 156 /*@ 157 DMPlexReverseCell - Give a mesh cell the opposite orientation 158 159 Input Parameters: 160 + dm - The DM 161 - cell - The cell number 162 163 Note: The modification of the DM is done in-place. 164 165 Level: advanced 166 167 .seealso: DMPlexOrient(), DMCreate(), DMPLEX 168 @*/ 169 PetscErrorCode DMPlexReverseCell(DM dm, PetscInt cell) 170 { 171 /* Note that the reverse orientation ro of a face with orientation o is: 172 173 ro = o >= 0 ? -(faceSize - o) : faceSize + o 174 175 where faceSize is the size of the cone for the face. 176 */ 177 const PetscInt *cone, *coneO, *support; 178 PetscInt *revcone, *revconeO; 179 PetscInt maxConeSize, coneSize, supportSize, faceSize, cp, sp; 180 PetscErrorCode ierr; 181 182 PetscFunctionBegin; 183 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr); 184 ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &revcone);CHKERRQ(ierr); 185 ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &revconeO);CHKERRQ(ierr); 186 /* Reverse cone, and reverse orientations of faces */ 187 ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 188 ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr); 189 ierr = DMPlexGetConeOrientation(dm, cell, &coneO);CHKERRQ(ierr); 190 for (cp = 0; cp < coneSize; ++cp) { 191 const PetscInt rcp = coneSize-cp-1; 192 193 ierr = DMPlexGetConeSize(dm, cone[rcp], &faceSize);CHKERRQ(ierr); 194 revcone[cp] = cone[rcp]; 195 revconeO[cp] = coneO[rcp] >= 0 ? -(faceSize-coneO[rcp]) : faceSize+coneO[rcp]; 196 } 197 ierr = DMPlexSetCone(dm, cell, revcone);CHKERRQ(ierr); 198 ierr = DMPlexSetConeOrientation(dm, cell, revconeO);CHKERRQ(ierr); 199 /* Reverse orientation of this cell in the support hypercells */ 200 faceSize = coneSize; 201 ierr = DMPlexGetSupportSize(dm, cell, &supportSize);CHKERRQ(ierr); 202 ierr = DMPlexGetSupport(dm, cell, &support);CHKERRQ(ierr); 203 for (sp = 0; sp < supportSize; ++sp) { 204 ierr = DMPlexGetConeSize(dm, support[sp], &coneSize);CHKERRQ(ierr); 205 ierr = DMPlexGetCone(dm, support[sp], &cone);CHKERRQ(ierr); 206 ierr = DMPlexGetConeOrientation(dm, support[sp], &coneO);CHKERRQ(ierr); 207 for (cp = 0; cp < coneSize; ++cp) { 208 if (cone[cp] != cell) continue; 209 ierr = DMPlexInsertConeOrientation(dm, support[sp], cp, coneO[cp] >= 0 ? -(faceSize-coneO[cp]) : faceSize+coneO[cp]);CHKERRQ(ierr); 210 } 211 } 212 ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &revcone);CHKERRQ(ierr); 213 ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &revconeO);CHKERRQ(ierr); 214 PetscFunctionReturn(0); 215 } 216 217 /* 218 - Checks face match 219 - Flips non-matching 220 - Inserts faces of support cells in FIFO 221 */ 222 static PetscErrorCode DMPlexCheckFace_Internal(DM dm, PetscInt *faceFIFO, PetscInt *fTop, PetscInt *fBottom, PetscInt cStart, PetscInt fStart, PetscInt fEnd, PetscBT seenCells, PetscBT flippedCells, PetscBT seenFaces) 223 { 224 const PetscInt *support, *coneA, *coneB, *coneOA, *coneOB; 225 PetscInt supportSize, coneSizeA, coneSizeB, posA = -1, posB = -1; 226 PetscInt face, dim, seenA, flippedA, seenB, flippedB, mismatch, c; 227 PetscErrorCode ierr; 228 229 PetscFunctionBegin; 230 face = faceFIFO[(*fTop)++]; 231 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 232 ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); 233 ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr); 234 if (supportSize < 2) PetscFunctionReturn(0); 235 if (supportSize != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %d", supportSize); 236 seenA = PetscBTLookup(seenCells, support[0]-cStart); 237 flippedA = PetscBTLookup(flippedCells, support[0]-cStart) ? 1 : 0; 238 seenB = PetscBTLookup(seenCells, support[1]-cStart); 239 flippedB = PetscBTLookup(flippedCells, support[1]-cStart) ? 1 : 0; 240 241 ierr = DMPlexGetConeSize(dm, support[0], &coneSizeA);CHKERRQ(ierr); 242 ierr = DMPlexGetConeSize(dm, support[1], &coneSizeB);CHKERRQ(ierr); 243 ierr = DMPlexGetCone(dm, support[0], &coneA);CHKERRQ(ierr); 244 ierr = DMPlexGetCone(dm, support[1], &coneB);CHKERRQ(ierr); 245 ierr = DMPlexGetConeOrientation(dm, support[0], &coneOA);CHKERRQ(ierr); 246 ierr = DMPlexGetConeOrientation(dm, support[1], &coneOB);CHKERRQ(ierr); 247 for (c = 0; c < coneSizeA; ++c) { 248 if (!PetscBTLookup(seenFaces, coneA[c]-fStart)) { 249 faceFIFO[(*fBottom)++] = coneA[c]; 250 ierr = PetscBTSet(seenFaces, coneA[c]-fStart);CHKERRQ(ierr); 251 } 252 if (coneA[c] == face) posA = c; 253 if (*fBottom > fEnd-fStart) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %d was pushed exceeding capacity %d > %d", coneA[c], *fBottom, fEnd-fStart); 254 } 255 if (posA < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d could not be located in cell %d", face, support[0]); 256 for (c = 0; c < coneSizeB; ++c) { 257 if (!PetscBTLookup(seenFaces, coneB[c]-fStart)) { 258 faceFIFO[(*fBottom)++] = coneB[c]; 259 ierr = PetscBTSet(seenFaces, coneB[c]-fStart);CHKERRQ(ierr); 260 } 261 if (coneB[c] == face) posB = c; 262 if (*fBottom > fEnd-fStart) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %d was pushed exceeding capacity %d > %d", coneA[c], *fBottom, fEnd-fStart); 263 } 264 if (posB < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d could not be located in cell %d", face, support[1]); 265 266 if (dim == 1) { 267 mismatch = posA == posB; 268 } else { 269 mismatch = coneOA[posA] == coneOB[posB]; 270 } 271 272 if (mismatch ^ (flippedA ^ flippedB)) { 273 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]); 274 if (!seenA && !flippedA) { 275 ierr = PetscBTSet(flippedCells, support[0]-cStart);CHKERRQ(ierr); 276 } else if (!seenB && !flippedB) { 277 ierr = PetscBTSet(flippedCells, support[1]-cStart);CHKERRQ(ierr); 278 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable"); 279 } else if (mismatch && flippedA && flippedB) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable"); 280 ierr = PetscBTSet(seenCells, support[0]-cStart);CHKERRQ(ierr); 281 ierr = PetscBTSet(seenCells, support[1]-cStart);CHKERRQ(ierr); 282 PetscFunctionReturn(0); 283 } 284 285 /*@ 286 DMPlexOrient - Give a consistent orientation to the input mesh 287 288 Input Parameters: 289 . dm - The DM 290 291 Note: The orientation data for the DM are change in-place. 292 $ This routine will fail for non-orientable surfaces, such as the Moebius strip. 293 294 Level: advanced 295 296 .seealso: DMCreate(), DMPLEX 297 @*/ 298 PetscErrorCode DMPlexOrient(DM dm) 299 { 300 MPI_Comm comm; 301 PetscSF sf; 302 const PetscInt *lpoints; 303 const PetscSFNode *rpoints; 304 PetscSFNode *rorntComp = NULL, *lorntComp = NULL; 305 PetscInt *numNeighbors, **neighbors; 306 PetscSFNode *nrankComp; 307 PetscBool *match, *flipped; 308 PetscBT seenCells, flippedCells, seenFaces; 309 PetscInt *faceFIFO, fTop, fBottom, *cellComp, *faceComp; 310 PetscInt numLeaves, numRoots, dim, h, cStart, cEnd, c, cell, fStart, fEnd, face, off, totNeighbors = 0; 311 PetscMPIInt rank, size, numComponents, comp = 0; 312 PetscBool flg, flg2; 313 PetscViewer viewer = NULL, selfviewer = NULL; 314 PetscErrorCode ierr; 315 316 PetscFunctionBegin; 317 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 318 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 319 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 320 ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view", &flg);CHKERRQ(ierr); 321 ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view_synchronized", &flg2);CHKERRQ(ierr); 322 ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 323 ierr = PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints);CHKERRQ(ierr); 324 /* Truth Table 325 mismatch flips do action mismatch flipA ^ flipB action 326 F 0 flips no F F F 327 F 1 flip yes F T T 328 F 2 flips no T F T 329 T 0 flips yes T T F 330 T 1 flip no 331 T 2 flips yes 332 */ 333 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 334 ierr = DMPlexGetVTKCellHeight(dm, &h);CHKERRQ(ierr); 335 ierr = DMPlexGetHeightStratum(dm, h, &cStart, &cEnd);CHKERRQ(ierr); 336 ierr = DMPlexGetHeightStratum(dm, h+1, &fStart, &fEnd);CHKERRQ(ierr); 337 ierr = PetscBTCreate(cEnd - cStart, &seenCells);CHKERRQ(ierr); 338 ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr); 339 ierr = PetscBTCreate(cEnd - cStart, &flippedCells);CHKERRQ(ierr); 340 ierr = PetscBTMemzero(cEnd - cStart, flippedCells);CHKERRQ(ierr); 341 ierr = PetscBTCreate(fEnd - fStart, &seenFaces);CHKERRQ(ierr); 342 ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr); 343 ierr = PetscCalloc3(fEnd - fStart, &faceFIFO, cEnd-cStart, &cellComp, fEnd-fStart, &faceComp);CHKERRQ(ierr); 344 /* 345 OLD STYLE 346 - Add an integer array over cells and faces (component) for connected component number 347 Foreach component 348 - Mark the initial cell as seen 349 - Process component as usual 350 - Set component for all seenCells 351 - Wipe seenCells and seenFaces (flippedCells can stay) 352 - Generate parallel adjacency for component using SF and seenFaces 353 - Collect numComponents adj data from each proc to 0 354 - Build same serial graph 355 - Use same solver 356 - Use Scatterv to to send back flipped flags for each component 357 - Negate flippedCells by component 358 359 NEW STYLE 360 - Create the adj on each process 361 - Bootstrap to complete graph on proc 0 362 */ 363 /* Loop over components */ 364 for (cell = cStart; cell < cEnd; ++cell) cellComp[cell-cStart] = -1; 365 do { 366 /* Look for first unmarked cell */ 367 for (cell = cStart; cell < cEnd; ++cell) if (cellComp[cell-cStart] < 0) break; 368 if (cell >= cEnd) break; 369 /* Initialize FIFO with first cell in component */ 370 { 371 const PetscInt *cone; 372 PetscInt coneSize; 373 374 fTop = fBottom = 0; 375 ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 376 ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr); 377 for (c = 0; c < coneSize; ++c) { 378 faceFIFO[fBottom++] = cone[c]; 379 ierr = PetscBTSet(seenFaces, cone[c]-fStart);CHKERRQ(ierr); 380 } 381 ierr = PetscBTSet(seenCells, cell-cStart);CHKERRQ(ierr); 382 } 383 /* Consider each face in FIFO */ 384 while (fTop < fBottom) { 385 ierr = DMPlexCheckFace_Internal(dm, faceFIFO, &fTop, &fBottom, cStart, fStart, fEnd, seenCells, flippedCells, seenFaces);CHKERRQ(ierr); 386 } 387 /* Set component for cells and faces */ 388 for (cell = 0; cell < cEnd-cStart; ++cell) { 389 if (PetscBTLookup(seenCells, cell)) cellComp[cell] = comp; 390 } 391 for (face = 0; face < fEnd-fStart; ++face) { 392 if (PetscBTLookup(seenFaces, face)) faceComp[face] = comp; 393 } 394 /* Wipe seenCells and seenFaces for next component */ 395 ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr); 396 ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr); 397 ++comp; 398 } while (1); 399 numComponents = comp; 400 if (flg) { 401 PetscViewer v; 402 403 ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr); 404 ierr = PetscViewerASCIIPushSynchronized(v);CHKERRQ(ierr); 405 ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank);CHKERRQ(ierr); 406 ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr); 407 ierr = PetscViewerFlush(v);CHKERRQ(ierr); 408 ierr = PetscViewerASCIIPopSynchronized(v);CHKERRQ(ierr); 409 } 410 /* Now all subdomains are oriented, but we need a consistent parallel orientation */ 411 if (numLeaves >= 0) { 412 /* Store orientations of boundary faces*/ 413 ierr = PetscCalloc2(numRoots,&rorntComp,numRoots,&lorntComp);CHKERRQ(ierr); 414 for (face = fStart; face < fEnd; ++face) { 415 const PetscInt *cone, *support, *ornt; 416 PetscInt coneSize, supportSize; 417 418 ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); 419 if (supportSize != 1) continue; 420 ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr); 421 422 ierr = DMPlexGetCone(dm, support[0], &cone);CHKERRQ(ierr); 423 ierr = DMPlexGetConeSize(dm, support[0], &coneSize);CHKERRQ(ierr); 424 ierr = DMPlexGetConeOrientation(dm, support[0], &ornt);CHKERRQ(ierr); 425 for (c = 0; c < coneSize; ++c) if (cone[c] == face) break; 426 if (dim == 1) { 427 /* Use cone position instead, shifted to -1 or 1 */ 428 if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = 1-c*2; 429 else rorntComp[face].rank = c*2-1; 430 } else { 431 if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 : 1; 432 else rorntComp[face].rank = ornt[c] < 0 ? 1 : -1; 433 } 434 rorntComp[face].index = faceComp[face-fStart]; 435 } 436 /* Communicate boundary edge orientations */ 437 ierr = PetscSFBcastBegin(sf, MPIU_2INT, rorntComp, lorntComp);CHKERRQ(ierr); 438 ierr = PetscSFBcastEnd(sf, MPIU_2INT, rorntComp, lorntComp);CHKERRQ(ierr); 439 } 440 /* Get process adjacency */ 441 ierr = PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors);CHKERRQ(ierr); 442 viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm)); 443 if (flg2) {ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);} 444 ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&selfviewer);CHKERRQ(ierr); 445 for (comp = 0; comp < numComponents; ++comp) { 446 PetscInt l, n; 447 448 numNeighbors[comp] = 0; 449 ierr = PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp]);CHKERRQ(ierr); 450 /* I know this is p^2 time in general, but for bounded degree its alright */ 451 for (l = 0; l < numLeaves; ++l) { 452 const PetscInt face = lpoints[l]; 453 454 /* Find a representative face (edge) separating pairs of procs */ 455 if ((face >= fStart) && (face < fEnd) && (faceComp[face-fStart] == comp)) { 456 const PetscInt rrank = rpoints[l].rank; 457 const PetscInt rcomp = lorntComp[face].index; 458 459 for (n = 0; n < numNeighbors[comp]; ++n) if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break; 460 if (n >= numNeighbors[comp]) { 461 PetscInt supportSize; 462 463 ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); 464 if (supportSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %d", supportSize); 465 if (flg) {ierr = PetscViewerASCIIPrintf(selfviewer, "[%d]: component %d, Found representative leaf %d (face %d) connecting to face %d on (%d, %d) with orientation %d\n", rank, comp, l, face, rpoints[l].index, rrank, rcomp, lorntComp[face].rank);CHKERRQ(ierr);} 466 neighbors[comp][numNeighbors[comp]++] = l; 467 } 468 } 469 } 470 totNeighbors += numNeighbors[comp]; 471 } 472 ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&selfviewer);CHKERRQ(ierr); 473 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 474 if (flg2) {ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);} 475 ierr = PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match);CHKERRQ(ierr); 476 for (comp = 0, off = 0; comp < numComponents; ++comp) { 477 PetscInt n; 478 479 for (n = 0; n < numNeighbors[comp]; ++n, ++off) { 480 const PetscInt face = lpoints[neighbors[comp][n]]; 481 const PetscInt o = rorntComp[face].rank*lorntComp[face].rank; 482 483 if (o < 0) match[off] = PETSC_TRUE; 484 else if (o > 0) match[off] = PETSC_FALSE; 485 else SETERRQ5(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid face %d (%d, %d) neighbor: %d comp: %d", face, rorntComp[face], lorntComp[face], neighbors[comp][n], comp); 486 nrankComp[off].rank = rpoints[neighbors[comp][n]].rank; 487 nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index; 488 } 489 ierr = PetscFree(neighbors[comp]);CHKERRQ(ierr); 490 } 491 /* Collect the graph on 0 */ 492 if (numLeaves >= 0) { 493 Mat G; 494 PetscBT seenProcs, flippedProcs; 495 PetscInt *procFIFO, pTop, pBottom; 496 PetscInt *N = NULL, *Noff; 497 PetscSFNode *adj = NULL; 498 PetscBool *val = NULL; 499 PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc, p, o; 500 PetscMPIInt size = 0; 501 502 ierr = PetscCalloc1(numComponents, &flipped);CHKERRQ(ierr); 503 if (!rank) {ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);} 504 ierr = PetscCalloc4(size, &recvcounts, size+1, &displs, size, &Nc, size+1, &Noff);CHKERRQ(ierr); 505 ierr = MPI_Gather(&numComponents, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 506 for (p = 0; p < size; ++p) { 507 displs[p+1] = displs[p] + Nc[p]; 508 } 509 if (!rank) {ierr = PetscMalloc1(displs[size],&N);CHKERRQ(ierr);} 510 ierr = MPI_Gatherv(numNeighbors, numComponents, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm);CHKERRQ(ierr); 511 for (p = 0, o = 0; p < size; ++p) { 512 recvcounts[p] = 0; 513 for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o]; 514 displs[p+1] = displs[p] + recvcounts[p]; 515 } 516 if (!rank) {ierr = PetscMalloc2(displs[size], &adj, displs[size], &val);CHKERRQ(ierr);} 517 ierr = MPI_Gatherv(nrankComp, totNeighbors, MPIU_2INT, adj, recvcounts, displs, MPIU_2INT, 0, comm);CHKERRQ(ierr); 518 ierr = MPI_Gatherv(match, totNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm);CHKERRQ(ierr); 519 ierr = PetscFree2(numNeighbors, neighbors);CHKERRQ(ierr); 520 if (!rank) { 521 for (p = 1; p <= size; ++p) {Noff[p] = Noff[p-1] + Nc[p-1];} 522 if (flg) { 523 PetscInt n; 524 525 for (p = 0, off = 0; p < size; ++p) { 526 for (c = 0; c < Nc[p]; ++c) { 527 ierr = PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %d:\n", p, c);CHKERRQ(ierr); 528 for (n = 0; n < N[Noff[p]+c]; ++n, ++off) { 529 ierr = PetscPrintf(PETSC_COMM_SELF, " edge (%d, %d) (%d):\n", adj[off].rank, adj[off].index, val[off]);CHKERRQ(ierr); 530 } 531 } 532 } 533 } 534 /* Symmetrize the graph */ 535 ierr = MatCreate(PETSC_COMM_SELF, &G);CHKERRQ(ierr); 536 ierr = MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size]);CHKERRQ(ierr); 537 ierr = MatSetUp(G);CHKERRQ(ierr); 538 for (p = 0, off = 0; p < size; ++p) { 539 for (c = 0; c < Nc[p]; ++c) { 540 const PetscInt r = Noff[p]+c; 541 PetscInt n; 542 543 for (n = 0; n < N[r]; ++n, ++off) { 544 const PetscInt q = Noff[adj[off].rank] + adj[off].index; 545 const PetscScalar o = val[off] ? 1.0 : 0.0; 546 547 ierr = MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES);CHKERRQ(ierr); 548 ierr = MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES);CHKERRQ(ierr); 549 } 550 } 551 } 552 ierr = MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 553 ierr = MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 554 555 ierr = PetscBTCreate(Noff[size], &seenProcs);CHKERRQ(ierr); 556 ierr = PetscBTMemzero(Noff[size], seenProcs);CHKERRQ(ierr); 557 ierr = PetscBTCreate(Noff[size], &flippedProcs);CHKERRQ(ierr); 558 ierr = PetscBTMemzero(Noff[size], flippedProcs);CHKERRQ(ierr); 559 ierr = PetscMalloc1(Noff[size], &procFIFO);CHKERRQ(ierr); 560 pTop = pBottom = 0; 561 for (p = 0; p < Noff[size]; ++p) { 562 if (PetscBTLookup(seenProcs, p)) continue; 563 /* Initialize FIFO with next proc */ 564 procFIFO[pBottom++] = p; 565 ierr = PetscBTSet(seenProcs, p);CHKERRQ(ierr); 566 /* Consider each proc in FIFO */ 567 while (pTop < pBottom) { 568 const PetscScalar *ornt; 569 const PetscInt *neighbors; 570 PetscInt proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n; 571 572 proc = procFIFO[pTop++]; 573 flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0; 574 ierr = MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt);CHKERRQ(ierr); 575 /* Loop over neighboring procs */ 576 for (n = 0; n < numNeighbors; ++n) { 577 nproc = neighbors[n]; 578 mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1; 579 seen = PetscBTLookup(seenProcs, nproc); 580 flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0; 581 582 if (mismatch ^ (flippedA ^ flippedB)) { 583 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); 584 if (!flippedB) { 585 ierr = PetscBTSet(flippedProcs, nproc);CHKERRQ(ierr); 586 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable"); 587 } else if (mismatch && flippedA && flippedB) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable"); 588 if (!seen) { 589 procFIFO[pBottom++] = nproc; 590 ierr = PetscBTSet(seenProcs, nproc);CHKERRQ(ierr); 591 } 592 } 593 } 594 } 595 ierr = PetscFree(procFIFO);CHKERRQ(ierr); 596 ierr = MatDestroy(&G);CHKERRQ(ierr); 597 ierr = PetscFree2(adj, val);CHKERRQ(ierr); 598 ierr = PetscBTDestroy(&seenProcs);CHKERRQ(ierr); 599 } 600 /* Scatter flip flags */ 601 { 602 PetscBool *flips = NULL; 603 604 if (!rank) { 605 ierr = PetscMalloc1(Noff[size], &flips);CHKERRQ(ierr); 606 for (p = 0; p < Noff[size]; ++p) { 607 flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE; 608 if (flg && flips[p]) {ierr = PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p);CHKERRQ(ierr);} 609 } 610 for (p = 0; p < size; ++p) { 611 displs[p+1] = displs[p] + Nc[p]; 612 } 613 } 614 ierr = MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, numComponents, MPIU_BOOL, 0, comm);CHKERRQ(ierr); 615 ierr = PetscFree(flips);CHKERRQ(ierr); 616 } 617 if (!rank) {ierr = PetscBTDestroy(&flippedProcs);CHKERRQ(ierr);} 618 ierr = PetscFree(N);CHKERRQ(ierr); 619 ierr = PetscFree4(recvcounts, displs, Nc, Noff);CHKERRQ(ierr); 620 ierr = PetscFree2(nrankComp, match);CHKERRQ(ierr); 621 622 /* Decide whether to flip cells in each component */ 623 for (c = 0; c < cEnd-cStart; ++c) {if (flipped[cellComp[c]]) {ierr = PetscBTNegate(flippedCells, c);CHKERRQ(ierr);}} 624 ierr = PetscFree(flipped);CHKERRQ(ierr); 625 } 626 if (flg) { 627 PetscViewer v; 628 629 ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr); 630 ierr = PetscViewerASCIIPushSynchronized(v);CHKERRQ(ierr); 631 ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank);CHKERRQ(ierr); 632 ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr); 633 ierr = PetscViewerFlush(v);CHKERRQ(ierr); 634 ierr = PetscViewerASCIIPopSynchronized(v);CHKERRQ(ierr); 635 } 636 /* Reverse flipped cells in the mesh */ 637 for (c = cStart; c < cEnd; ++c) { 638 if (PetscBTLookup(flippedCells, c-cStart)) { 639 ierr = DMPlexReverseCell(dm, c);CHKERRQ(ierr); 640 } 641 } 642 ierr = PetscBTDestroy(&seenCells);CHKERRQ(ierr); 643 ierr = PetscBTDestroy(&flippedCells);CHKERRQ(ierr); 644 ierr = PetscBTDestroy(&seenFaces);CHKERRQ(ierr); 645 ierr = PetscFree2(numNeighbors, neighbors);CHKERRQ(ierr); 646 ierr = PetscFree2(rorntComp, lorntComp);CHKERRQ(ierr); 647 ierr = PetscFree3(faceFIFO, cellComp, faceComp);CHKERRQ(ierr); 648 PetscFunctionReturn(0); 649 } 650