Actual source code: plexorient.c
petsc-3.7.3 2016-08-01
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)->options,((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: PetscViewerASCIIPushSynchronized(v);
256: PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank);
257: PetscBTView(cEnd-cStart, flippedCells, v);
258: PetscViewerFlush(v);
259: PetscViewerASCIIPopSynchronized(v);
260: }
261: /* Now all subdomains are oriented, but we need a consistent parallel orientation */
262: if (numLeaves >= 0) {
263: /* Store orientations of boundary faces*/
264: PetscCalloc2(numRoots,&rorntComp,numRoots,&lorntComp);
265: for (face = fStart; face < fEnd; ++face) {
266: const PetscInt *cone, *support, *ornt;
267: PetscInt coneSize, supportSize;
269: DMPlexGetSupportSize(dm, face, &supportSize);
270: if (supportSize != 1) continue;
271: DMPlexGetSupport(dm, face, &support);
273: DMPlexGetCone(dm, support[0], &cone);
274: DMPlexGetConeSize(dm, support[0], &coneSize);
275: DMPlexGetConeOrientation(dm, support[0], &ornt);
276: for (c = 0; c < coneSize; ++c) if (cone[c] == face) break;
277: if (dim == 1) {
278: /* Use cone position instead, shifted to -1 or 1 */
279: if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = 1-c*2;
280: else rorntComp[face].rank = c*2-1;
281: } else {
282: if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 : 1;
283: else rorntComp[face].rank = ornt[c] < 0 ? 1 : -1;
284: }
285: rorntComp[face].index = faceComp[face-fStart];
286: }
287: /* Communicate boundary edge orientations */
288: PetscSFBcastBegin(sf, MPIU_2INT, rorntComp, lorntComp);
289: PetscSFBcastEnd(sf, MPIU_2INT, rorntComp, lorntComp);
290: }
291: /* Get process adjacency */
292: PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors);
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) {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);}
314: neighbors[comp][numNeighbors[comp]++] = l;
315: }
316: }
317: }
318: totNeighbors += numNeighbors[comp];
319: }
320: PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match);
321: for (comp = 0, off = 0; comp < numComponents; ++comp) {
322: PetscInt n;
324: for (n = 0; n < numNeighbors[comp]; ++n, ++off) {
325: const PetscInt face = lpoints[neighbors[comp][n]];
326: const PetscInt o = rorntComp[face].rank*lorntComp[face].rank;
328: if (o < 0) match[off] = PETSC_TRUE;
329: else if (o > 0) match[off] = PETSC_FALSE;
330: 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);
331: nrankComp[off].rank = rpoints[neighbors[comp][n]].rank;
332: nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index;
333: }
334: PetscFree(neighbors[comp]);
335: }
336: /* Collect the graph on 0 */
337: if (numLeaves >= 0) {
338: Mat G;
339: PetscBT seenProcs, flippedProcs;
340: PetscInt *procFIFO, pTop, pBottom;
341: PetscInt *N = NULL, *Noff;
342: PetscSFNode *adj = NULL;
343: PetscBool *val = NULL;
344: PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc, p, o;
345: PetscMPIInt numProcs = 0;
347: PetscCalloc1(numComponents, &flipped);
348: if (!rank) {MPI_Comm_size(comm, &numProcs);}
349: PetscCalloc4(numProcs, &recvcounts, numProcs+1, &displs, numProcs, &Nc, numProcs+1, &Noff);
350: MPI_Gather(&numComponents, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm);
351: for (p = 0; p < numProcs; ++p) {
352: displs[p+1] = displs[p] + Nc[p];
353: }
354: if (!rank) {PetscMalloc1(displs[numProcs],&N);}
355: MPI_Gatherv(numNeighbors, numComponents, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm);
356: for (p = 0, o = 0; p < numProcs; ++p) {
357: recvcounts[p] = 0;
358: for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o];
359: displs[p+1] = displs[p] + recvcounts[p];
360: }
361: if (!rank) {PetscMalloc2(displs[numProcs], &adj, displs[numProcs], &val);}
362: MPI_Gatherv(nrankComp, totNeighbors, MPIU_2INT, adj, recvcounts, displs, MPIU_2INT, 0, comm);
363: MPI_Gatherv(match, totNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm);
364: PetscFree2(numNeighbors, neighbors);
365: if (!rank) {
366: for (p = 1; p <= numProcs; ++p) {Noff[p] = Noff[p-1] + Nc[p-1];}
367: if (flg) {
368: PetscInt n;
370: for (p = 0, off = 0; p < numProcs; ++p) {
371: for (c = 0; c < Nc[p]; ++c) {
372: PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %d:\n", p, c);
373: for (n = 0; n < N[Noff[p]+c]; ++n, ++off) {
374: PetscPrintf(PETSC_COMM_SELF, " edge (%d, %d) (%d):\n", adj[off].rank, adj[off].index, val[off]);
375: }
376: }
377: }
378: }
379: /* Symmetrize the graph */
380: MatCreate(PETSC_COMM_SELF, &G);
381: MatSetSizes(G, Noff[numProcs], Noff[numProcs], Noff[numProcs], Noff[numProcs]);
382: MatSetUp(G);
383: for (p = 0, off = 0; p < numProcs; ++p) {
384: for (c = 0; c < Nc[p]; ++c) {
385: const PetscInt r = Noff[p]+c;
386: PetscInt n;
388: for (n = 0; n < N[r]; ++n, ++off) {
389: const PetscInt q = Noff[adj[off].rank] + adj[off].index;
390: const PetscScalar o = val[off] ? 1.0 : 0.0;
392: MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES);
393: MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES);
394: }
395: }
396: }
397: MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);
398: MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);
400: PetscBTCreate(Noff[numProcs], &seenProcs);
401: PetscBTMemzero(Noff[numProcs], seenProcs);
402: PetscBTCreate(Noff[numProcs], &flippedProcs);
403: PetscBTMemzero(Noff[numProcs], flippedProcs);
404: PetscMalloc1(Noff[numProcs], &procFIFO);
405: pTop = pBottom = 0;
406: for (p = 0; p < Noff[numProcs]; ++p) {
407: if (PetscBTLookup(seenProcs, p)) continue;
408: /* Initialize FIFO with next proc */
409: procFIFO[pBottom++] = p;
410: PetscBTSet(seenProcs, p);
411: /* Consider each proc in FIFO */
412: while (pTop < pBottom) {
413: const PetscScalar *ornt;
414: const PetscInt *neighbors;
415: PetscInt proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n;
417: proc = procFIFO[pTop++];
418: flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0;
419: MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt);
420: /* Loop over neighboring procs */
421: for (n = 0; n < numNeighbors; ++n) {
422: nproc = neighbors[n];
423: mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1;
424: seen = PetscBTLookup(seenProcs, nproc);
425: flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0;
427: if (mismatch ^ (flippedA ^ flippedB)) {
428: 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);
429: if (!flippedB) {
430: PetscBTSet(flippedProcs, nproc);
431: } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
432: } else if (mismatch && flippedA && flippedB) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
433: if (!seen) {
434: procFIFO[pBottom++] = nproc;
435: PetscBTSet(seenProcs, nproc);
436: }
437: }
438: }
439: }
440: PetscFree(procFIFO);
441: MatDestroy(&G);
442: PetscFree2(adj, val);
443: PetscBTDestroy(&seenProcs);
444: }
445: /* Scatter flip flags */
446: {
447: PetscBool *flips = NULL;
449: if (!rank) {
450: PetscMalloc1(Noff[numProcs], &flips);
451: for (p = 0; p < Noff[numProcs]; ++p) {
452: flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE;
453: if (flg && flips[p]) {PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p);}
454: }
455: for (p = 0; p < numProcs; ++p) {
456: displs[p+1] = displs[p] + Nc[p];
457: }
458: }
459: MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, numComponents, MPIU_BOOL, 0, comm);
460: PetscFree(flips);
461: }
462: if (!rank) {PetscBTDestroy(&flippedProcs);}
463: PetscFree(N);
464: PetscFree4(recvcounts, displs, Nc, Noff);
465: PetscFree2(nrankComp, match);
467: /* Decide whether to flip cells in each component */
468: for (c = 0; c < cEnd-cStart; ++c) {if (flipped[cellComp[c]]) {PetscBTNegate(flippedCells, c);}}
469: PetscFree(flipped);
470: }
471: if (flg) {
472: PetscViewer v;
474: PetscViewerASCIIGetStdout(comm, &v);
475: PetscViewerASCIIPushSynchronized(v);
476: PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank);
477: PetscBTView(cEnd-cStart, flippedCells, v);
478: PetscViewerFlush(v);
479: PetscViewerASCIIPopSynchronized(v);
480: }
481: /* Reverse flipped cells in the mesh */
482: for (c = cStart; c < cEnd; ++c) {
483: if (PetscBTLookup(flippedCells, c-cStart)) {DMPlexReverseCell(dm, c);}
484: }
485: PetscBTDestroy(&seenCells);
486: PetscBTDestroy(&flippedCells);
487: PetscBTDestroy(&seenFaces);
488: PetscFree2(numNeighbors, neighbors);
489: PetscFree2(rorntComp, lorntComp);
490: PetscFree3(faceFIFO, cellComp, faceComp);
491: return(0);
492: }