Actual source code: plexorient.c

  1: #include <petsc/private/dmpleximpl.h>
  2: #include <petscsf.h>

  4: /*@
  5:   DMPlexOrientPoint - Act with the given orientation on the cone points of this mesh point, and update its use in the mesh.

  7:   Not Collective

  9:   Input Parameters:
 10: + dm - The DM
 11: . p  - The mesh point
 12: - o  - The orientation

 14:   Level: intermediate

 16: .seealso: `DMPlexOrient()`, `DMPlexGetCone()`, `DMPlexGetConeOrientation()`, `DMPlexInterpolate()`, `DMPlexGetChart()`
 17: @*/
 18: PetscErrorCode DMPlexOrientPoint(DM dm, PetscInt p, PetscInt o)
 19: {
 20:   DMPolytopeType  ct;
 21:   const PetscInt *arr, *cone, *ornt, *support;
 22:   PetscInt       *newcone, *newornt;
 23:   PetscInt        coneSize, c, supportSize, s;

 26:   DMPlexGetCellType(dm, p, &ct);
 27:   arr = DMPolytopeTypeGetArrangment(ct, o);
 28:   DMPlexGetConeSize(dm, p, &coneSize);
 29:   DMPlexGetCone(dm, p, &cone);
 30:   DMPlexGetConeOrientation(dm, p, &ornt);
 31:   DMGetWorkArray(dm, coneSize, MPIU_INT, &newcone);
 32:   DMGetWorkArray(dm, coneSize, MPIU_INT, &newornt);
 33:   for (c = 0; c < coneSize; ++c) {
 34:     DMPolytopeType ft;
 35:     PetscInt       nO;

 37:     DMPlexGetCellType(dm, cone[c], &ft);
 38:     nO         = DMPolytopeTypeGetNumArrangments(ft) / 2;
 39:     newcone[c] = cone[arr[c * 2 + 0]];
 40:     newornt[c] = DMPolytopeTypeComposeOrientation(ft, arr[c * 2 + 1], ornt[arr[c * 2 + 0]]);
 42:   }
 43:   DMPlexSetCone(dm, p, newcone);
 44:   DMPlexSetConeOrientation(dm, p, newornt);
 45:   DMRestoreWorkArray(dm, coneSize, MPIU_INT, &newcone);
 46:   DMRestoreWorkArray(dm, coneSize, MPIU_INT, &newornt);
 47:   /* Update orientation of this point in the support points */
 48:   DMPlexGetSupportSize(dm, p, &supportSize);
 49:   DMPlexGetSupport(dm, p, &support);
 50:   for (s = 0; s < supportSize; ++s) {
 51:     DMPlexGetConeSize(dm, support[s], &coneSize);
 52:     DMPlexGetCone(dm, support[s], &cone);
 53:     DMPlexGetConeOrientation(dm, support[s], &ornt);
 54:     for (c = 0; c < coneSize; ++c) {
 55:       PetscInt po;

 57:       if (cone[c] != p) continue;
 58:       /* ornt[c] * 0 = target = po * o so that po = ornt[c] * o^{-1} */
 59:       po = DMPolytopeTypeComposeOrientationInv(ct, ornt[c], o);
 60:       DMPlexInsertConeOrientation(dm, support[s], c, po);
 61:     }
 62:   }
 63:   return 0;
 64: }

 66: /*
 67:   - Checks face match
 68:     - Flips non-matching
 69:   - Inserts faces of support cells in FIFO
 70: */
 71: static PetscErrorCode DMPlexCheckFace_Internal(DM dm, PetscInt *faceFIFO, PetscInt *fTop, PetscInt *fBottom, PetscInt cStart, PetscInt fStart, PetscInt fEnd, PetscBT seenCells, PetscBT flippedCells, PetscBT seenFaces)
 72: {
 73:   const PetscInt *support, *coneA, *coneB, *coneOA, *coneOB;
 74:   PetscInt        supportSize, coneSizeA, coneSizeB, posA = -1, posB = -1;
 75:   PetscInt        face, dim, seenA, flippedA, seenB, flippedB, mismatch, c;

 77:   face = faceFIFO[(*fTop)++];
 78:   DMGetDimension(dm, &dim);
 79:   DMPlexGetSupportSize(dm, face, &supportSize);
 80:   DMPlexGetSupport(dm, face, &support);
 81:   if (supportSize < 2) return 0;
 83:   seenA    = PetscBTLookup(seenCells, support[0] - cStart);
 84:   flippedA = PetscBTLookup(flippedCells, support[0] - cStart) ? 1 : 0;
 85:   seenB    = PetscBTLookup(seenCells, support[1] - cStart);
 86:   flippedB = PetscBTLookup(flippedCells, support[1] - cStart) ? 1 : 0;

 88:   DMPlexGetConeSize(dm, support[0], &coneSizeA);
 89:   DMPlexGetConeSize(dm, support[1], &coneSizeB);
 90:   DMPlexGetCone(dm, support[0], &coneA);
 91:   DMPlexGetCone(dm, support[1], &coneB);
 92:   DMPlexGetConeOrientation(dm, support[0], &coneOA);
 93:   DMPlexGetConeOrientation(dm, support[1], &coneOB);
 94:   for (c = 0; c < coneSizeA; ++c) {
 95:     if (!PetscBTLookup(seenFaces, coneA[c] - fStart)) {
 96:       faceFIFO[(*fBottom)++] = coneA[c];
 97:       PetscBTSet(seenFaces, coneA[c] - fStart);
 98:     }
 99:     if (coneA[c] == face) posA = c;
101:   }
103:   for (c = 0; c < coneSizeB; ++c) {
104:     if (!PetscBTLookup(seenFaces, coneB[c] - fStart)) {
105:       faceFIFO[(*fBottom)++] = coneB[c];
106:       PetscBTSet(seenFaces, coneB[c] - fStart);
107:     }
108:     if (coneB[c] == face) posB = c;
110:   }

113:   if (dim == 1) {
114:     mismatch = posA == posB;
115:   } else {
116:     mismatch = coneOA[posA] == coneOB[posB];
117:   }

119:   if (mismatch ^ (flippedA ^ flippedB)) {
121:     if (!seenA && !flippedA) {
122:       PetscBTSet(flippedCells, support[0] - cStart);
123:     } else if (!seenB && !flippedB) {
124:       PetscBTSet(flippedCells, support[1] - cStart);
125:     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
127:   PetscBTSet(seenCells, support[0] - cStart);
128:   PetscBTSet(seenCells, support[1] - cStart);
129:   return 0;
130: }

132: /*@
133:   DMPlexOrient - Give a consistent orientation to the input mesh

135:   Input Parameters:
136: . dm - The DM

138:   Note: The orientation data for the DM are change in-place.
139: $ This routine will fail for non-orientable surfaces, such as the Moebius strip.

141:   Level: advanced

143: .seealso: `DMCreate()`, `DMPLEX`
144: @*/
145: PetscErrorCode DMPlexOrient(DM dm)
146: {
147:   MPI_Comm           comm;
148:   PetscSF            sf;
149:   const PetscInt    *lpoints;
150:   const PetscSFNode *rpoints;
151:   PetscSFNode       *rorntComp = NULL, *lorntComp = NULL;
152:   PetscInt          *numNeighbors, **neighbors, *locSupport = NULL;
153:   PetscSFNode       *nrankComp;
154:   PetscBool         *match, *flipped;
155:   PetscBT            seenCells, flippedCells, seenFaces;
156:   PetscInt          *faceFIFO, fTop, fBottom, *cellComp, *faceComp;
157:   PetscInt           numLeaves, numRoots, dim, h, cStart, cEnd, c, cell, fStart, fEnd, face, off, totNeighbors = 0;
158:   PetscMPIInt        rank, size, numComponents, comp = 0;
159:   PetscBool          flg, flg2;
160:   PetscViewer        viewer = NULL, selfviewer = NULL;

162:   PetscObjectGetComm((PetscObject)dm, &comm);
163:   MPI_Comm_rank(comm, &rank);
164:   MPI_Comm_size(comm, &size);
165:   PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-orientation_view", &flg);
166:   PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-orientation_view_synchronized", &flg2);
167:   DMGetPointSF(dm, &sf);
168:   PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints);
169:   /* Truth Table
170:      mismatch    flips   do action   mismatch   flipA ^ flipB   action
171:          F       0 flips     no         F             F           F
172:          F       1 flip      yes        F             T           T
173:          F       2 flips     no         T             F           T
174:          T       0 flips     yes        T             T           F
175:          T       1 flip      no
176:          T       2 flips     yes
177:   */
178:   DMGetDimension(dm, &dim);
179:   DMPlexGetVTKCellHeight(dm, &h);
180:   DMPlexGetHeightStratum(dm, h, &cStart, &cEnd);
181:   DMPlexGetHeightStratum(dm, h + 1, &fStart, &fEnd);
182:   PetscBTCreate(cEnd - cStart, &seenCells);
183:   PetscBTMemzero(cEnd - cStart, seenCells);
184:   PetscBTCreate(cEnd - cStart, &flippedCells);
185:   PetscBTMemzero(cEnd - cStart, flippedCells);
186:   PetscBTCreate(fEnd - fStart, &seenFaces);
187:   PetscBTMemzero(fEnd - fStart, seenFaces);
188:   PetscCalloc3(fEnd - fStart, &faceFIFO, cEnd - cStart, &cellComp, fEnd - fStart, &faceComp);
189:   /*
190:    OLD STYLE
191:    - Add an integer array over cells and faces (component) for connected component number
192:    Foreach component
193:      - Mark the initial cell as seen
194:      - Process component as usual
195:      - Set component for all seenCells
196:      - Wipe seenCells and seenFaces (flippedCells can stay)
197:    - Generate parallel adjacency for component using SF and seenFaces
198:    - Collect numComponents adj data from each proc to 0
199:    - Build same serial graph
200:    - Use same solver
201:    - Use Scatterv to to send back flipped flags for each component
202:    - Negate flippedCells by component

204:    NEW STYLE
205:    - Create the adj on each process
206:    - Bootstrap to complete graph on proc 0
207:   */
208:   /* Loop over components */
209:   for (cell = cStart; cell < cEnd; ++cell) cellComp[cell - cStart] = -1;
210:   do {
211:     /* Look for first unmarked cell */
212:     for (cell = cStart; cell < cEnd; ++cell)
213:       if (cellComp[cell - cStart] < 0) break;
214:     if (cell >= cEnd) break;
215:     /* Initialize FIFO with first cell in component */
216:     {
217:       const PetscInt *cone;
218:       PetscInt        coneSize;

220:       fTop = fBottom = 0;
221:       DMPlexGetConeSize(dm, cell, &coneSize);
222:       DMPlexGetCone(dm, cell, &cone);
223:       for (c = 0; c < coneSize; ++c) {
224:         faceFIFO[fBottom++] = cone[c];
225:         PetscBTSet(seenFaces, cone[c] - fStart);
226:       }
227:       PetscBTSet(seenCells, cell - cStart);
228:     }
229:     /* Consider each face in FIFO */
230:     while (fTop < fBottom) DMPlexCheckFace_Internal(dm, faceFIFO, &fTop, &fBottom, cStart, fStart, fEnd, seenCells, flippedCells, seenFaces);
231:     /* Set component for cells and faces */
232:     for (cell = 0; cell < cEnd - cStart; ++cell) {
233:       if (PetscBTLookup(seenCells, cell)) cellComp[cell] = comp;
234:     }
235:     for (face = 0; face < fEnd - fStart; ++face) {
236:       if (PetscBTLookup(seenFaces, face)) faceComp[face] = comp;
237:     }
238:     /* Wipe seenCells and seenFaces for next component */
239:     PetscBTMemzero(fEnd - fStart, seenFaces);
240:     PetscBTMemzero(cEnd - cStart, seenCells);
241:     ++comp;
242:   } while (1);
243:   numComponents = comp;
244:   if (flg) {
245:     PetscViewer v;

247:     PetscViewerASCIIGetStdout(comm, &v);
248:     PetscViewerASCIIPushSynchronized(v);
249:     PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank);
250:     PetscBTView(cEnd - cStart, flippedCells, v);
251:     PetscViewerFlush(v);
252:     PetscViewerASCIIPopSynchronized(v);
253:   }
254:   /* Now all subdomains are oriented, but we need a consistent parallel orientation */
255:   if (numLeaves >= 0) {
256:     PetscInt maxSupportSize, neighbor;

258:     /* Store orientations of boundary faces*/
259:     DMPlexGetMaxSizes(dm, NULL, &maxSupportSize);
260:     PetscCalloc3(numRoots, &rorntComp, numRoots, &lorntComp, maxSupportSize, &locSupport);
261:     for (face = fStart; face < fEnd; ++face) {
262:       const PetscInt *cone, *support, *ornt;
263:       PetscInt        coneSize, supportSize, Ns = 0, s, l;

265:       DMPlexGetSupportSize(dm, face, &supportSize);
266:       /* Ignore overlapping cells */
267:       DMPlexGetSupport(dm, face, &support);
268:       for (s = 0; s < supportSize; ++s) {
269:         PetscFindInt(support[s], numLeaves, lpoints, &l);
270:         if (l >= 0) continue;
271:         locSupport[Ns++] = support[s];
272:       }
273:       if (Ns != 1) continue;
274:       neighbor = locSupport[0];
275:       DMPlexGetCone(dm, neighbor, &cone);
276:       DMPlexGetConeSize(dm, neighbor, &coneSize);
277:       DMPlexGetConeOrientation(dm, neighbor, &ornt);
278:       for (c = 0; c < coneSize; ++c)
279:         if (cone[c] == face) break;
280:       if (dim == 1) {
281:         /* Use cone position instead, shifted to -1 or 1 */
282:         if (PetscBTLookup(flippedCells, neighbor - cStart)) rorntComp[face].rank = 1 - c * 2;
283:         else rorntComp[face].rank = c * 2 - 1;
284:       } else {
285:         if (PetscBTLookup(flippedCells, neighbor - cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 : 1;
286:         else rorntComp[face].rank = ornt[c] < 0 ? 1 : -1;
287:       }
288:       rorntComp[face].index = faceComp[face - fStart];
289:     }
290:     /* Communicate boundary edge orientations */
291:     PetscSFBcastBegin(sf, MPIU_2INT, rorntComp, lorntComp, MPI_REPLACE);
292:     PetscSFBcastEnd(sf, MPIU_2INT, rorntComp, lorntComp, MPI_REPLACE);
293:   }
294:   /* Get process adjacency */
295:   PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors);
296:   viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm));
297:   if (flg2) PetscViewerASCIIPushSynchronized(viewer);
298:   PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &selfviewer);
299:   for (comp = 0; comp < numComponents; ++comp) {
300:     PetscInt l, n;

302:     numNeighbors[comp] = 0;
303:     PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp]);
304:     /* I know this is p^2 time in general, but for bounded degree its alright */
305:     for (l = 0; l < numLeaves; ++l) {
306:       const PetscInt face = lpoints[l];

308:       /* Find a representative face (edge) separating pairs of procs */
309:       if ((face >= fStart) && (face < fEnd) && (faceComp[face - fStart] == comp) && rorntComp[face].rank) {
310:         const PetscInt rrank = rpoints[l].rank;
311:         const PetscInt rcomp = lorntComp[face].index;

313:         for (n = 0; n < numNeighbors[comp]; ++n)
314:           if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break;
315:         if (n >= numNeighbors[comp]) {
316:           PetscInt supportSize;

318:           DMPlexGetSupportSize(dm, face, &supportSize);
320:           if (flg)
321:             PetscCall(PetscViewerASCIIPrintf(selfviewer, "[%d]: component %d, Found representative leaf %" PetscInt_FMT " (face %" PetscInt_FMT ") connecting to face %" PetscInt_FMT " on (%" PetscInt_FMT ", %" PetscInt_FMT ") with orientation %" PetscInt_FMT "\n", rank, comp, l, face,
322:                                              rpoints[l].index, rrank, rcomp, lorntComp[face].rank));
323:           neighbors[comp][numNeighbors[comp]++] = l;
324:         }
325:       }
326:     }
327:     totNeighbors += numNeighbors[comp];
328:   }
329:   PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &selfviewer);
330:   PetscViewerFlush(viewer);
331:   if (flg2) PetscViewerASCIIPopSynchronized(viewer);
332:   PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match);
333:   for (comp = 0, off = 0; comp < numComponents; ++comp) {
334:     PetscInt n;

336:     for (n = 0; n < numNeighbors[comp]; ++n, ++off) {
337:       const PetscInt face = lpoints[neighbors[comp][n]];
338:       const PetscInt o    = rorntComp[face].rank * lorntComp[face].rank;

340:       if (o < 0) match[off] = PETSC_TRUE;
341:       else if (o > 0) match[off] = PETSC_FALSE;
342:       else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid face %" PetscInt_FMT " (%" PetscInt_FMT ", %" PetscInt_FMT ") neighbor: %" PetscInt_FMT " comp: %d", face, rorntComp[face].rank, lorntComp[face].rank, neighbors[comp][n], comp);
343:       nrankComp[off].rank  = rpoints[neighbors[comp][n]].rank;
344:       nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index;
345:     }
346:     PetscFree(neighbors[comp]);
347:   }
348:   /* Collect the graph on 0 */
349:   if (numLeaves >= 0) {
350:     Mat          G;
351:     PetscBT      seenProcs, flippedProcs;
352:     PetscInt    *procFIFO, pTop, pBottom;
353:     PetscInt    *N          = NULL, *Noff;
354:     PetscSFNode *adj        = NULL;
355:     PetscBool   *val        = NULL;
356:     PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc, p, o;
357:     PetscMPIInt  size = 0;

359:     PetscCalloc1(numComponents, &flipped);
360:     if (rank == 0) MPI_Comm_size(comm, &size);
361:     PetscCalloc4(size, &recvcounts, size + 1, &displs, size, &Nc, size + 1, &Noff);
362:     MPI_Gather(&numComponents, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm);
363:     for (p = 0; p < size; ++p) displs[p + 1] = displs[p] + Nc[p];
364:     if (rank == 0) PetscMalloc1(displs[size], &N);
365:     MPI_Gatherv(numNeighbors, numComponents, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm);
366:     for (p = 0, o = 0; p < size; ++p) {
367:       recvcounts[p] = 0;
368:       for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o];
369:       displs[p + 1] = displs[p] + recvcounts[p];
370:     }
371:     if (rank == 0) PetscMalloc2(displs[size], &adj, displs[size], &val);
372:     MPI_Gatherv(nrankComp, totNeighbors, MPIU_2INT, adj, recvcounts, displs, MPIU_2INT, 0, comm);
373:     MPI_Gatherv(match, totNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm);
374:     PetscFree2(numNeighbors, neighbors);
375:     if (rank == 0) {
376:       for (p = 1; p <= size; ++p) Noff[p] = Noff[p - 1] + Nc[p - 1];
377:       if (flg) {
378:         PetscInt n;

380:         for (p = 0, off = 0; p < size; ++p) {
381:           for (c = 0; c < Nc[p]; ++c) {
382:             PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %" PetscInt_FMT ":\n", p, c);
383:             for (n = 0; n < N[Noff[p] + c]; ++n, ++off) PetscPrintf(PETSC_COMM_SELF, "  edge (%" PetscInt_FMT ", %" PetscInt_FMT ") (%s):\n", adj[off].rank, adj[off].index, PetscBools[val[off]]);
384:           }
385:         }
386:       }
387:       /* Symmetrize the graph */
388:       MatCreate(PETSC_COMM_SELF, &G);
389:       MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size]);
390:       MatSetUp(G);
391:       for (p = 0, off = 0; p < size; ++p) {
392:         for (c = 0; c < Nc[p]; ++c) {
393:           const PetscInt r = Noff[p] + c;
394:           PetscInt       n;

396:           for (n = 0; n < N[r]; ++n, ++off) {
397:             const PetscInt    q = Noff[adj[off].rank] + adj[off].index;
398:             const PetscScalar o = val[off] ? 1.0 : 0.0;

400:             MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES);
401:             MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES);
402:           }
403:         }
404:       }
405:       MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);
406:       MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);

408:       PetscBTCreate(Noff[size], &seenProcs);
409:       PetscBTMemzero(Noff[size], seenProcs);
410:       PetscBTCreate(Noff[size], &flippedProcs);
411:       PetscBTMemzero(Noff[size], flippedProcs);
412:       PetscMalloc1(Noff[size], &procFIFO);
413:       pTop = pBottom = 0;
414:       for (p = 0; p < Noff[size]; ++p) {
415:         if (PetscBTLookup(seenProcs, p)) continue;
416:         /* Initialize FIFO with next proc */
417:         procFIFO[pBottom++] = p;
418:         PetscBTSet(seenProcs, p);
419:         /* Consider each proc in FIFO */
420:         while (pTop < pBottom) {
421:           const PetscScalar *ornt;
422:           const PetscInt    *neighbors;
423:           PetscInt           proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n;

425:           proc     = procFIFO[pTop++];
426:           flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0;
427:           MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt);
428:           /* Loop over neighboring procs */
429:           for (n = 0; n < numNeighbors; ++n) {
430:             nproc    = neighbors[n];
431:             mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1;
432:             seen     = PetscBTLookup(seenProcs, nproc);
433:             flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0;

435:             if (mismatch ^ (flippedA ^ flippedB)) {
437:               if (!flippedB) {
438:                 PetscBTSet(flippedProcs, nproc);
439:               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
441:             if (!seen) {
442:               procFIFO[pBottom++] = nproc;
443:               PetscBTSet(seenProcs, nproc);
444:             }
445:           }
446:         }
447:       }
448:       PetscFree(procFIFO);
449:       MatDestroy(&G);
450:       PetscFree2(adj, val);
451:       PetscBTDestroy(&seenProcs);
452:     }
453:     /* Scatter flip flags */
454:     {
455:       PetscBool *flips = NULL;

457:       if (rank == 0) {
458:         PetscMalloc1(Noff[size], &flips);
459:         for (p = 0; p < Noff[size]; ++p) {
460:           flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE;
461:           if (flg && flips[p]) PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p);
462:         }
463:         for (p = 0; p < size; ++p) displs[p + 1] = displs[p] + Nc[p];
464:       }
465:       MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, numComponents, MPIU_BOOL, 0, comm);
466:       PetscFree(flips);
467:     }
468:     if (rank == 0) PetscBTDestroy(&flippedProcs);
469:     PetscFree(N);
470:     PetscFree4(recvcounts, displs, Nc, Noff);
471:     PetscFree2(nrankComp, match);

473:     /* Decide whether to flip cells in each component */
474:     for (c = 0; c < cEnd - cStart; ++c) {
475:       if (flipped[cellComp[c]]) PetscBTNegate(flippedCells, c);
476:     }
477:     PetscFree(flipped);
478:   }
479:   if (flg) {
480:     PetscViewer v;

482:     PetscViewerASCIIGetStdout(comm, &v);
483:     PetscViewerASCIIPushSynchronized(v);
484:     PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank);
485:     PetscBTView(cEnd - cStart, flippedCells, v);
486:     PetscViewerFlush(v);
487:     PetscViewerASCIIPopSynchronized(v);
488:   }
489:   /* Reverse flipped cells in the mesh */
490:   for (c = cStart; c < cEnd; ++c) {
491:     if (PetscBTLookup(flippedCells, c - cStart)) DMPlexOrientPoint(dm, c, -1);
492:   }
493:   PetscBTDestroy(&seenCells);
494:   PetscBTDestroy(&flippedCells);
495:   PetscBTDestroy(&seenFaces);
496:   PetscFree2(numNeighbors, neighbors);
497:   PetscFree3(rorntComp, lorntComp, locSupport);
498:   PetscFree3(faceFIFO, cellComp, faceComp);
499:   return 0;
500: }