Actual source code: plexorient.c
petsc-3.8.4 2018-03-24
1: #include <petsc/private/dmpleximpl.h>
2: #include <petscsf.h>
4: /*@
5: DMPlexReverseCell - Give a mesh cell the opposite orientation
7: Input Parameters:
8: + dm - The DM
9: - cell - The cell number
11: Note: The modification of the DM is done in-place.
13: Level: advanced
15: .seealso: DMPlexOrient(), DMCreate(), DMPLEX
16: @*/
17: PetscErrorCode DMPlexReverseCell(DM dm, PetscInt cell)
18: {
19: /* Note that the reverse orientation ro of a face with orientation o is:
21: ro = o >= 0 ? -(faceSize - o) : faceSize + o
23: where faceSize is the size of the cone for the face.
24: */
25: const PetscInt *cone, *coneO, *support;
26: PetscInt *revcone, *revconeO;
27: PetscInt maxConeSize, coneSize, supportSize, faceSize, cp, sp;
28: PetscErrorCode ierr;
31: DMPlexGetMaxSizes(dm, &maxConeSize, NULL);
32: DMGetWorkArray(dm, maxConeSize, PETSC_INT, &revcone);
33: DMGetWorkArray(dm, maxConeSize, PETSC_INT, &revconeO);
34: /* Reverse cone, and reverse orientations of faces */
35: DMPlexGetConeSize(dm, cell, &coneSize);
36: DMPlexGetCone(dm, cell, &cone);
37: DMPlexGetConeOrientation(dm, cell, &coneO);
38: for (cp = 0; cp < coneSize; ++cp) {
39: const PetscInt rcp = coneSize-cp-1;
41: DMPlexGetConeSize(dm, cone[rcp], &faceSize);
42: revcone[cp] = cone[rcp];
43: revconeO[cp] = coneO[rcp] >= 0 ? -(faceSize-coneO[rcp]) : faceSize+coneO[rcp];
44: }
45: DMPlexSetCone(dm, cell, revcone);
46: DMPlexSetConeOrientation(dm, cell, revconeO);
47: /* Reverse orientation of this cell in the support hypercells */
48: faceSize = coneSize;
49: DMPlexGetSupportSize(dm, cell, &supportSize);
50: DMPlexGetSupport(dm, cell, &support);
51: for (sp = 0; sp < supportSize; ++sp) {
52: DMPlexGetConeSize(dm, support[sp], &coneSize);
53: DMPlexGetCone(dm, support[sp], &cone);
54: DMPlexGetConeOrientation(dm, support[sp], &coneO);
55: for (cp = 0; cp < coneSize; ++cp) {
56: if (cone[cp] != cell) continue;
57: DMPlexInsertConeOrientation(dm, support[sp], cp, coneO[cp] >= 0 ? -(faceSize-coneO[cp]) : faceSize+coneO[cp]);
58: }
59: }
60: DMRestoreWorkArray(dm, maxConeSize, PETSC_INT, &revcone);
61: DMRestoreWorkArray(dm, maxConeSize, PETSC_INT, &revconeO);
62: return(0);
63: }
65: /*
66: - Checks face match
67: - Flips non-matching
68: - Inserts faces of support cells in FIFO
69: */
70: static PetscErrorCode DMPlexCheckFace_Internal(DM dm, PetscInt *faceFIFO, PetscInt *fTop, PetscInt *fBottom, PetscInt cStart, PetscInt fStart, PetscInt fEnd, PetscBT seenCells, PetscBT flippedCells, PetscBT seenFaces)
71: {
72: const PetscInt *support, *coneA, *coneB, *coneOA, *coneOB;
73: PetscInt supportSize, coneSizeA, coneSizeB, posA = -1, posB = -1;
74: PetscInt face, dim, seenA, flippedA, seenB, flippedB, mismatch, c;
75: PetscErrorCode ierr;
78: face = faceFIFO[(*fTop)++];
79: DMGetDimension(dm, &dim);
80: DMPlexGetSupportSize(dm, face, &supportSize);
81: DMPlexGetSupport(dm, face, &support);
82: if (supportSize < 2) return(0);
83: if (supportSize != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %d", supportSize);
84: seenA = PetscBTLookup(seenCells, support[0]-cStart);
85: flippedA = PetscBTLookup(flippedCells, support[0]-cStart) ? 1 : 0;
86: seenB = PetscBTLookup(seenCells, support[1]-cStart);
87: flippedB = PetscBTLookup(flippedCells, support[1]-cStart) ? 1 : 0;
89: DMPlexGetConeSize(dm, support[0], &coneSizeA);
90: DMPlexGetConeSize(dm, support[1], &coneSizeB);
91: DMPlexGetCone(dm, support[0], &coneA);
92: DMPlexGetCone(dm, support[1], &coneB);
93: DMPlexGetConeOrientation(dm, support[0], &coneOA);
94: DMPlexGetConeOrientation(dm, support[1], &coneOB);
95: for (c = 0; c < coneSizeA; ++c) {
96: if (!PetscBTLookup(seenFaces, coneA[c]-fStart)) {
97: faceFIFO[(*fBottom)++] = coneA[c];
98: PetscBTSet(seenFaces, coneA[c]-fStart);
99: }
100: if (coneA[c] == face) posA = c;
101: if (*fBottom > fEnd-fStart) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %d was pushed exceeding capacity %d > %d", coneA[c], *fBottom, fEnd-fStart);
102: }
103: if (posA < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d could not be located in cell %d", face, support[0]);
104: for (c = 0; c < coneSizeB; ++c) {
105: if (!PetscBTLookup(seenFaces, coneB[c]-fStart)) {
106: faceFIFO[(*fBottom)++] = coneB[c];
107: PetscBTSet(seenFaces, coneB[c]-fStart);
108: }
109: if (coneB[c] == face) posB = c;
110: if (*fBottom > fEnd-fStart) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %d was pushed exceeding capacity %d > %d", coneA[c], *fBottom, fEnd-fStart);
111: }
112: if (posB < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d could not be located in cell %d", face, support[1]);
114: if (dim == 1) {
115: mismatch = posA == posB;
116: } else {
117: mismatch = coneOA[posA] == coneOB[posB];
118: }
120: if (mismatch ^ (flippedA ^ flippedB)) {
121: 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]);
122: if (!seenA && !flippedA) {
123: PetscBTSet(flippedCells, support[0]-cStart);
124: } else if (!seenB && !flippedB) {
125: PetscBTSet(flippedCells, support[1]-cStart);
126: } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
127: } else if (mismatch && flippedA && flippedB) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
128: PetscBTSet(seenCells, support[0]-cStart);
129: PetscBTSet(seenCells, support[1]-cStart);
130: return(0);
131: }
133: /*@
134: DMPlexOrient - Give a consistent orientation to the input mesh
136: Input Parameters:
137: . dm - The DM
139: Note: The orientation data for the DM are change in-place.
140: $ This routine will fail for non-orientable surfaces, such as the Moebius strip.
142: Level: advanced
144: .seealso: DMCreate(), DMPLEX
145: @*/
146: PetscErrorCode DMPlexOrient(DM dm)
147: {
148: MPI_Comm comm;
149: PetscSF sf;
150: const PetscInt *lpoints;
151: const PetscSFNode *rpoints;
152: PetscSFNode *rorntComp = NULL, *lorntComp = NULL;
153: PetscInt *numNeighbors, **neighbors;
154: PetscSFNode *nrankComp;
155: PetscBool *match, *flipped;
156: PetscBT seenCells, flippedCells, seenFaces;
157: PetscInt *faceFIFO, fTop, fBottom, *cellComp, *faceComp;
158: PetscInt numLeaves, numRoots, dim, h, cStart, cEnd, c, cell, fStart, fEnd, face, off, totNeighbors = 0;
159: PetscMPIInt rank, size, numComponents, comp = 0;
160: PetscBool flg, flg2;
161: PetscViewer viewer = NULL, selfviewer = NULL;
162: PetscErrorCode ierr;
165: PetscObjectGetComm((PetscObject) dm, &comm);
166: MPI_Comm_rank(comm, &rank);
167: MPI_Comm_size(comm, &size);
168: PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view", &flg);
169: PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view_synchronized", &flg2);
170: DMGetPointSF(dm, &sf);
171: 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: DMGetDimension(dm, &dim);
182: DMPlexGetVTKCellHeight(dm, &h);
183: DMPlexGetHeightStratum(dm, h, &cStart, &cEnd);
184: DMPlexGetHeightStratum(dm, h+1, &fStart, &fEnd);
185: PetscBTCreate(cEnd - cStart, &seenCells);
186: PetscBTMemzero(cEnd - cStart, seenCells);
187: PetscBTCreate(cEnd - cStart, &flippedCells);
188: PetscBTMemzero(cEnd - cStart, flippedCells);
189: PetscBTCreate(fEnd - fStart, &seenFaces);
190: PetscBTMemzero(fEnd - fStart, seenFaces);
191: 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
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;
222: fTop = fBottom = 0;
223: DMPlexGetConeSize(dm, cell, &coneSize);
224: DMPlexGetCone(dm, cell, &cone);
225: for (c = 0; c < coneSize; ++c) {
226: faceFIFO[fBottom++] = cone[c];
227: PetscBTSet(seenFaces, cone[c]-fStart);
228: }
229: PetscBTSet(seenCells, cell-cStart);
230: }
231: /* Consider each face in FIFO */
232: while (fTop < fBottom) {
233: 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: PetscBTMemzero(fEnd - fStart, seenFaces);
244: PetscBTMemzero(cEnd - cStart, seenCells);
245: ++comp;
246: } while (1);
247: numComponents = comp;
248: if (flg) {
249: PetscViewer v;
251: PetscViewerASCIIGetStdout(comm, &v);
252: PetscViewerASCIIPushSynchronized(v);
253: PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank);
254: PetscBTView(cEnd-cStart, flippedCells, v);
255: PetscViewerFlush(v);
256: PetscViewerASCIIPopSynchronized(v);
257: }
258: /* Now all subdomains are oriented, but we need a consistent parallel orientation */
259: if (numLeaves >= 0) {
260: /* Store orientations of boundary faces*/
261: PetscCalloc2(numRoots,&rorntComp,numRoots,&lorntComp);
262: for (face = fStart; face < fEnd; ++face) {
263: const PetscInt *cone, *support, *ornt;
264: PetscInt coneSize, supportSize;
266: DMPlexGetSupportSize(dm, face, &supportSize);
267: if (supportSize != 1) continue;
268: DMPlexGetSupport(dm, face, &support);
270: DMPlexGetCone(dm, support[0], &cone);
271: DMPlexGetConeSize(dm, support[0], &coneSize);
272: DMPlexGetConeOrientation(dm, support[0], &ornt);
273: for (c = 0; c < coneSize; ++c) if (cone[c] == face) break;
274: if (dim == 1) {
275: /* Use cone position instead, shifted to -1 or 1 */
276: if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = 1-c*2;
277: else rorntComp[face].rank = c*2-1;
278: } else {
279: if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 : 1;
280: else rorntComp[face].rank = ornt[c] < 0 ? 1 : -1;
281: }
282: rorntComp[face].index = faceComp[face-fStart];
283: }
284: /* Communicate boundary edge orientations */
285: PetscSFBcastBegin(sf, MPIU_2INT, rorntComp, lorntComp);
286: PetscSFBcastEnd(sf, MPIU_2INT, rorntComp, lorntComp);
287: }
288: /* Get process adjacency */
289: PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors);
290: viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm));
291: if (flg2) {PetscViewerASCIIPushSynchronized(viewer);}
292: PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&selfviewer);
293: for (comp = 0; comp < numComponents; ++comp) {
294: PetscInt l, n;
296: numNeighbors[comp] = 0;
297: PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp]);
298: /* I know this is p^2 time in general, but for bounded degree its alright */
299: for (l = 0; l < numLeaves; ++l) {
300: const PetscInt face = lpoints[l];
302: /* Find a representative face (edge) separating pairs of procs */
303: if ((face >= fStart) && (face < fEnd) && (faceComp[face-fStart] == comp)) {
304: const PetscInt rrank = rpoints[l].rank;
305: const PetscInt rcomp = lorntComp[face].index;
307: for (n = 0; n < numNeighbors[comp]; ++n) if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break;
308: if (n >= numNeighbors[comp]) {
309: PetscInt supportSize;
311: DMPlexGetSupportSize(dm, face, &supportSize);
312: if (supportSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %d", supportSize);
313: 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);}
314: neighbors[comp][numNeighbors[comp]++] = l;
315: }
316: }
317: }
318: totNeighbors += numNeighbors[comp];
319: }
320: PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&selfviewer);
321: PetscViewerFlush(viewer);
322: if (flg2) {PetscViewerASCIIPopSynchronized(viewer);}
323: PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match);
324: for (comp = 0, off = 0; comp < numComponents; ++comp) {
325: PetscInt n;
327: for (n = 0; n < numNeighbors[comp]; ++n, ++off) {
328: const PetscInt face = lpoints[neighbors[comp][n]];
329: const PetscInt o = rorntComp[face].rank*lorntComp[face].rank;
331: if (o < 0) match[off] = PETSC_TRUE;
332: else if (o > 0) match[off] = PETSC_FALSE;
333: 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);
334: nrankComp[off].rank = rpoints[neighbors[comp][n]].rank;
335: nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index;
336: }
337: PetscFree(neighbors[comp]);
338: }
339: /* Collect the graph on 0 */
340: if (numLeaves >= 0) {
341: Mat G;
342: PetscBT seenProcs, flippedProcs;
343: PetscInt *procFIFO, pTop, pBottom;
344: PetscInt *N = NULL, *Noff;
345: PetscSFNode *adj = NULL;
346: PetscBool *val = NULL;
347: PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc, p, o;
348: PetscMPIInt size = 0;
350: PetscCalloc1(numComponents, &flipped);
351: if (!rank) {MPI_Comm_size(comm, &size);}
352: PetscCalloc4(size, &recvcounts, size+1, &displs, size, &Nc, size+1, &Noff);
353: MPI_Gather(&numComponents, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm);
354: for (p = 0; p < size; ++p) {
355: displs[p+1] = displs[p] + Nc[p];
356: }
357: if (!rank) {PetscMalloc1(displs[size],&N);}
358: MPI_Gatherv(numNeighbors, numComponents, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm);
359: for (p = 0, o = 0; p < size; ++p) {
360: recvcounts[p] = 0;
361: for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o];
362: displs[p+1] = displs[p] + recvcounts[p];
363: }
364: if (!rank) {PetscMalloc2(displs[size], &adj, displs[size], &val);}
365: MPI_Gatherv(nrankComp, totNeighbors, MPIU_2INT, adj, recvcounts, displs, MPIU_2INT, 0, comm);
366: MPI_Gatherv(match, totNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm);
367: PetscFree2(numNeighbors, neighbors);
368: if (!rank) {
369: for (p = 1; p <= size; ++p) {Noff[p] = Noff[p-1] + Nc[p-1];}
370: if (flg) {
371: PetscInt n;
373: for (p = 0, off = 0; p < size; ++p) {
374: for (c = 0; c < Nc[p]; ++c) {
375: PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %d:\n", p, c);
376: for (n = 0; n < N[Noff[p]+c]; ++n, ++off) {
377: PetscPrintf(PETSC_COMM_SELF, " edge (%d, %d) (%d):\n", adj[off].rank, adj[off].index, val[off]);
378: }
379: }
380: }
381: }
382: /* Symmetrize the graph */
383: MatCreate(PETSC_COMM_SELF, &G);
384: MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size]);
385: MatSetUp(G);
386: for (p = 0, off = 0; p < size; ++p) {
387: for (c = 0; c < Nc[p]; ++c) {
388: const PetscInt r = Noff[p]+c;
389: PetscInt n;
391: for (n = 0; n < N[r]; ++n, ++off) {
392: const PetscInt q = Noff[adj[off].rank] + adj[off].index;
393: const PetscScalar o = val[off] ? 1.0 : 0.0;
395: MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES);
396: MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES);
397: }
398: }
399: }
400: MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);
401: MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);
403: PetscBTCreate(Noff[size], &seenProcs);
404: PetscBTMemzero(Noff[size], seenProcs);
405: PetscBTCreate(Noff[size], &flippedProcs);
406: PetscBTMemzero(Noff[size], flippedProcs);
407: PetscMalloc1(Noff[size], &procFIFO);
408: pTop = pBottom = 0;
409: for (p = 0; p < Noff[size]; ++p) {
410: if (PetscBTLookup(seenProcs, p)) continue;
411: /* Initialize FIFO with next proc */
412: procFIFO[pBottom++] = p;
413: PetscBTSet(seenProcs, p);
414: /* Consider each proc in FIFO */
415: while (pTop < pBottom) {
416: const PetscScalar *ornt;
417: const PetscInt *neighbors;
418: PetscInt proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n;
420: proc = procFIFO[pTop++];
421: flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0;
422: MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt);
423: /* Loop over neighboring procs */
424: for (n = 0; n < numNeighbors; ++n) {
425: nproc = neighbors[n];
426: mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1;
427: seen = PetscBTLookup(seenProcs, nproc);
428: flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0;
430: if (mismatch ^ (flippedA ^ flippedB)) {
431: 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);
432: if (!flippedB) {
433: PetscBTSet(flippedProcs, nproc);
434: } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
435: } else if (mismatch && flippedA && flippedB) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
436: if (!seen) {
437: procFIFO[pBottom++] = nproc;
438: PetscBTSet(seenProcs, nproc);
439: }
440: }
441: }
442: }
443: PetscFree(procFIFO);
444: MatDestroy(&G);
445: PetscFree2(adj, val);
446: PetscBTDestroy(&seenProcs);
447: }
448: /* Scatter flip flags */
449: {
450: PetscBool *flips = NULL;
452: if (!rank) {
453: PetscMalloc1(Noff[size], &flips);
454: for (p = 0; p < Noff[size]; ++p) {
455: flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE;
456: if (flg && flips[p]) {PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p);}
457: }
458: for (p = 0; p < size; ++p) {
459: displs[p+1] = displs[p] + Nc[p];
460: }
461: }
462: MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, numComponents, MPIU_BOOL, 0, comm);
463: PetscFree(flips);
464: }
465: if (!rank) {PetscBTDestroy(&flippedProcs);}
466: PetscFree(N);
467: PetscFree4(recvcounts, displs, Nc, Noff);
468: PetscFree2(nrankComp, match);
470: /* Decide whether to flip cells in each component */
471: for (c = 0; c < cEnd-cStart; ++c) {if (flipped[cellComp[c]]) {PetscBTNegate(flippedCells, c);}}
472: PetscFree(flipped);
473: }
474: if (flg) {
475: PetscViewer v;
477: PetscViewerASCIIGetStdout(comm, &v);
478: PetscViewerASCIIPushSynchronized(v);
479: PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank);
480: PetscBTView(cEnd-cStart, flippedCells, v);
481: PetscViewerFlush(v);
482: PetscViewerASCIIPopSynchronized(v);
483: }
484: /* Reverse flipped cells in the mesh */
485: for (c = cStart; c < cEnd; ++c) {
486: if (PetscBTLookup(flippedCells, c-cStart)) {
487: DMPlexReverseCell(dm, c);
488: }
489: }
490: PetscBTDestroy(&seenCells);
491: PetscBTDestroy(&flippedCells);
492: PetscBTDestroy(&seenFaces);
493: PetscFree2(numNeighbors, neighbors);
494: PetscFree2(rorntComp, lorntComp);
495: PetscFree3(faceFIFO, cellComp, faceComp);
496: return(0);
497: }