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