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