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: [](ch_unstructured), `DM`, `DMPLEX`, `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;

 25:   PetscFunctionBegin;
 27:   PetscCall(DMPlexGetCellType(dm, p, &ct));
 28:   arr = DMPolytopeTypeGetArrangement(ct, o);
 29:   if (!arr) PetscFunctionReturn(PETSC_SUCCESS);
 30:   PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
 31:   PetscCall(DMPlexGetCone(dm, p, &cone));
 32:   PetscCall(DMPlexGetConeOrientation(dm, p, &ornt));
 33:   PetscCall(DMGetWorkArray(dm, coneSize, MPIU_INT, &newcone));
 34:   PetscCall(DMGetWorkArray(dm, coneSize, MPIU_INT, &newornt));
 35:   for (c = 0; c < coneSize; ++c) {
 36:     DMPolytopeType ft;
 37:     PetscInt       nO;

 39:     PetscCall(DMPlexGetCellType(dm, cone[c], &ft));
 40:     nO         = DMPolytopeTypeGetNumArrangements(ft) / 2;
 41:     newcone[c] = cone[arr[c * 2 + 0]];
 42:     newornt[c] = DMPolytopeTypeComposeOrientation(ft, arr[c * 2 + 1], ornt[arr[c * 2 + 0]]);
 43:     PetscCheck(!newornt[c] || !(newornt[c] >= nO || newornt[c] < -nO), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid orientation %" PetscInt_FMT " not in [%" PetscInt_FMT ",%" PetscInt_FMT ") for %s %" PetscInt_FMT, newornt[c], -nO, nO, DMPolytopeTypes[ft], cone[c]);
 44:   }
 45:   PetscCall(DMPlexSetCone(dm, p, newcone));
 46:   PetscCall(DMPlexSetConeOrientation(dm, p, newornt));
 47:   PetscCall(DMRestoreWorkArray(dm, coneSize, MPIU_INT, &newcone));
 48:   PetscCall(DMRestoreWorkArray(dm, coneSize, MPIU_INT, &newornt));
 49:   /* Update orientation of this point in the support points */
 50:   PetscCall(DMPlexGetSupportSize(dm, p, &supportSize));
 51:   PetscCall(DMPlexGetSupport(dm, p, &support));
 52:   for (s = 0; s < supportSize; ++s) {
 53:     PetscCall(DMPlexGetConeSize(dm, support[s], &coneSize));
 54:     PetscCall(DMPlexGetCone(dm, support[s], &cone));
 55:     PetscCall(DMPlexGetConeOrientation(dm, support[s], &ornt));
 56:     for (c = 0; c < coneSize; ++c) {
 57:       PetscInt po;

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

 68: static PetscInt GetPointIndex(PetscInt point, PetscInt pStart, PetscInt pEnd, const PetscInt points[])
 69: {
 70:   if (points) {
 71:     PetscInt loc;

 73:     PetscCallAbort(PETSC_COMM_SELF, PetscFindInt(point, pEnd - pStart, points, &loc));
 74:     if (loc >= 0) return loc;
 75:   } else {
 76:     if (point >= pStart && point < pEnd) return point - pStart;
 77:   }
 78:   return -1;
 79: }

 81: /*
 82:   - Checks face match
 83:     - Flips non-matching
 84:   - Inserts faces of support cells in FIFO
 85: */
 86: static PetscErrorCode DMPlexCheckFace_Internal(DM dm, PetscInt *faceFIFO, PetscInt *fTop, PetscInt *fBottom, IS cellIS, IS faceIS, PetscBT seenCells, PetscBT flippedCells, PetscBT seenFaces)
 87: {
 88:   const PetscInt *supp, *coneA, *coneB, *coneOA, *coneOB;
 89:   PetscInt        suppSize, Ns = 0, coneSizeA, coneSizeB, posA = -1, posB = -1;
 90:   PetscInt        face, dim, indC[3], indS[3], seenA, flippedA, seenB, flippedB, mismatch;
 91:   const PetscInt *cells, *faces;
 92:   PetscInt        cStart, cEnd, fStart, fEnd;

 94:   PetscFunctionBegin;
 95:   face = faceFIFO[(*fTop)++];
 96:   PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
 97:   PetscCall(ISGetPointRange(faceIS, &fStart, &fEnd, &faces));
 98:   PetscCall(DMPlexGetPointDepth(dm, cells ? cells[cStart] : cStart, &dim));
 99:   PetscCall(DMPlexGetSupportSize(dm, face, &suppSize));
100:   PetscCall(DMPlexGetSupport(dm, face, &supp));
101:   // Filter the support
102:   for (PetscInt s = 0; s < suppSize; ++s) {
103:     // Filter support
104:     indC[Ns] = GetPointIndex(supp[s], cStart, cEnd, cells);
105:     indS[Ns] = s;
106:     if (indC[Ns] >= 0) ++Ns;
107:   }
108:   if (Ns < 2) PetscFunctionReturn(PETSC_SUCCESS);
109:   PetscCheck(Ns == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %" PetscInt_FMT, Ns);
110:   PetscCheck(indC[0] >= 0 && indC[1] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Support cells %" PetscInt_FMT " (%" PetscInt_FMT ") and %" PetscInt_FMT " (%" PetscInt_FMT ") are not both valid", supp[0], indC[0], supp[1], indC[1]);
111:   seenA    = PetscBTLookup(seenCells, indC[0]);
112:   flippedA = PetscBTLookup(flippedCells, indC[0]) ? 1 : 0;
113:   seenB    = PetscBTLookup(seenCells, indC[1]);
114:   flippedB = PetscBTLookup(flippedCells, indC[1]) ? 1 : 0;

116:   PetscCall(DMPlexGetConeSize(dm, supp[indS[0]], &coneSizeA));
117:   PetscCall(DMPlexGetConeSize(dm, supp[indS[1]], &coneSizeB));
118:   PetscCall(DMPlexGetCone(dm, supp[indS[0]], &coneA));
119:   PetscCall(DMPlexGetCone(dm, supp[indS[1]], &coneB));
120:   PetscCall(DMPlexGetConeOrientation(dm, supp[indS[0]], &coneOA));
121:   PetscCall(DMPlexGetConeOrientation(dm, supp[indS[1]], &coneOB));
122:   for (PetscInt c = 0; c < coneSizeA; ++c) {
123:     const PetscInt indF = GetPointIndex(coneA[c], fStart, fEnd, faces);

125:     // Filter cone
126:     if (indF < 0) continue;
127:     if (!PetscBTLookup(seenFaces, indF)) {
128:       faceFIFO[(*fBottom)++] = coneA[c];
129:       PetscCall(PetscBTSet(seenFaces, indF));
130:     }
131:     if (coneA[c] == face) posA = c;
132:     PetscCheck(*fBottom <= fEnd - fStart, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %" PetscInt_FMT " was pushed exceeding capacity %" PetscInt_FMT " > %" PetscInt_FMT, coneA[c], *fBottom, fEnd - fStart);
133:   }
134:   PetscCheck(posA >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " could not be located in cell %" PetscInt_FMT, face, supp[indS[0]]);
135:   for (PetscInt c = 0; c < coneSizeB; ++c) {
136:     const PetscInt indF = GetPointIndex(coneB[c], fStart, fEnd, faces);

138:     // Filter cone
139:     if (indF < 0) continue;
140:     if (!PetscBTLookup(seenFaces, indF)) {
141:       faceFIFO[(*fBottom)++] = coneB[c];
142:       PetscCall(PetscBTSet(seenFaces, indF));
143:     }
144:     if (coneB[c] == face) posB = c;
145:     PetscCheck(*fBottom <= fEnd - fStart, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %" PetscInt_FMT " was pushed exceeding capacity %" PetscInt_FMT " > %" PetscInt_FMT, coneA[c], *fBottom, fEnd - fStart);
146:   }
147:   PetscCheck(posB >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " could not be located in cell %" PetscInt_FMT, face, supp[indS[1]]);

149:   if (dim == 1) {
150:     mismatch = posA == posB;
151:   } else {
152:     mismatch = coneOA[posA] == coneOB[posB];
153:   }

155:   if (mismatch ^ (flippedA ^ flippedB)) {
156:     PetscCheck(!seenA || !seenB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Previously seen cells %" PetscInt_FMT " and %" PetscInt_FMT " do not match: Fault mesh is non-orientable", supp[indS[0]], supp[indS[1]]);
157:     if (!seenA && !flippedA) {
158:       PetscCall(PetscBTSet(flippedCells, indC[0]));
159:     } else if (!seenB && !flippedB) {
160:       PetscCall(PetscBTSet(flippedCells, indC[1]));
161:     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
162:   } else PetscCheck(!mismatch || !flippedA || !flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
163:   PetscCall(PetscBTSet(seenCells, indC[0]));
164:   PetscCall(PetscBTSet(seenCells, indC[1]));
165:   PetscFunctionReturn(PETSC_SUCCESS);
166: }

168: static PetscErrorCode DMPlexCheckFace_Old_Internal(DM dm, PetscInt *faceFIFO, PetscInt *fTop, PetscInt *fBottom, PetscInt cStart, PetscInt fStart, PetscInt fEnd, PetscBT seenCells, PetscBT flippedCells, PetscBT seenFaces)
169: {
170:   const PetscInt *support, *coneA, *coneB, *coneOA, *coneOB;
171:   PetscInt        supportSize, coneSizeA, coneSizeB, posA = -1, posB = -1;
172:   PetscInt        face, dim, seenA, flippedA, seenB, flippedB, mismatch, c;

174:   PetscFunctionBegin;
175:   face = faceFIFO[(*fTop)++];
176:   PetscCall(DMGetDimension(dm, &dim));
177:   PetscCall(DMPlexGetSupportSize(dm, face, &supportSize));
178:   PetscCall(DMPlexGetSupport(dm, face, &support));
179:   if (supportSize < 2) PetscFunctionReturn(PETSC_SUCCESS);
180:   PetscCheck(supportSize == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %" PetscInt_FMT, supportSize);
181:   seenA    = PetscBTLookup(seenCells, support[0] - cStart);
182:   flippedA = PetscBTLookup(flippedCells, support[0] - cStart) ? 1 : 0;
183:   seenB    = PetscBTLookup(seenCells, support[1] - cStart);
184:   flippedB = PetscBTLookup(flippedCells, support[1] - cStart) ? 1 : 0;

186:   PetscCall(DMPlexGetConeSize(dm, support[0], &coneSizeA));
187:   PetscCall(DMPlexGetConeSize(dm, support[1], &coneSizeB));
188:   PetscCall(DMPlexGetCone(dm, support[0], &coneA));
189:   PetscCall(DMPlexGetCone(dm, support[1], &coneB));
190:   PetscCall(DMPlexGetConeOrientation(dm, support[0], &coneOA));
191:   PetscCall(DMPlexGetConeOrientation(dm, support[1], &coneOB));
192:   for (c = 0; c < coneSizeA; ++c) {
193:     if (!PetscBTLookup(seenFaces, coneA[c] - fStart)) {
194:       faceFIFO[(*fBottom)++] = coneA[c];
195:       PetscCall(PetscBTSet(seenFaces, coneA[c] - fStart));
196:     }
197:     if (coneA[c] == face) posA = c;
198:     PetscCheck(*fBottom <= fEnd - fStart, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %" PetscInt_FMT " was pushed exceeding capacity %" PetscInt_FMT " > %" PetscInt_FMT, coneA[c], *fBottom, fEnd - fStart);
199:   }
200:   PetscCheck(posA >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " could not be located in cell %" PetscInt_FMT, face, support[0]);
201:   for (c = 0; c < coneSizeB; ++c) {
202:     if (!PetscBTLookup(seenFaces, coneB[c] - fStart)) {
203:       faceFIFO[(*fBottom)++] = coneB[c];
204:       PetscCall(PetscBTSet(seenFaces, coneB[c] - fStart));
205:     }
206:     if (coneB[c] == face) posB = c;
207:     PetscCheck(*fBottom <= fEnd - fStart, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %" PetscInt_FMT " was pushed exceeding capacity %" PetscInt_FMT " > %" PetscInt_FMT, coneA[c], *fBottom, fEnd - fStart);
208:   }
209:   PetscCheck(posB >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " could not be located in cell %" PetscInt_FMT, face, support[1]);

211:   if (dim == 1) {
212:     mismatch = posA == posB;
213:   } else {
214:     mismatch = coneOA[posA] == coneOB[posB];
215:   }

217:   if (mismatch ^ (flippedA ^ flippedB)) {
218:     PetscCheck(!seenA || !seenB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Previously seen cells %" PetscInt_FMT " and %" PetscInt_FMT " do not match: Fault mesh is non-orientable", support[0], support[1]);
219:     if (!seenA && !flippedA) {
220:       PetscCall(PetscBTSet(flippedCells, support[0] - cStart));
221:     } else if (!seenB && !flippedB) {
222:       PetscCall(PetscBTSet(flippedCells, support[1] - cStart));
223:     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
224:   } else PetscCheck(!mismatch || !flippedA || !flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
225:   PetscCall(PetscBTSet(seenCells, support[0] - cStart));
226:   PetscCall(PetscBTSet(seenCells, support[1] - cStart));
227:   PetscFunctionReturn(PETSC_SUCCESS);
228: }

230: /*
231:   DMPlexOrient_Serial - Compute valid orientation for local connected components

233:   Not collective

235:   Input Parameters:
236:   + dm         - The `DM`
237:   - cellHeight - The height of k-cells to be oriented

239:   Output Parameters:
240:   + Ncomp        - The number of connected component
241:   . cellComp     - The connected component for each local cell
242:   . faceComp     - The connected component for each local face
243:   - flippedCells - Marked cells should be inverted

245:   Level: developer

247: .seealso: `DMPlexOrient()`
248: */
249: static PetscErrorCode DMPlexOrient_Serial(DM dm, IS cellIS, IS faceIS, PetscInt *Ncomp, PetscInt cellComp[], PetscInt faceComp[], PetscBT flippedCells)
250: {
251:   PetscBT         seenCells, seenFaces;
252:   PetscInt       *faceFIFO;
253:   const PetscInt *cells = NULL, *faces = NULL;
254:   PetscInt        cStart = 0, cEnd = 0, fStart = 0, fEnd = 0;

256:   PetscFunctionBegin;
257:   if (cellIS) PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
258:   if (faceIS) PetscCall(ISGetPointRange(faceIS, &fStart, &fEnd, &faces));
259:   PetscCall(PetscBTCreate(cEnd - cStart, &seenCells));
260:   PetscCall(PetscBTMemzero(cEnd - cStart, seenCells));
261:   PetscCall(PetscBTCreate(fEnd - fStart, &seenFaces));
262:   PetscCall(PetscBTMemzero(fEnd - fStart, seenFaces));
263:   PetscCall(PetscMalloc1(fEnd - fStart, &faceFIFO));
264:   *Ncomp = 0;
265:   for (PetscInt c = 0; c < cEnd - cStart; ++c) cellComp[c] = -1;
266:   do {
267:     PetscInt cc, fTop, fBottom;

269:     // Look for first unmarked cell
270:     for (cc = cStart; cc < cEnd; ++cc)
271:       if (cellComp[cc - cStart] < 0) break;
272:     if (cc >= cEnd) break;
273:     // Initialize FIFO with first cell in component
274:     {
275:       const PetscInt  cell = cells ? cells[cc] : cc;
276:       const PetscInt *cone;
277:       PetscInt        coneSize;

279:       fTop = fBottom = 0;
280:       PetscCall(DMPlexGetConeSize(dm, cell, &coneSize));
281:       PetscCall(DMPlexGetCone(dm, cell, &cone));
282:       for (PetscInt c = 0; c < coneSize; ++c) {
283:         const PetscInt idx = GetPointIndex(cone[c], fStart, fEnd, faces);

285:         // Cell faces are guaranteed to be in the face set
286:         PetscCheck(idx >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " of cell %" PetscInt_FMT " is not present in the label", cone[c], cell);
287:         faceFIFO[fBottom++] = cone[c];
288:         PetscCall(PetscBTSet(seenFaces, idx));
289:       }
290:       PetscCall(PetscBTSet(seenCells, cc - cStart));
291:     }
292:     // Consider each face in FIFO
293:     while (fTop < fBottom) PetscCall(DMPlexCheckFace_Internal(dm, faceFIFO, &fTop, &fBottom, cellIS, faceIS, seenCells, flippedCells, seenFaces));
294:     // Set component for cells and faces
295:     for (PetscInt c = 0; c < cEnd - cStart; ++c) {
296:       if (PetscBTLookup(seenCells, c)) cellComp[c] = *Ncomp;
297:     }
298:     for (PetscInt f = 0; f < fEnd - fStart; ++f) {
299:       if (PetscBTLookup(seenFaces, f)) faceComp[f] = *Ncomp;
300:     }
301:     // Wipe seenCells and seenFaces for next component
302:     PetscCall(PetscBTMemzero(fEnd - fStart, seenFaces));
303:     PetscCall(PetscBTMemzero(cEnd - cStart, seenCells));
304:     ++(*Ncomp);
305:   } while (1);
306:   PetscCall(PetscBTDestroy(&seenCells));
307:   PetscCall(PetscBTDestroy(&seenFaces));
308:   PetscCall(PetscFree(faceFIFO));
309:   PetscFunctionReturn(PETSC_SUCCESS);
310: }

312: /*@
313:   DMPlexOrient - Give a consistent orientation to the input mesh

315:   Input Parameter:
316: . dm - The `DM`

318:   Note:
319:   The orientation data for the `DM` are change in-place.

321:   This routine will fail for non-orientable surfaces, such as the Moebius strip.

323:   Level: advanced

325: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMCreate()`
326: @*/
327: PetscErrorCode DMPlexOrient(DM dm)
328: {
329: #if 0
330:   IS cellIS, faceIS;

332:   PetscFunctionBegin;
333:   PetscCall(DMPlexGetAllCells_Internal(dm, &cellIS));
334:   PetscCall(DMPlexGetAllFaces_Internal(dm, &faceIS));
335:   PetscCall(DMPlexOrientCells_Internal(dm, cellIS, faceIS));
336:   PetscCall(ISDestroy(&cellIS));
337:   PetscCall(ISDestroy(&faceIS));
338:   PetscFunctionReturn(PETSC_SUCCESS);
339: #else
340:   MPI_Comm           comm;
341:   PetscSF            sf;
342:   const PetscInt    *lpoints;
343:   const PetscSFNode *rpoints;
344:   PetscSFNode       *rorntComp = NULL, *lorntComp = NULL;
345:   PetscInt          *numNeighbors, **neighbors, *locSupport = NULL;
346:   PetscSFNode       *nrankComp;
347:   PetscBool         *match, *flipped;
348:   PetscBT            seenCells, flippedCells, seenFaces;
349:   PetscInt          *faceFIFO, fTop, fBottom, *cellComp, *faceComp;
350:   PetscInt           numLeaves, numRoots, dim, h, cStart, cEnd, c, cell, fStart, fEnd, face, off, totNeighbors = 0;
351:   PetscMPIInt        rank, size, numComponents, comp = 0;
352:   PetscBool          flg, flg2;
353:   PetscViewer        viewer = NULL, selfviewer = NULL;

355:   PetscFunctionBegin;
356:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
357:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
358:   PetscCallMPI(MPI_Comm_size(comm, &size));
359:   PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-orientation_view", &flg));
360:   PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-orientation_view_synchronized", &flg2));
361:   PetscCall(DMGetPointSF(dm, &sf));
362:   PetscCall(PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints));
363:   /* Truth Table
364:      mismatch    flips   do action   mismatch   flipA ^ flipB   action
365:          F       0 flips     no         F             F           F
366:          F       1 flip      yes        F             T           T
367:          F       2 flips     no         T             F           T
368:          T       0 flips     yes        T             T           F
369:          T       1 flip      no
370:          T       2 flips     yes
371:   */
372:   PetscCall(DMGetDimension(dm, &dim));
373:   PetscCall(DMPlexGetVTKCellHeight(dm, &h));
374:   PetscCall(DMPlexGetHeightStratum(dm, h, &cStart, &cEnd));
375:   PetscCall(DMPlexGetHeightStratum(dm, h + 1, &fStart, &fEnd));
376:   PetscCall(PetscBTCreate(cEnd - cStart, &seenCells));
377:   PetscCall(PetscBTMemzero(cEnd - cStart, seenCells));
378:   PetscCall(PetscBTCreate(cEnd - cStart, &flippedCells));
379:   PetscCall(PetscBTMemzero(cEnd - cStart, flippedCells));
380:   PetscCall(PetscBTCreate(fEnd - fStart, &seenFaces));
381:   PetscCall(PetscBTMemzero(fEnd - fStart, seenFaces));
382:   PetscCall(PetscCalloc3(fEnd - fStart, &faceFIFO, cEnd - cStart, &cellComp, fEnd - fStart, &faceComp));
383:   /*
384:    OLD STYLE
385:    - Add an integer array over cells and faces (component) for connected component number
386:    Foreach component
387:      - Mark the initial cell as seen
388:      - Process component as usual
389:      - Set component for all seenCells
390:      - Wipe seenCells and seenFaces (flippedCells can stay)
391:    - Generate parallel adjacency for component using SF and seenFaces
392:    - Collect numComponents adj data from each proc to 0
393:    - Build same serial graph
394:    - Use same solver
395:    - Use Scatterv to send back flipped flags for each component
396:    - Negate flippedCells by component

398:    NEW STYLE
399:    - Create the adj on each process
400:    - Bootstrap to complete graph on proc 0
401:   */
402:   /* Loop over components */
403:   for (cell = cStart; cell < cEnd; ++cell) cellComp[cell - cStart] = -1;
404:   do {
405:     /* Look for first unmarked cell */
406:     for (cell = cStart; cell < cEnd; ++cell)
407:       if (cellComp[cell - cStart] < 0) break;
408:     if (cell >= cEnd) break;
409:     /* Initialize FIFO with first cell in component */
410:     {
411:       const PetscInt *cone;
412:       PetscInt        coneSize;

414:       fTop = fBottom = 0;
415:       PetscCall(DMPlexGetConeSize(dm, cell, &coneSize));
416:       PetscCall(DMPlexGetCone(dm, cell, &cone));
417:       for (c = 0; c < coneSize; ++c) {
418:         faceFIFO[fBottom++] = cone[c];
419:         PetscCall(PetscBTSet(seenFaces, cone[c] - fStart));
420:       }
421:       PetscCall(PetscBTSet(seenCells, cell - cStart));
422:     }
423:     /* Consider each face in FIFO */
424:     while (fTop < fBottom) PetscCall(DMPlexCheckFace_Old_Internal(dm, faceFIFO, &fTop, &fBottom, cStart, fStart, fEnd, seenCells, flippedCells, seenFaces));
425:     /* Set component for cells and faces */
426:     for (cell = 0; cell < cEnd - cStart; ++cell) {
427:       if (PetscBTLookup(seenCells, cell)) cellComp[cell] = comp;
428:     }
429:     for (face = 0; face < fEnd - fStart; ++face) {
430:       if (PetscBTLookup(seenFaces, face)) faceComp[face] = comp;
431:     }
432:     /* Wipe seenCells and seenFaces for next component */
433:     PetscCall(PetscBTMemzero(fEnd - fStart, seenFaces));
434:     PetscCall(PetscBTMemzero(cEnd - cStart, seenCells));
435:     ++comp;
436:   } while (1);
437:   numComponents = comp;
438:   if (flg) {
439:     PetscViewer v;

441:     PetscCall(PetscViewerASCIIGetStdout(comm, &v));
442:     PetscCall(PetscViewerASCIIPushSynchronized(v));
443:     PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank));
444:     PetscCall(PetscBTView(cEnd - cStart, flippedCells, v));
445:     PetscCall(PetscViewerFlush(v));
446:     PetscCall(PetscViewerASCIIPopSynchronized(v));
447:   }
448:   /* Now all subdomains are oriented, but we need a consistent parallel orientation */
449:   if (numLeaves >= 0) {
450:     PetscInt maxSupportSize, neighbor;

452:     /* Store orientations of boundary faces*/
453:     PetscCall(DMPlexGetMaxSizes(dm, NULL, &maxSupportSize));
454:     PetscCall(PetscCalloc3(numRoots, &rorntComp, numRoots, &lorntComp, maxSupportSize, &locSupport));
455:     for (face = fStart; face < fEnd; ++face) {
456:       const PetscInt *cone, *support, *ornt;
457:       PetscInt        coneSize, supportSize, Ns = 0, s, l;

459:       PetscCall(DMPlexGetSupportSize(dm, face, &supportSize));
460:       /* Ignore overlapping cells */
461:       PetscCall(DMPlexGetSupport(dm, face, &support));
462:       for (s = 0; s < supportSize; ++s) {
463:         PetscCall(PetscFindInt(support[s], numLeaves, lpoints, &l));
464:         if (l >= 0) continue;
465:         locSupport[Ns++] = support[s];
466:       }
467:       if (Ns != 1) continue;
468:       neighbor = locSupport[0];
469:       PetscCall(DMPlexGetCone(dm, neighbor, &cone));
470:       PetscCall(DMPlexGetConeSize(dm, neighbor, &coneSize));
471:       PetscCall(DMPlexGetConeOrientation(dm, neighbor, &ornt));
472:       for (c = 0; c < coneSize; ++c)
473:         if (cone[c] == face) break;
474:       if (dim == 1) {
475:         /* Use cone position instead, shifted to -1 or 1 */
476:         if (PetscBTLookup(flippedCells, neighbor - cStart)) rorntComp[face].rank = 1 - c * 2;
477:         else rorntComp[face].rank = c * 2 - 1;
478:       } else {
479:         if (PetscBTLookup(flippedCells, neighbor - cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 : 1;
480:         else rorntComp[face].rank = ornt[c] < 0 ? 1 : -1;
481:       }
482:       rorntComp[face].index = faceComp[face - fStart];
483:     }
484:     /* Communicate boundary edge orientations */
485:     PetscCall(PetscSFBcastBegin(sf, MPIU_SF_NODE, rorntComp, lorntComp, MPI_REPLACE));
486:     PetscCall(PetscSFBcastEnd(sf, MPIU_SF_NODE, rorntComp, lorntComp, MPI_REPLACE));
487:   }
488:   /* Get process adjacency */
489:   PetscCall(PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors));
490:   viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm));
491:   if (flg2) PetscCall(PetscViewerASCIIPushSynchronized(viewer));
492:   PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &selfviewer));
493:   for (comp = 0; comp < numComponents; ++comp) {
494:     PetscInt l, n;

496:     numNeighbors[comp] = 0;
497:     PetscCall(PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp]));
498:     /* I know this is p^2 time in general, but for bounded degree its alright */
499:     for (l = 0; l < numLeaves; ++l) {
500:       const PetscInt face = lpoints[l];

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

507:         for (n = 0; n < numNeighbors[comp]; ++n)
508:           if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break;
509:         if (n >= numNeighbors[comp]) {
510:           PetscInt supportSize;

512:           PetscCall(DMPlexGetSupportSize(dm, face, &supportSize));
513:           PetscCheck(supportSize == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %" PetscInt_FMT, supportSize);
514:           if (flg)
515:             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,
516:                                              rpoints[l].index, rrank, rcomp, lorntComp[face].rank));
517:           neighbors[comp][numNeighbors[comp]++] = l;
518:         }
519:       }
520:     }
521:     totNeighbors += numNeighbors[comp];
522:   }
523:   PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &selfviewer));
524:   if (flg2) PetscCall(PetscViewerASCIIPopSynchronized(viewer));
525:   PetscCall(PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match));
526:   for (comp = 0, off = 0; comp < numComponents; ++comp) {
527:     PetscInt n;

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

533:       if (o < 0) match[off] = PETSC_TRUE;
534:       else if (o > 0) match[off] = PETSC_FALSE;
535:       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);
536:       nrankComp[off].rank  = rpoints[neighbors[comp][n]].rank;
537:       nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index;
538:     }
539:     PetscCall(PetscFree(neighbors[comp]));
540:   }
541:   /* Collect the graph on 0 */
542:   if (numLeaves >= 0) {
543:     Mat          G;
544:     PetscBT      seenProcs, flippedProcs;
545:     PetscInt    *procFIFO, pTop, pBottom;
546:     PetscInt    *N          = NULL, *Noff;
547:     PetscSFNode *adj        = NULL;
548:     PetscBool   *val        = NULL;
549:     PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc, p, o, itotNeighbors;
550:     PetscMPIInt  size = 0;

552:     PetscCall(PetscCalloc1(numComponents, &flipped));
553:     if (rank == 0) PetscCallMPI(MPI_Comm_size(comm, &size));
554:     PetscCall(PetscCalloc4(size, &recvcounts, size + 1, &displs, size, &Nc, size + 1, &Noff));
555:     PetscCallMPI(MPI_Gather(&numComponents, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm));
556:     for (p = 0; p < size; ++p) displs[p + 1] = displs[p] + Nc[p];
557:     if (rank == 0) PetscCall(PetscMalloc1(displs[size], &N));
558:     PetscCallMPI(MPI_Gatherv(numNeighbors, numComponents, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm));
559:     for (p = 0, o = 0; p < size; ++p) {
560:       recvcounts[p] = 0;
561:       for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o];
562:       displs[p + 1] = displs[p] + recvcounts[p];
563:     }
564:     if (rank == 0) PetscCall(PetscMalloc2(displs[size], &adj, displs[size], &val));
565:     PetscCall(PetscMPIIntCast(totNeighbors, &itotNeighbors));
566:     PetscCallMPI(MPI_Gatherv(nrankComp, itotNeighbors, MPIU_SF_NODE, adj, recvcounts, displs, MPIU_SF_NODE, 0, comm));
567:     PetscCallMPI(MPI_Gatherv(match, itotNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm));
568:     PetscCall(PetscFree2(numNeighbors, neighbors));
569:     if (rank == 0) {
570:       for (p = 1; p <= size; ++p) Noff[p] = Noff[p - 1] + Nc[p - 1];
571:       if (flg) {
572:         PetscInt n;

574:         for (p = 0, off = 0; p < size; ++p) {
575:           for (c = 0; c < Nc[p]; ++c) {
576:             PetscCall(PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %" PetscInt_FMT ":\n", p, c));
577:             for (n = 0; n < N[Noff[p] + c]; ++n, ++off) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  edge (%" PetscInt_FMT ", %" PetscInt_FMT ") (%s):\n", adj[off].rank, adj[off].index, PetscBools[val[off]]));
578:           }
579:         }
580:       }
581:       /* Symmetrize the graph */
582:       PetscCall(MatCreate(PETSC_COMM_SELF, &G));
583:       PetscCall(MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size]));
584:       PetscCall(MatSetUp(G));
585:       for (p = 0, off = 0; p < size; ++p) {
586:         for (c = 0; c < Nc[p]; ++c) {
587:           const PetscInt r = Noff[p] + c;
588:           PetscInt       n;

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

594:             PetscCall(MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES));
595:             PetscCall(MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES));
596:           }
597:         }
598:       }
599:       PetscCall(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY));
600:       PetscCall(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY));

602:       PetscCall(PetscBTCreate(Noff[size], &seenProcs));
603:       PetscCall(PetscBTMemzero(Noff[size], seenProcs));
604:       PetscCall(PetscBTCreate(Noff[size], &flippedProcs));
605:       PetscCall(PetscBTMemzero(Noff[size], flippedProcs));
606:       PetscCall(PetscMalloc1(Noff[size], &procFIFO));
607:       pTop = pBottom = 0;
608:       for (p = 0; p < Noff[size]; ++p) {
609:         if (PetscBTLookup(seenProcs, p)) continue;
610:         /* Initialize FIFO with next proc */
611:         procFIFO[pBottom++] = p;
612:         PetscCall(PetscBTSet(seenProcs, p));
613:         /* Consider each proc in FIFO */
614:         while (pTop < pBottom) {
615:           const PetscScalar *ornt;
616:           const PetscInt    *neighbors;
617:           PetscInt           proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n;

619:           proc     = procFIFO[pTop++];
620:           flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0;
621:           PetscCall(MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt));
622:           /* Loop over neighboring procs */
623:           for (n = 0; n < numNeighbors; ++n) {
624:             nproc    = neighbors[n];
625:             mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1;
626:             seen     = PetscBTLookup(seenProcs, nproc);
627:             flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0;

629:             if (mismatch ^ (flippedA ^ flippedB)) {
630:               PetscCheck(!seen, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Previously seen procs %" PetscInt_FMT " and %" PetscInt_FMT " do not match: Fault mesh is non-orientable", proc, nproc);
631:               if (!flippedB) {
632:                 PetscCall(PetscBTSet(flippedProcs, nproc));
633:               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
634:             } else PetscCheck(!mismatch || !flippedA || !flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
635:             if (!seen) {
636:               procFIFO[pBottom++] = nproc;
637:               PetscCall(PetscBTSet(seenProcs, nproc));
638:             }
639:           }
640:         }
641:       }
642:       PetscCall(PetscFree(procFIFO));
643:       PetscCall(MatDestroy(&G));
644:       PetscCall(PetscFree2(adj, val));
645:       PetscCall(PetscBTDestroy(&seenProcs));
646:     }
647:     /* Scatter flip flags */
648:     {
649:       PetscBool *flips = NULL;

651:       if (rank == 0) {
652:         PetscCall(PetscMalloc1(Noff[size], &flips));
653:         for (p = 0; p < Noff[size]; ++p) {
654:           flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE;
655:           if (flg && flips[p]) PetscCall(PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p));
656:         }
657:         for (p = 0; p < size; ++p) displs[p + 1] = displs[p] + Nc[p];
658:       }
659:       PetscCallMPI(MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, numComponents, MPIU_BOOL, 0, comm));
660:       PetscCall(PetscFree(flips));
661:     }
662:     if (rank == 0) PetscCall(PetscBTDestroy(&flippedProcs));
663:     PetscCall(PetscFree(N));
664:     PetscCall(PetscFree4(recvcounts, displs, Nc, Noff));
665:     PetscCall(PetscFree2(nrankComp, match));

667:     /* Decide whether to flip cells in each component */
668:     for (c = 0; c < cEnd - cStart; ++c) {
669:       if (flipped[cellComp[c]]) PetscCall(PetscBTNegate(flippedCells, c));
670:     }
671:     PetscCall(PetscFree(flipped));
672:   }
673:   if (flg) {
674:     PetscViewer v;

676:     PetscCall(PetscViewerASCIIGetStdout(comm, &v));
677:     PetscCall(PetscViewerASCIIPushSynchronized(v));
678:     PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank));
679:     PetscCall(PetscBTView(cEnd - cStart, flippedCells, v));
680:     PetscCall(PetscViewerFlush(v));
681:     PetscCall(PetscViewerASCIIPopSynchronized(v));
682:   }
683:   /* Reverse flipped cells in the mesh */
684:   for (c = cStart; c < cEnd; ++c) {
685:     if (PetscBTLookup(flippedCells, c - cStart)) PetscCall(DMPlexOrientPoint(dm, c, -1));
686:   }
687:   PetscCall(PetscBTDestroy(&seenCells));
688:   PetscCall(PetscBTDestroy(&flippedCells));
689:   PetscCall(PetscBTDestroy(&seenFaces));
690:   PetscCall(PetscFree2(numNeighbors, neighbors));
691:   PetscCall(PetscFree3(rorntComp, lorntComp, locSupport));
692:   PetscCall(PetscFree3(faceFIFO, cellComp, faceComp));
693:   PetscFunctionReturn(PETSC_SUCCESS);
694: #endif
695: }

697: static PetscErrorCode CreateCellAndFaceIS_Private(DM dm, DMLabel label, IS *cellIS, IS *faceIS)
698: {
699:   IS              valueIS;
700:   const PetscInt *values;
701:   PetscInt        Nv, depth = 0;

703:   PetscFunctionBegin;
704:   PetscCall(DMLabelGetValueIS(label, &valueIS));
705:   PetscCall(ISGetLocalSize(valueIS, &Nv));
706:   PetscCall(ISGetIndices(valueIS, &values));
707:   for (PetscInt v = 0; v < Nv; ++v) {
708:     const PetscInt val = values[v] < 0 || values[v] >= 100 ? 0 : values[v];
709:     PetscInt       n;

711:     PetscCall(DMLabelGetStratumSize(label, val, &n));
712:     if (!n) continue;
713:     depth = PetscMax(val, depth);
714:   }
715:   PetscCall(ISDestroy(&valueIS));
716:   PetscCheck(depth >= 1 || !Nv, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Depth for interface must be at least 1, not %" PetscInt_FMT, depth);
717:   PetscCall(DMLabelGetStratumIS(label, depth, cellIS));
718:   PetscCall(DMLabelGetStratumIS(label, depth - 1, faceIS));
719:   if (!*cellIS) PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, cellIS));
720:   if (!*faceIS) PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, faceIS));
721:   PetscFunctionReturn(PETSC_SUCCESS);
722: }

724: PetscErrorCode DMPlexOrientLabel(DM dm, DMLabel label)
725: {
726:   IS cellIS, faceIS;

728:   PetscFunctionBegin;
729:   PetscCall(CreateCellAndFaceIS_Private(dm, label, &cellIS, &faceIS));
730:   PetscCall(DMPlexOrientCells_Internal(dm, cellIS, faceIS));
731:   PetscCall(ISDestroy(&cellIS));
732:   PetscCall(ISDestroy(&faceIS));
733:   PetscFunctionReturn(PETSC_SUCCESS);
734: }

736: PetscErrorCode DMPlexOrientCells_Internal(DM dm, IS cellIS, IS faceIS)
737: {
738:   MPI_Comm           comm;
739:   PetscSF            sf;
740:   const PetscInt    *lpoints;
741:   const PetscSFNode *rpoints;
742:   PetscSFNode       *rorntComp = NULL, *lorntComp = NULL;
743:   PetscInt          *numNeighbors, **neighbors, *locSupp = NULL;
744:   PetscSFNode       *nrankComp;
745:   PetscBool         *match, *flipped;
746:   PetscBT            flippedCells;
747:   PetscInt          *cellComp, *faceComp;
748:   const PetscInt    *cells = NULL, *faces = NULL;
749:   PetscInt           cStart = 0, cEnd = 0, fStart = 0, fEnd = 0;
750:   PetscInt           numLeaves, numRoots, dim, Ncomp, totNeighbors = 0;
751:   PetscMPIInt        rank, size;
752:   PetscBool          view, viewSync;
753:   PetscViewer        viewer = NULL, selfviewer = NULL;

755:   PetscFunctionBegin;
756:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
757:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
758:   PetscCallMPI(MPI_Comm_size(comm, &size));
759:   PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-orientation_view", &view));
760:   PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-orientation_view_synchronized", &viewSync));

762:   if (cellIS) PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
763:   if (faceIS) PetscCall(ISGetPointRange(faceIS, &fStart, &fEnd, &faces));
764:   PetscCall(DMGetPointSF(dm, &sf));
765:   PetscCall(PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints));
766:   /* Truth Table
767:      mismatch    flips   do action   mismatch   flipA ^ flipB   action
768:          F       0 flips     no         F             F           F
769:          F       1 flip      yes        F             T           T
770:          F       2 flips     no         T             F           T
771:          T       0 flips     yes        T             T           F
772:          T       1 flip      no
773:          T       2 flips     yes
774:   */
775:   PetscCall(DMGetDimension(dm, &dim));
776:   PetscCall(PetscBTCreate(cEnd - cStart, &flippedCells));
777:   PetscCall(PetscBTMemzero(cEnd - cStart, flippedCells));
778:   PetscCall(PetscCalloc2(cEnd - cStart, &cellComp, fEnd - fStart, &faceComp));
779:   /*
780:    OLD STYLE
781:    - Add an integer array over cells and faces (component) for connected component number
782:    Foreach component
783:      - Mark the initial cell as seen
784:      - Process component as usual
785:      - Set component for all seenCells
786:      - Wipe seenCells and seenFaces (flippedCells can stay)
787:    - Generate parallel adjacency for component using SF and seenFaces
788:    - Collect Ncomp adj data from each proc to 0
789:    - Build same serial graph
790:    - Use same solver
791:    - Use Scatterv to send back flipped flags for each component
792:    - Negate flippedCells by component

794:    NEW STYLE
795:    - Create the adj on each process
796:    - Bootstrap to complete graph on proc 0
797:   */
798:   PetscCall(DMPlexOrient_Serial(dm, cellIS, faceIS, &Ncomp, cellComp, faceComp, flippedCells));
799:   if (view) {
800:     PetscViewer v;
801:     PetscInt    cdepth = -1;

803:     PetscCall(PetscViewerASCIIGetStdout(comm, &v));
804:     PetscCall(PetscViewerASCIIPushSynchronized(v));
805:     if (cEnd > cStart) PetscCall(DMPlexGetPointDepth(dm, cells ? cells[cStart] : cStart, &cdepth));
806:     PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]New Orientation %" PetscInt_FMT " cells (depth %" PetscInt_FMT ") and %" PetscInt_FMT " faces\n", rank, cEnd - cStart, cdepth, fEnd - fStart));
807:     PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank));
808:     PetscCall(PetscBTView(cEnd - cStart, flippedCells, v));
809:     PetscCall(PetscViewerFlush(v));
810:     PetscCall(PetscViewerASCIIPopSynchronized(v));
811:   }
812:   /* Now all subdomains are oriented, but we need a consistent parallel orientation */
813:   // TODO: This all has to be rewritten to filter cones/supports to the ISes
814:   if (numLeaves >= 0) {
815:     PetscInt maxSuppSize, neighbor;

817:     // Store orientations of boundary faces
818:     PetscCall(DMPlexGetMaxSizes(dm, NULL, &maxSuppSize));
819:     PetscCall(PetscCalloc3(numRoots, &rorntComp, numRoots, &lorntComp, maxSuppSize, &locSupp));
820:     for (PetscInt f = fStart; f < fEnd; ++f) {
821:       const PetscInt  face = faces ? faces[f] : f;
822:       const PetscInt *cone, *supp, *ornt;
823:       PetscInt        coneSize, suppSize, nind, c, Ns = 0;

825:       PetscCall(DMPlexGetSupportSize(dm, face, &suppSize));
826:       PetscCall(DMPlexGetSupport(dm, face, &supp));
827:       for (PetscInt s = 0; s < suppSize; ++s) {
828:         PetscInt ind, l;

830:         // Filter support
831:         ind = GetPointIndex(supp[s], cStart, cEnd, cells);
832:         if (ind < 0) continue;
833:         // Ignore overlapping cells
834:         PetscCall(PetscFindInt(supp[s], numLeaves, lpoints, &l));
835:         if (l >= 0) continue;
836:         locSupp[Ns++] = supp[s];
837:       }
838:       PetscCheck(Ns < maxSuppSize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Index %" PetscInt_FMT " exceeds array size %" PetscInt_FMT, Ns, maxSuppSize);
839:       if (Ns != 1) continue;
840:       neighbor = locSupp[0];
841:       nind     = GetPointIndex(neighbor, cStart, cEnd, cells);
842:       PetscCall(DMPlexGetCone(dm, neighbor, &cone));
843:       PetscCall(DMPlexGetConeSize(dm, neighbor, &coneSize));
844:       PetscCall(DMPlexGetConeOrientation(dm, neighbor, &ornt));
845:       for (c = 0; c < coneSize; ++c)
846:         if (cone[c] == face) break;
847:       if (dim == 1) {
848:         /* Use cone position instead, shifted to -1 or 1 */
849:         if (PetscBTLookup(flippedCells, nind)) rorntComp[face].rank = 1 - c * 2;
850:         else rorntComp[face].rank = c * 2 - 1;
851:       } else {
852:         if (PetscBTLookup(flippedCells, nind)) rorntComp[face].rank = ornt[c] < 0 ? -1 : 1;
853:         else rorntComp[face].rank = ornt[c] < 0 ? 1 : -1;
854:       }
855:       rorntComp[face].index = faceComp[GetPointIndex(face, fStart, fEnd, faces)];
856:     }
857:     // Communicate boundary edge orientations
858:     PetscCall(PetscSFBcastBegin(sf, MPIU_SF_NODE, rorntComp, lorntComp, MPI_REPLACE));
859:     PetscCall(PetscSFBcastEnd(sf, MPIU_SF_NODE, rorntComp, lorntComp, MPI_REPLACE));
860:   }
861:   /* Get process adjacency */
862:   PetscCall(PetscMalloc2(Ncomp, &numNeighbors, Ncomp, &neighbors));
863:   viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm));
864:   if (viewSync) PetscCall(PetscViewerASCIIPushSynchronized(viewer));
865:   PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &selfviewer));
866:   for (PetscInt comp = 0; comp < Ncomp; ++comp) {
867:     PetscInt n;

869:     numNeighbors[comp] = 0;
870:     PetscCall(PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp]));
871:     /* I know this is p^2 time in general, but for bounded degree its alright */
872:     for (PetscInt l = 0; l < numLeaves; ++l) {
873:       const PetscInt face = lpoints[l];
874:       PetscInt       find;

876:       /* Find a representative face (edge) separating pairs of procs */
877:       find = GetPointIndex(face, fStart, fEnd, faces);
878:       if ((find >= 0) && (faceComp[find] == comp) && rorntComp[face].rank) {
879:         const PetscInt rrank = rpoints[l].rank;
880:         const PetscInt rcomp = lorntComp[face].index;

882:         for (n = 0; n < numNeighbors[comp]; ++n)
883:           if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break;
884:         if (n >= numNeighbors[comp]) {
885:           const PetscInt *supp;
886:           PetscInt        suppSize, Ns = 0;

888:           PetscCall(DMPlexGetSupport(dm, face, &supp));
889:           PetscCall(DMPlexGetSupportSize(dm, face, &suppSize));
890:           for (PetscInt s = 0; s < suppSize; ++s) {
891:             // Filter support
892:             if (GetPointIndex(supp[s], cStart, cEnd, cells) >= 0) ++Ns;
893:           }
894:           PetscCheck(Ns == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary face %" PetscInt_FMT " should see one cell, not %" PetscInt_FMT, face, Ns);
895:           if (view)
896:             PetscCall(PetscViewerASCIIPrintf(selfviewer, "[%d]: component %" PetscInt_FMT ", 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,
897:                                              rpoints[l].index, rrank, rcomp, lorntComp[face].rank));
898:           neighbors[comp][numNeighbors[comp]++] = l;
899:         }
900:       }
901:     }
902:     totNeighbors += numNeighbors[comp];
903:   }
904:   PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &selfviewer));
905:   if (viewSync) PetscCall(PetscViewerASCIIPopSynchronized(viewer));
906:   PetscCall(PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match));
907:   for (PetscInt comp = 0, off = 0; comp < Ncomp; ++comp) {
908:     for (PetscInt n = 0; n < numNeighbors[comp]; ++n, ++off) {
909:       const PetscInt face = lpoints[neighbors[comp][n]];
910:       const PetscInt o    = rorntComp[face].rank * lorntComp[face].rank;

912:       if (o < 0) match[off] = PETSC_TRUE;
913:       else if (o > 0) match[off] = PETSC_FALSE;
914:       else
915:         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid face %" PetscInt_FMT " (%" PetscInt_FMT ", %" PetscInt_FMT ") neighbor: %" PetscInt_FMT " comp: %" PetscInt_FMT, face, rorntComp[face].rank, lorntComp[face].rank, neighbors[comp][n], comp);
916:       nrankComp[off].rank  = rpoints[neighbors[comp][n]].rank;
917:       nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index;
918:     }
919:     PetscCall(PetscFree(neighbors[comp]));
920:   }
921:   /* Collect the graph on 0 */
922:   if (numLeaves >= 0) {
923:     Mat          G;
924:     PetscBT      seenProcs, flippedProcs;
925:     PetscInt    *procFIFO, pTop, pBottom;
926:     PetscInt    *N          = NULL, *Noff;
927:     PetscSFNode *adj        = NULL;
928:     PetscBool   *val        = NULL;
929:     PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc;
930:     PetscMPIInt  size = 0, iNcomp, itotNeighbors;

932:     PetscCall(PetscCalloc1(Ncomp, &flipped));
933:     if (rank == 0) PetscCallMPI(MPI_Comm_size(comm, &size));
934:     PetscCall(PetscCalloc4(size, &recvcounts, size + 1, &displs, size, &Nc, size + 1, &Noff));
935:     PetscCallMPI(MPI_Gather(&Ncomp, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm));
936:     for (PetscInt p = 0; p < size; ++p) displs[p + 1] = displs[p] + Nc[p];
937:     if (rank == 0) PetscCall(PetscMalloc1(displs[size], &N));
938:     PetscCall(PetscMPIIntCast(Ncomp, &iNcomp));
939:     PetscCallMPI(MPI_Gatherv(numNeighbors, iNcomp, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm));
940:     for (PetscInt p = 0, o = 0; p < size; ++p) {
941:       recvcounts[p] = 0;
942:       for (PetscInt c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o];
943:       displs[p + 1] = displs[p] + recvcounts[p];
944:     }
945:     if (rank == 0) PetscCall(PetscMalloc2(displs[size], &adj, displs[size], &val));
946:     PetscCall(PetscMPIIntCast(totNeighbors, &itotNeighbors));
947:     PetscCallMPI(MPI_Gatherv(nrankComp, itotNeighbors, MPIU_SF_NODE, adj, recvcounts, displs, MPIU_SF_NODE, 0, comm));
948:     PetscCallMPI(MPI_Gatherv(match, itotNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm));
949:     PetscCall(PetscFree2(numNeighbors, neighbors));
950:     if (rank == 0) {
951:       for (PetscInt p = 1; p <= size; ++p) Noff[p] = Noff[p - 1] + Nc[p - 1];
952:       if (view) {
953:         for (PetscInt p = 0, off = 0; p < size; ++p) {
954:           for (PetscInt c = 0; c < Nc[p]; ++c) {
955:             PetscCall(PetscPrintf(PETSC_COMM_SELF, "Proc %" PetscInt_FMT " Comp %" PetscInt_FMT ":\n", p, c));
956:             for (PetscInt n = 0; n < N[Noff[p] + c]; ++n, ++off) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  edge (%" PetscInt_FMT ", %" PetscInt_FMT ") (%s):\n", adj[off].rank, adj[off].index, PetscBools[val[off]]));
957:           }
958:         }
959:       }
960:       /* Symmetrize the graph */
961:       PetscCall(MatCreate(PETSC_COMM_SELF, &G));
962:       PetscCall(MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size]));
963:       PetscCall(MatSetUp(G));
964:       for (PetscInt p = 0, off = 0; p < size; ++p) {
965:         for (PetscInt c = 0; c < Nc[p]; ++c) {
966:           const PetscInt r = Noff[p] + c;

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

972:             PetscCall(MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES));
973:             PetscCall(MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES));
974:           }
975:         }
976:       }
977:       PetscCall(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY));
978:       PetscCall(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY));

980:       PetscCall(PetscBTCreate(Noff[size], &seenProcs));
981:       PetscCall(PetscBTMemzero(Noff[size], seenProcs));
982:       PetscCall(PetscBTCreate(Noff[size], &flippedProcs));
983:       PetscCall(PetscBTMemzero(Noff[size], flippedProcs));
984:       PetscCall(PetscMalloc1(Noff[size], &procFIFO));
985:       pTop = pBottom = 0;
986:       for (PetscInt p = 0; p < Noff[size]; ++p) {
987:         if (PetscBTLookup(seenProcs, p)) continue;
988:         /* Initialize FIFO with next proc */
989:         procFIFO[pBottom++] = p;
990:         PetscCall(PetscBTSet(seenProcs, p));
991:         /* Consider each proc in FIFO */
992:         while (pTop < pBottom) {
993:           const PetscScalar *ornt;
994:           const PetscInt    *neighbors;
995:           PetscInt           proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors;

997:           proc     = procFIFO[pTop++];
998:           flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0;
999:           PetscCall(MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt));
1000:           /* Loop over neighboring procs */
1001:           for (PetscInt n = 0; n < numNeighbors; ++n) {
1002:             nproc    = neighbors[n];
1003:             mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1;
1004:             seen     = PetscBTLookup(seenProcs, nproc);
1005:             flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0;

1007:             if (mismatch ^ (flippedA ^ flippedB)) {
1008:               PetscCheck(!seen, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Previously seen procs %" PetscInt_FMT " and %" PetscInt_FMT " do not match: Fault mesh is non-orientable", proc, nproc);
1009:               if (!flippedB) {
1010:                 PetscCall(PetscBTSet(flippedProcs, nproc));
1011:               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
1012:             } else PetscCheck(!mismatch || !flippedA || !flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
1013:             if (!seen) {
1014:               procFIFO[pBottom++] = nproc;
1015:               PetscCall(PetscBTSet(seenProcs, nproc));
1016:             }
1017:           }
1018:         }
1019:       }
1020:       PetscCall(PetscFree(procFIFO));
1021:       PetscCall(MatDestroy(&G));
1022:       PetscCall(PetscFree2(adj, val));
1023:       PetscCall(PetscBTDestroy(&seenProcs));
1024:     }
1025:     /* Scatter flip flags */
1026:     {
1027:       PetscBool *flips = NULL;

1029:       if (rank == 0) {
1030:         PetscCall(PetscMalloc1(Noff[size], &flips));
1031:         for (PetscInt p = 0; p < Noff[size]; ++p) {
1032:           flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE;
1033:           if (view && flips[p]) PetscCall(PetscPrintf(comm, "Flipping Proc+Comp %" PetscInt_FMT ":\n", p));
1034:         }
1035:         for (PetscInt p = 0; p < size; ++p) displs[p + 1] = displs[p] + Nc[p];
1036:       }
1037:       PetscCall(PetscMPIIntCast(Ncomp, &iNcomp));
1038:       PetscCallMPI(MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, iNcomp, MPIU_BOOL, 0, comm));
1039:       PetscCall(PetscFree(flips));
1040:     }
1041:     if (rank == 0) PetscCall(PetscBTDestroy(&flippedProcs));
1042:     PetscCall(PetscFree(N));
1043:     PetscCall(PetscFree4(recvcounts, displs, Nc, Noff));
1044:     PetscCall(PetscFree2(nrankComp, match));

1046:     /* Decide whether to flip cells in each component */
1047:     for (PetscInt c = 0; c < cEnd - cStart; ++c) {
1048:       if (flipped[cellComp[c]]) PetscCall(PetscBTNegate(flippedCells, c));
1049:     }
1050:     PetscCall(PetscFree(flipped));
1051:   }
1052:   if (view) {
1053:     PetscViewer v;

1055:     PetscCall(PetscViewerASCIIGetStdout(comm, &v));
1056:     PetscCall(PetscViewerASCIIPushSynchronized(v));
1057:     PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank));
1058:     PetscCall(PetscBTView(cEnd - cStart, flippedCells, v));
1059:     PetscCall(PetscViewerFlush(v));
1060:     PetscCall(PetscViewerASCIIPopSynchronized(v));
1061:   }
1062:   // Reverse flipped cells in the mesh
1063:   PetscViewer     v;
1064:   const PetscInt *degree = NULL;
1065:   PetscInt       *points;
1066:   PetscInt        pStart, pEnd;

1068:   if (view) {
1069:     PetscCall(PetscViewerASCIIGetStdout(comm, &v));
1070:     PetscCall(PetscViewerASCIIPushSynchronized(v));
1071:   }
1072:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
1073:   if (numRoots >= 0) {
1074:     PetscCall(PetscSFComputeDegreeBegin(sf, &degree));
1075:     PetscCall(PetscSFComputeDegreeEnd(sf, &degree));
1076:   }
1077:   PetscCall(PetscCalloc1(pEnd - pStart, &points));
1078:   for (PetscInt c = cStart; c < cEnd; ++c) {
1079:     if (PetscBTLookup(flippedCells, c - cStart)) {
1080:       const PetscInt cell = cells ? cells[c] : c;

1082:       PetscCall(DMPlexOrientPoint(dm, cell, -1));
1083:       if (degree && degree[cell]) points[cell] = 1;
1084:       if (view) PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]Flipping cell %" PetscInt_FMT "%s\n", rank, cell, degree && degree[cell] ? " and sending to overlap" : ""));
1085:     }
1086:   }
1087:   // Must propagate flips for cells in the overlap
1088:   if (numRoots >= 0) {
1089:     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, points, points, MPI_SUM));
1090:     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, points, points, MPI_SUM));
1091:   }
1092:   for (PetscInt c = cStart; c < cEnd; ++c) {
1093:     const PetscInt cell = cells ? cells[c] : c;

1095:     if (points[cell] && !PetscBTLookup(flippedCells, c - cStart)) {
1096:       PetscCall(DMPlexOrientPoint(dm, cell, -1));
1097:       if (view) PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]Flipping cell %" PetscInt_FMT " through overlap\n", rank, cell));
1098:     }
1099:   }
1100:   if (view) {
1101:     PetscCall(PetscViewerFlush(v));
1102:     PetscCall(PetscViewerASCIIPopSynchronized(v));
1103:   }
1104:   PetscCall(PetscFree(points));
1105:   PetscCall(PetscBTDestroy(&flippedCells));
1106:   PetscCall(PetscFree2(numNeighbors, neighbors));
1107:   PetscCall(PetscFree3(rorntComp, lorntComp, locSupp));
1108:   PetscCall(PetscFree2(cellComp, faceComp));
1109:   PetscFunctionReturn(PETSC_SUCCESS);
1110: }