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) PetscCall(PetscBTSet(flippedCells, indC[0]));
158:     else {
159:       PetscCheck(!seenB && !flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
160:       PetscCall(PetscBTSet(flippedCells, indC[1]));
161:     }
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:         if (lpoints) PetscCall(PetscFindInt(support[s], numLeaves, lpoints, &l));
464:         else {
465:           if (support[s] >= 0 && support[s] < numLeaves) l = support[s];
466:           else l = -1;
467:         }
468:         if (l >= 0) continue;
469:         locSupport[Ns++] = support[s];
470:       }
471:       if (Ns != 1) continue;
472:       neighbor = locSupport[0];
473:       PetscCall(DMPlexGetCone(dm, neighbor, &cone));
474:       PetscCall(DMPlexGetConeSize(dm, neighbor, &coneSize));
475:       PetscCall(DMPlexGetConeOrientation(dm, neighbor, &ornt));
476:       for (c = 0; c < coneSize; ++c)
477:         if (cone[c] == face) break;
478:       if (dim == 1) {
479:         /* Use cone position instead, shifted to -1 or 1 */
480:         if (PetscBTLookup(flippedCells, neighbor - cStart)) rorntComp[face].rank = 1 - c * 2;
481:         else rorntComp[face].rank = c * 2 - 1;
482:       } else {
483:         if (PetscBTLookup(flippedCells, neighbor - cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 : 1;
484:         else rorntComp[face].rank = ornt[c] < 0 ? 1 : -1;
485:       }
486:       rorntComp[face].index = faceComp[face - fStart];
487:     }
488:     /* Communicate boundary edge orientations */
489:     PetscCall(PetscSFBcastBegin(sf, MPIU_SF_NODE, rorntComp, lorntComp, MPI_REPLACE));
490:     PetscCall(PetscSFBcastEnd(sf, MPIU_SF_NODE, rorntComp, lorntComp, MPI_REPLACE));
491:   }
492:   /* Get process adjacency */
493:   PetscCall(PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors));
494:   viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm));
495:   if (flg2) PetscCall(PetscViewerASCIIPushSynchronized(viewer));
496:   PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &selfviewer));
497:   for (comp = 0; comp < numComponents; ++comp) {
498:     PetscInt n;

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

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

511:         for (n = 0; n < numNeighbors[comp]; ++n)
512:           if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break;
513:         if (n >= numNeighbors[comp]) {
514:           PetscInt supportSize;

516:           PetscCall(DMPlexGetSupportSize(dm, face, &supportSize));
517:           // We can have internal faces in the SF if we have cells in the SF
518:           if (supportSize > 1) continue;
519:           if (flg)
520:             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,
521:                                              rpoints[l].index, rrank, rcomp, lorntComp[face].rank));
522:           neighbors[comp][numNeighbors[comp]++] = l;
523:         }
524:       }
525:     }
526:     totNeighbors += numNeighbors[comp];
527:   }
528:   PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &selfviewer));
529:   if (flg2) PetscCall(PetscViewerASCIIPopSynchronized(viewer));
530:   PetscCall(PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match));
531:   for (comp = 0, off = 0; comp < numComponents; ++comp) {
532:     PetscInt n;

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

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

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

579:         for (p = 0, off = 0; p < size; ++p) {
580:           for (c = 0; c < Nc[p]; ++c) {
581:             PetscCall(PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %" PetscInt_FMT ":\n", p, c));
582:             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]]));
583:           }
584:         }
585:       }
586:       /* Symmetrize the graph */
587:       PetscCall(MatCreate(PETSC_COMM_SELF, &G));
588:       PetscCall(MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size]));
589:       PetscCall(MatSetUp(G));
590:       for (p = 0, off = 0; p < size; ++p) {
591:         for (c = 0; c < Nc[p]; ++c) {
592:           const PetscInt r = Noff[p] + c;

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

598:             PetscCall(MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES));
599:             PetscCall(MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES));
600:           }
601:         }
602:       }
603:       PetscCall(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY));
604:       PetscCall(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY));

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

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

633:             if (mismatch ^ (flippedA ^ flippedB)) {
634:               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);
635:               PetscCheck(!flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
636:               PetscCall(PetscBTSet(flippedProcs, nproc));
637:             } else PetscCheck(!mismatch || !flippedA || !flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
638:             if (!seen) {
639:               procFIFO[pBottom++] = nproc;
640:               PetscCall(PetscBTSet(seenProcs, nproc));
641:             }
642:           }
643:         }
644:       }
645:       PetscCall(PetscFree(procFIFO));
646:       PetscCall(MatDestroy(&G));
647:       PetscCall(PetscFree2(adj, val));
648:       PetscCall(PetscBTDestroy(&seenProcs));
649:     }
650:     /* Scatter flip flags */
651:     {
652:       PetscBool *flips = NULL;

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

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

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

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

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

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

727: PetscErrorCode DMPlexOrientLabel(DM dm, DMLabel label)
728: {
729:   IS cellIS, faceIS;

731:   PetscFunctionBegin;
732:   PetscCall(CreateCellAndFaceIS_Private(dm, label, &cellIS, &faceIS));
733:   PetscCall(DMPlexOrientCells_Internal(dm, cellIS, faceIS));
734:   PetscCall(ISDestroy(&cellIS));
735:   PetscCall(ISDestroy(&faceIS));
736:   PetscFunctionReturn(PETSC_SUCCESS);
737: }

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

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

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

798:    NEW STYLE
799:    - Create the adj on each process
800:    - Bootstrap to complete graph on proc 0
801:   */
802:   viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm));
803:   if (viewSync) PetscCall(PetscViewerASCIIPushSynchronized(viewer));
804:   PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &selfviewer));
805:   if (faceIsVertex && lpoints) {
806:     // Need to first flip cells which hit parallel boundary on wrong face
807:     for (PetscInt c = cStart; c < cEnd; ++c) {
808:       const PetscInt  cell = cells ? cells[c] : c;
809:       const PetscInt *cone;
810:       PetscInt        cS, ls, le;

812:       PetscCall(DMPlexGetCone(dm, cell, &cone));
813:       PetscCall(DMPlexGetConeSize(dm, cell, &cS));
814:       PetscCheck(cS == 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Size of cone for edge %" PetscInt_FMT " must be 2, not %" PetscInt_FMT, cell, cS);
815:       PetscCall(PetscFindInt(cone[0], numLeaves, lpoints, &ls));
816:       PetscCall(PetscFindInt(cone[1], numLeaves, lpoints, &le));
817:       if (ls >= 0 && le < 0) {
818:         PetscCall(PetscBTSet(flippedCells, c - cStart));
819:         if (view) PetscCall(PetscViewerASCIIPrintf(selfviewer, "[%d]: Flipped cell %" PetscInt_FMT " to meet parallel boundary on shared face %" PetscInt_FMT "\n", rank, cell, cone[0]));
820:       }
821:     }
822:   }
823:   PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &selfviewer));
824:   if (viewSync) PetscCall(PetscViewerASCIIPopSynchronized(viewer));
825:   PetscCall(DMPlexOrient_Serial(dm, cellIS, faceIS, &Ncomp, cellComp, faceComp, flippedCells));
826:   if (view) {
827:     PetscViewer v;
828:     PetscInt    cdepth = -1;

830:     PetscCall(PetscViewerASCIIGetStdout(comm, &v));
831:     PetscCall(PetscViewerASCIIPushSynchronized(v));
832:     if (cEnd > cStart) PetscCall(DMPlexGetPointDepth(dm, cells ? cells[cStart] : cStart, &cdepth));
833:     PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]New Orientation %" PetscInt_FMT " cells (depth %" PetscInt_FMT ") and %" PetscInt_FMT " faces\n", rank, cEnd - cStart, cdepth, fEnd - fStart));
834:     PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank));
835:     PetscCall(PetscBTView(cEnd - cStart, flippedCells, v));
836:     PetscCall(PetscViewerFlush(v));
837:     PetscCall(PetscViewerASCIIPopSynchronized(v));
838:   }
839:   if (viewSync) PetscCall(PetscViewerASCIIPushSynchronized(viewer));
840:   PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &selfviewer));
841:   /* Now all subdomains are oriented, but we need a consistent parallel orientation */
842:   // TODO: This all has to be rewritten to filter cones/supports to the ISes
843:   if (numLeaves >= 0) {
844:     PetscInt maxSuppSize, neighbor;

846:     // Store orientations of boundary faces
847:     PetscCall(DMPlexGetMaxSizes(dm, NULL, &maxSuppSize));
848:     PetscCall(PetscCalloc3(numRoots, &rorntComp, numRoots, &lorntComp, maxSuppSize, &locSupp));
849:     for (PetscInt f = fStart; f < fEnd; ++f) {
850:       const PetscInt  face = faces ? faces[f] : f;
851:       const PetscInt *cone, *supp, *ornt;
852:       PetscInt        coneSize, suppSize, nind, c, Ns = 0;

854:       PetscCall(DMPlexGetSupportSize(dm, face, &suppSize));
855:       PetscCall(DMPlexGetSupport(dm, face, &supp));
856:       for (PetscInt s = 0; s < suppSize; ++s) {
857:         PetscInt ind, l;

859:         // Filter support
860:         ind = GetPointIndex(supp[s], cStart, cEnd, cells);
861:         if (ind < 0) continue;
862:         // Ignore overlapping cells
863:         PetscCall(PetscFindInt(supp[s], numLeaves, lpoints, &l));
864:         if (l >= 0) continue;
865:         locSupp[Ns++] = supp[s];
866:       }
867:       PetscCheck(Ns < maxSuppSize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Index %" PetscInt_FMT " exceeds array size %" PetscInt_FMT, Ns, maxSuppSize);
868:       if (Ns != 1) continue;
869:       neighbor = locSupp[0];
870:       nind     = GetPointIndex(neighbor, cStart, cEnd, cells);
871:       PetscCall(DMPlexGetCone(dm, neighbor, &cone));
872:       PetscCall(DMPlexGetConeSize(dm, neighbor, &coneSize));
873:       PetscCall(DMPlexGetConeOrientation(dm, neighbor, &ornt));
874:       for (c = 0; c < coneSize; ++c)
875:         if (cone[c] == face) break;
876:       if (faceIsVertex) {
877:         /* Use cone position instead, shifted to -1 or 1 */
878:         if (PetscBTLookup(flippedCells, nind)) rorntComp[face].rank = 1 - c * 2;
879:         else rorntComp[face].rank = c * 2 - 1;
880:       } else {
881:         if (PetscBTLookup(flippedCells, nind)) rorntComp[face].rank = ornt[c] < 0 ? -1 : 1;
882:         else rorntComp[face].rank = ornt[c] < 0 ? 1 : -1;
883:       }
884:       rorntComp[face].index = faceComp[GetPointIndex(face, fStart, fEnd, faces)];
885:       if (view) PetscCall(PetscViewerASCIIPrintf(selfviewer, "[%d]: Boundary face %" PetscInt_FMT " component %" PetscInt_FMT " orientation %" PetscInt_FMT "\n", rank, face, rorntComp[face].index, rorntComp[face].rank));
886:     }
887:     // Communicate boundary edge orientations
888:     PetscCall(PetscSFBcastBegin(sf, MPIU_SF_NODE, rorntComp, lorntComp, MPI_REPLACE));
889:     PetscCall(PetscSFBcastEnd(sf, MPIU_SF_NODE, rorntComp, lorntComp, MPI_REPLACE));
890:   }
891:   /* Get process adjacency */
892:   PetscCall(PetscMalloc2(Ncomp, &numNeighbors, Ncomp, &neighbors));
893:   for (PetscInt comp = 0; comp < Ncomp; ++comp) {
894:     PetscInt n;

896:     numNeighbors[comp] = 0;
897:     PetscCall(PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp]));
898:     /* I know this is p^2 time in general, but for bounded degree its alright */
899:     for (PetscInt l = 0; l < numLeaves; ++l) {
900:       const PetscInt face = lpoints[l];
901:       PetscInt       find;

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

909:         for (n = 0; n < numNeighbors[comp]; ++n)
910:           if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break;
911:         if (n >= numNeighbors[comp]) {
912:           const PetscInt *supp;
913:           PetscInt        suppSize, Ns = 0;

915:           PetscCall(DMPlexGetSupport(dm, face, &supp));
916:           PetscCall(DMPlexGetSupportSize(dm, face, &suppSize));
917:           for (PetscInt s = 0; s < suppSize; ++s) {
918:             // Filter support
919:             if (GetPointIndex(supp[s], cStart, cEnd, cells) >= 0) ++Ns;
920:           }
921:           // Skip if the face is in the interior
922:           if (Ns > 1) continue;
923:           // Skip if there is a lone vertex
924:           if (!lorntComp[face].rank) continue;
925:           if (view)
926:             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,
927:                                              rpoints[l].index, rrank, rcomp, lorntComp[face].rank));
928:           neighbors[comp][numNeighbors[comp]++] = l;
929:         }
930:       }
931:     }
932:     totNeighbors += numNeighbors[comp];
933:   }
934:   PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &selfviewer));
935:   if (viewSync) PetscCall(PetscViewerASCIIPopSynchronized(viewer));
936:   PetscCall(PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match));
937:   for (PetscInt comp = 0, off = 0; comp < Ncomp; ++comp) {
938:     for (PetscInt n = 0; n < numNeighbors[comp]; ++n, ++off) {
939:       const PetscInt face = lpoints[neighbors[comp][n]];
940:       const PetscInt o    = rorntComp[face].rank * lorntComp[face].rank;

942:       if (faceIsVertex) {
943:         match[off] = rorntComp[face].rank != lorntComp[face].rank ? PETSC_TRUE : PETSC_FALSE;
944:       } else {
945:         if (o < 0) match[off] = PETSC_TRUE;
946:         else if (o > 0) match[off] = PETSC_FALSE;
947:         else
948:           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);
949:       }
950:       nrankComp[off].rank  = rpoints[neighbors[comp][n]].rank;
951:       nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index;
952:     }
953:     PetscCall(PetscFree(neighbors[comp]));
954:   }
955:   /* Collect the graph on 0 */
956:   if (numLeaves >= 0) {
957:     Mat          G;
958:     PetscBT      seenProcs, flippedProcs;
959:     PetscInt    *procFIFO, pTop, pBottom;
960:     PetscInt    *N          = NULL, *Noff;
961:     PetscSFNode *adj        = NULL;
962:     PetscBool   *val        = NULL;
963:     PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc;
964:     PetscMPIInt  size = 0, iNcomp, itotNeighbors;

966:     PetscCall(PetscCalloc1(Ncomp, &flipped));
967:     if (rank == 0) PetscCallMPI(MPI_Comm_size(comm, &size));
968:     PetscCall(PetscCalloc4(size, &recvcounts, size + 1, &displs, size, &Nc, size + 1, &Noff));
969:     PetscCallMPI(MPI_Gather(&Ncomp, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm));
970:     for (PetscInt p = 0; p < size; ++p) displs[p + 1] = displs[p] + Nc[p];
971:     if (rank == 0) PetscCall(PetscMalloc1(displs[size], &N));
972:     PetscCall(PetscMPIIntCast(Ncomp, &iNcomp));
973:     PetscCallMPI(MPI_Gatherv(numNeighbors, iNcomp, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm));
974:     for (PetscInt p = 0, o = 0; p < size; ++p) {
975:       recvcounts[p] = 0;
976:       for (PetscInt c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o];
977:       displs[p + 1] = displs[p] + recvcounts[p];
978:     }
979:     if (rank == 0) PetscCall(PetscMalloc2(displs[size], &adj, displs[size], &val));
980:     PetscCall(PetscMPIIntCast(totNeighbors, &itotNeighbors));
981:     PetscCallMPI(MPI_Gatherv(nrankComp, itotNeighbors, MPIU_SF_NODE, adj, recvcounts, displs, MPIU_SF_NODE, 0, comm));
982:     PetscCallMPI(MPI_Gatherv(match, itotNeighbors, MPI_C_BOOL, val, recvcounts, displs, MPI_C_BOOL, 0, comm));
983:     PetscCall(PetscFree2(numNeighbors, neighbors));
984:     if (rank == 0) {
985:       for (PetscInt p = 1; p <= size; ++p) Noff[p] = Noff[p - 1] + Nc[p - 1];
986:       if (view) {
987:         for (PetscInt p = 0, off = 0; p < size; ++p) {
988:           for (PetscInt c = 0; c < Nc[p]; ++c) {
989:             PetscCall(PetscPrintf(PETSC_COMM_SELF, "Proc %" PetscInt_FMT " Comp %" PetscInt_FMT ":\n", p, c));
990:             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]]));
991:           }
992:         }
993:       }
994:       /* Symmetrize the graph */
995:       PetscCall(MatCreate(PETSC_COMM_SELF, &G));
996:       PetscCall(MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size]));
997:       PetscCall(MatSetUp(G));
998:       for (PetscInt p = 0, off = 0; p < size; ++p) {
999:         for (PetscInt c = 0; c < Nc[p]; ++c) {
1000:           const PetscInt r = Noff[p] + c;

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

1006:             // Do not set values for processes that have no face to orient
1007:             if (!Nc[adj[off].rank]) continue;
1008:             PetscCall(MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES));
1009:             PetscCall(MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES));
1010:           }
1011:         }
1012:       }
1013:       PetscCall(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY));
1014:       PetscCall(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY));

1016:       PetscCall(PetscBTCreate(Noff[size], &seenProcs));
1017:       PetscCall(PetscBTMemzero(Noff[size], seenProcs));
1018:       PetscCall(PetscBTCreate(Noff[size], &flippedProcs));
1019:       PetscCall(PetscBTMemzero(Noff[size], flippedProcs));
1020:       PetscCall(PetscMalloc1(Noff[size], &procFIFO));
1021:       pTop = pBottom = 0;
1022:       for (PetscInt p = 0; p < Noff[size]; ++p) {
1023:         if (PetscBTLookup(seenProcs, p)) continue;
1024:         /* Initialize FIFO with next proc */
1025:         procFIFO[pBottom++] = p;
1026:         PetscCall(PetscBTSet(seenProcs, p));
1027:         /* Consider each proc in FIFO */
1028:         while (pTop < pBottom) {
1029:           const PetscScalar *ornt;
1030:           const PetscInt    *neighbors;
1031:           PetscInt           proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors;

1033:           proc     = procFIFO[pTop++];
1034:           flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0;
1035:           PetscCall(MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt));
1036:           /* Loop over neighboring procs */
1037:           for (PetscInt n = 0; n < numNeighbors; ++n) {
1038:             nproc    = neighbors[n];
1039:             mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1;
1040:             seen     = PetscBTLookup(seenProcs, nproc);
1041:             flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0;

1043:             if (mismatch ^ (flippedA ^ flippedB)) {
1044:               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);
1045:               PetscCheck(!flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable");
1046:               PetscCall(PetscBTSet(flippedProcs, nproc));
1047:             } else PetscCheck(!mismatch || !flippedA || !flippedB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable");
1048:             if (!seen) {
1049:               procFIFO[pBottom++] = nproc;
1050:               PetscCall(PetscBTSet(seenProcs, nproc));
1051:             }
1052:           }
1053:         }
1054:       }
1055:       PetscCall(PetscFree(procFIFO));
1056:       PetscCall(MatDestroy(&G));
1057:       PetscCall(PetscFree2(adj, val));
1058:       PetscCall(PetscBTDestroy(&seenProcs));
1059:     }
1060:     /* Scatter flip flags */
1061:     {
1062:       PetscBool *flips = NULL;

1064:       if (rank == 0) {
1065:         PetscCall(PetscMalloc1(Noff[size], &flips));
1066:         for (PetscInt p = 0; p < Noff[size]; ++p) {
1067:           flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE;
1068:           if (view && flips[p]) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Flipping Proc+Comp %" PetscInt_FMT ":\n", p));
1069:         }
1070:         for (PetscInt p = 0; p < size; ++p) displs[p + 1] = displs[p] + Nc[p];
1071:       }
1072:       PetscCall(PetscMPIIntCast(Ncomp, &iNcomp));
1073:       PetscCallMPI(MPI_Scatterv(flips, Nc, displs, MPI_C_BOOL, flipped, iNcomp, MPI_C_BOOL, 0, comm));
1074:       PetscCall(PetscFree(flips));
1075:     }
1076:     if (rank == 0) PetscCall(PetscBTDestroy(&flippedProcs));
1077:     PetscCall(PetscFree(N));
1078:     PetscCall(PetscFree4(recvcounts, displs, Nc, Noff));
1079:     PetscCall(PetscFree2(nrankComp, match));

1081:     /* Decide whether to flip cells in each component */
1082:     for (PetscInt c = 0; c < cEnd - cStart; ++c) {
1083:       if (flipped[cellComp[c]]) PetscCall(PetscBTNegate(flippedCells, c));
1084:     }
1085:     PetscCall(PetscFree(flipped));
1086:   }
1087:   if (view) {
1088:     PetscViewer v;

1090:     PetscCall(PetscViewerASCIIGetStdout(comm, &v));
1091:     PetscCall(PetscViewerASCIIPushSynchronized(v));
1092:     PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank));
1093:     PetscCall(PetscBTView(cEnd - cStart, flippedCells, v));
1094:     PetscCall(PetscViewerFlush(v));
1095:     PetscCall(PetscViewerASCIIPopSynchronized(v));
1096:   }
1097:   // Reverse flipped cells in the mesh
1098:   PetscViewer     v;
1099:   const PetscInt *degree = NULL;
1100:   PetscInt       *points;
1101:   PetscInt        pStart, pEnd;

1103:   if (view) {
1104:     PetscCall(PetscViewerASCIIGetStdout(comm, &v));
1105:     PetscCall(PetscViewerASCIIPushSynchronized(v));
1106:   }
1107:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
1108:   if (numRoots >= 0) {
1109:     PetscCall(PetscSFComputeDegreeBegin(sf, &degree));
1110:     PetscCall(PetscSFComputeDegreeEnd(sf, &degree));
1111:   }
1112:   PetscCall(PetscCalloc1(pEnd - pStart, &points));
1113:   for (PetscInt c = cStart; c < cEnd; ++c) {
1114:     if (PetscBTLookup(flippedCells, c - cStart)) {
1115:       const PetscInt cell = cells ? cells[c] : c;

1117:       PetscCall(DMPlexOrientPoint(dm, cell, -1));
1118:       if (degree && degree[cell]) points[cell] = 1;
1119:       if (view) PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]Flipping cell %" PetscInt_FMT "%s\n", rank, cell, degree && degree[cell] ? " and sending to overlap" : ""));
1120:     }
1121:   }
1122:   // Must propagate flips for cells in the overlap
1123:   if (numRoots >= 0) {
1124:     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, points, points, MPI_SUM));
1125:     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, points, points, MPI_SUM));
1126:   }
1127:   for (PetscInt c = cStart; c < cEnd; ++c) {
1128:     const PetscInt cell = cells ? cells[c] : c;

1130:     if (points[cell] && !PetscBTLookup(flippedCells, c - cStart)) {
1131:       PetscCall(DMPlexOrientPoint(dm, cell, -1));
1132:       if (view) PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d]Flipping cell %" PetscInt_FMT " through overlap\n", rank, cell));
1133:     }
1134:   }
1135:   if (view) {
1136:     PetscCall(PetscViewerFlush(v));
1137:     PetscCall(PetscViewerASCIIPopSynchronized(v));
1138:   }
1139:   PetscCall(PetscFree(points));
1140:   PetscCall(PetscBTDestroy(&flippedCells));
1141:   PetscCall(PetscFree2(numNeighbors, neighbors));
1142:   PetscCall(PetscFree3(rorntComp, lorntComp, locSupp));
1143:   PetscCall(PetscFree2(cellComp, faceComp));
1144:   PetscFunctionReturn(PETSC_SUCCESS);
1145: }