Actual source code: isdiff.c


  2: #include <petsc/private/isimpl.h>
  3: #include <petsc/private/sectionimpl.h>
  4: #include <petscbt.h>

  6: /*@
  7:    ISDifference - Computes the difference between two index sets.

  9:    Collective on is1

 11:    Input Parameters:
 12: +  is1 - first index, to have items removed from it
 13: -  is2 - index values to be removed

 15:    Output Parameters:
 16: .  isout - is1 - is2

 18:    Level: intermediate

 20:    Notes:
 21:    Negative values are removed from the lists. is2 may have values
 22:    that are not in is1. This requires O(imax-imin) memory and O(imax-imin)
 23:    work, where imin and imax are the bounds on the indices in is1.

 25:    If is2 is NULL, the result is the same as for an empty IS, i.e., a duplicate of is1.

 27:    The difference is computed separately on each MPI rank

 29: .seealso: [](sec_scatter), `IS`, `ISDestroy()`, `ISView()`, `ISSum()`, `ISExpand()`
 30: @*/
 31: PetscErrorCode ISDifference(IS is1, IS is2, IS *isout)
 32: {
 33:   PetscInt        i, n1, n2, imin, imax, nout, *iout;
 34:   const PetscInt *i1, *i2;
 35:   PetscBT         mask;
 36:   MPI_Comm        comm;

 38:   PetscFunctionBegin;
 41:   if (!is2) {
 42:     PetscCall(ISDuplicate(is1, isout));
 43:     PetscFunctionReturn(PETSC_SUCCESS);
 44:   }

 47:   PetscCall(ISGetIndices(is1, &i1));
 48:   PetscCall(ISGetLocalSize(is1, &n1));

 50:   /* Create a bit mask array to contain required values */
 51:   if (n1) {
 52:     imin = PETSC_MAX_INT;
 53:     imax = 0;
 54:     for (i = 0; i < n1; i++) {
 55:       if (i1[i] < 0) continue;
 56:       imin = PetscMin(imin, i1[i]);
 57:       imax = PetscMax(imax, i1[i]);
 58:     }
 59:   } else imin = imax = 0;

 61:   PetscCall(PetscBTCreate(imax - imin, &mask));
 62:   /* Put the values from is1 */
 63:   for (i = 0; i < n1; i++) {
 64:     if (i1[i] < 0) continue;
 65:     PetscCall(PetscBTSet(mask, i1[i] - imin));
 66:   }
 67:   PetscCall(ISRestoreIndices(is1, &i1));
 68:   /* Remove the values from is2 */
 69:   PetscCall(ISGetIndices(is2, &i2));
 70:   PetscCall(ISGetLocalSize(is2, &n2));
 71:   for (i = 0; i < n2; i++) {
 72:     if (i2[i] < imin || i2[i] > imax) continue;
 73:     PetscCall(PetscBTClear(mask, i2[i] - imin));
 74:   }
 75:   PetscCall(ISRestoreIndices(is2, &i2));

 77:   /* Count the number in the difference */
 78:   nout = 0;
 79:   for (i = 0; i < imax - imin + 1; i++) {
 80:     if (PetscBTLookup(mask, i)) nout++;
 81:   }

 83:   /* create the new IS containing the difference */
 84:   PetscCall(PetscMalloc1(nout, &iout));
 85:   nout = 0;
 86:   for (i = 0; i < imax - imin + 1; i++) {
 87:     if (PetscBTLookup(mask, i)) iout[nout++] = i + imin;
 88:   }
 89:   PetscCall(PetscObjectGetComm((PetscObject)is1, &comm));
 90:   PetscCall(ISCreateGeneral(comm, nout, iout, PETSC_OWN_POINTER, isout));

 92:   PetscCall(PetscBTDestroy(&mask));
 93:   PetscFunctionReturn(PETSC_SUCCESS);
 94: }

 96: /*@
 97:    ISSum - Computes the sum (union) of two index sets.

 99:    Only sequential version (at the moment)

101:    Input Parameters:
102: +  is1 - index set to be extended
103: -  is2 - index values to be added

105:    Output Parameter:
106: .   is3 - the sum; this can not be is1 or is2

108:    Level: intermediate

110:    Notes:
111:    If n1 and n2 are the sizes of the sets, this takes O(n1+n2) time;

113:    Both index sets need to be sorted on input.

115:    The sum is computed separately on each MPI rank

117: .seealso: [](sec_scatter), `IS`, `ISDestroy()`, `ISView()`, `ISDifference()`, `ISExpand()`
118: @*/
119: PetscErrorCode ISSum(IS is1, IS is2, IS *is3)
120: {
121:   MPI_Comm        comm;
122:   PetscBool       f;
123:   PetscMPIInt     size;
124:   const PetscInt *i1, *i2;
125:   PetscInt        n1, n2, n3, p1, p2, *iout;

127:   PetscFunctionBegin;
130:   PetscCall(PetscObjectGetComm((PetscObject)(is1), &comm));
131:   PetscCallMPI(MPI_Comm_size(comm, &size));
132:   PetscCheck(size <= 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Currently only for uni-processor IS");

134:   PetscCall(ISSorted(is1, &f));
135:   PetscCheck(f, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Arg 1 is not sorted");
136:   PetscCall(ISSorted(is2, &f));
137:   PetscCheck(f, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Arg 2 is not sorted");

139:   PetscCall(ISGetLocalSize(is1, &n1));
140:   PetscCall(ISGetLocalSize(is2, &n2));
141:   if (!n2) {
142:     PetscCall(ISDuplicate(is1, is3));
143:     PetscFunctionReturn(PETSC_SUCCESS);
144:   }
145:   PetscCall(ISGetIndices(is1, &i1));
146:   PetscCall(ISGetIndices(is2, &i2));

148:   p1 = 0;
149:   p2 = 0;
150:   n3 = 0;
151:   do {
152:     if (p1 == n1) { /* cleanup for is2 */
153:       n3 += n2 - p2;
154:       break;
155:     } else {
156:       while (p2 < n2 && i2[p2] < i1[p1]) {
157:         n3++;
158:         p2++;
159:       }
160:       if (p2 == n2) {
161:         /* cleanup for is1 */
162:         n3 += n1 - p1;
163:         break;
164:       } else {
165:         if (i2[p2] == i1[p1]) {
166:           n3++;
167:           p1++;
168:           p2++;
169:         }
170:       }
171:     }
172:     if (p2 == n2) {
173:       /* cleanup for is1 */
174:       n3 += n1 - p1;
175:       break;
176:     } else {
177:       while (p1 < n1 && i1[p1] < i2[p2]) {
178:         n3++;
179:         p1++;
180:       }
181:       if (p1 == n1) {
182:         /* clean up for is2 */
183:         n3 += n2 - p2;
184:         break;
185:       } else {
186:         if (i1[p1] == i2[p2]) {
187:           n3++;
188:           p1++;
189:           p2++;
190:         }
191:       }
192:     }
193:   } while (p1 < n1 || p2 < n2);

195:   if (n3 == n1) { /* no new elements to be added */
196:     PetscCall(ISRestoreIndices(is1, &i1));
197:     PetscCall(ISRestoreIndices(is2, &i2));
198:     PetscCall(ISDuplicate(is1, is3));
199:     PetscFunctionReturn(PETSC_SUCCESS);
200:   }
201:   PetscCall(PetscMalloc1(n3, &iout));

203:   p1 = 0;
204:   p2 = 0;
205:   n3 = 0;
206:   do {
207:     if (p1 == n1) { /* cleanup for is2 */
208:       while (p2 < n2) iout[n3++] = i2[p2++];
209:       break;
210:     } else {
211:       while (p2 < n2 && i2[p2] < i1[p1]) iout[n3++] = i2[p2++];
212:       if (p2 == n2) { /* cleanup for is1 */
213:         while (p1 < n1) iout[n3++] = i1[p1++];
214:         break;
215:       } else {
216:         if (i2[p2] == i1[p1]) {
217:           iout[n3++] = i1[p1++];
218:           p2++;
219:         }
220:       }
221:     }
222:     if (p2 == n2) { /* cleanup for is1 */
223:       while (p1 < n1) iout[n3++] = i1[p1++];
224:       break;
225:     } else {
226:       while (p1 < n1 && i1[p1] < i2[p2]) iout[n3++] = i1[p1++];
227:       if (p1 == n1) { /* clean up for is2 */
228:         while (p2 < n2) iout[n3++] = i2[p2++];
229:         break;
230:       } else {
231:         if (i1[p1] == i2[p2]) {
232:           iout[n3++] = i1[p1++];
233:           p2++;
234:         }
235:       }
236:     }
237:   } while (p1 < n1 || p2 < n2);

239:   PetscCall(ISRestoreIndices(is1, &i1));
240:   PetscCall(ISRestoreIndices(is2, &i2));
241:   PetscCall(ISCreateGeneral(comm, n3, iout, PETSC_OWN_POINTER, is3));
242:   PetscFunctionReturn(PETSC_SUCCESS);
243: }

245: /*@
246:    ISExpand - Computes the union of two index sets, by concatenating 2 lists and
247:    removing duplicates.

249:    Collective on is1

251:    Input Parameters:
252: +  is1 - first index set
253: -  is2 - index values to be added

255:    Output Parameters:
256: .  isout - is1 + is2 The index set is2 is appended to is1 removing duplicates

258:    Level: intermediate

260:    Notes:
261:    Negative values are removed from the lists. This requires O(imax-imin)
262:    memory and O(imax-imin) work, where imin and imax are the bounds on the
263:    indices in is1 and is2.

265:    The IS's do not need to be sorted.

267:    The operations are performed separately on each MPI rank

269: .seealso: [](sec_scatter), `IS`, `ISDestroy()`, `ISView()`, `ISDifference()`, `ISSum()`, `ISIntersect()`
270: @*/
271: PetscErrorCode ISExpand(IS is1, IS is2, IS *isout)
272: {
273:   PetscInt        i, n1, n2, imin, imax, nout, *iout;
274:   const PetscInt *i1, *i2;
275:   PetscBT         mask;
276:   MPI_Comm        comm;

278:   PetscFunctionBegin;

283:   PetscCheck(is1 || is2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Both arguments cannot be NULL");
284:   if (!is1) {
285:     PetscCall(ISDuplicate(is2, isout));
286:     PetscFunctionReturn(PETSC_SUCCESS);
287:   }
288:   if (!is2) {
289:     PetscCall(ISDuplicate(is1, isout));
290:     PetscFunctionReturn(PETSC_SUCCESS);
291:   }
292:   PetscCall(ISGetIndices(is1, &i1));
293:   PetscCall(ISGetLocalSize(is1, &n1));
294:   PetscCall(ISGetIndices(is2, &i2));
295:   PetscCall(ISGetLocalSize(is2, &n2));

297:   /* Create a bit mask array to contain required values */
298:   if (n1 || n2) {
299:     imin = PETSC_MAX_INT;
300:     imax = 0;
301:     for (i = 0; i < n1; i++) {
302:       if (i1[i] < 0) continue;
303:       imin = PetscMin(imin, i1[i]);
304:       imax = PetscMax(imax, i1[i]);
305:     }
306:     for (i = 0; i < n2; i++) {
307:       if (i2[i] < 0) continue;
308:       imin = PetscMin(imin, i2[i]);
309:       imax = PetscMax(imax, i2[i]);
310:     }
311:   } else imin = imax = 0;

313:   PetscCall(PetscMalloc1(n1 + n2, &iout));
314:   nout = 0;
315:   PetscCall(PetscBTCreate(imax - imin, &mask));
316:   /* Put the values from is1 */
317:   for (i = 0; i < n1; i++) {
318:     if (i1[i] < 0) continue;
319:     if (!PetscBTLookupSet(mask, i1[i] - imin)) iout[nout++] = i1[i];
320:   }
321:   PetscCall(ISRestoreIndices(is1, &i1));
322:   /* Put the values from is2 */
323:   for (i = 0; i < n2; i++) {
324:     if (i2[i] < 0) continue;
325:     if (!PetscBTLookupSet(mask, i2[i] - imin)) iout[nout++] = i2[i];
326:   }
327:   PetscCall(ISRestoreIndices(is2, &i2));

329:   /* create the new IS containing the sum */
330:   PetscCall(PetscObjectGetComm((PetscObject)is1, &comm));
331:   PetscCall(ISCreateGeneral(comm, nout, iout, PETSC_OWN_POINTER, isout));

333:   PetscCall(PetscBTDestroy(&mask));
334:   PetscFunctionReturn(PETSC_SUCCESS);
335: }

337: /*@
338:    ISIntersect - Computes the intersection of two index sets, by sorting and comparing.

340:    Collective on is1

342:    Input Parameters:
343: +  is1 - first index set
344: -  is2 - second index set

346:    Output Parameters:
347: .  isout - the sorted intersection of is1 and is2

349:    Level: intermediate

351:    Notes:
352:    Negative values are removed from the lists. This requires O(min(is1,is2))
353:    memory and O(max(is1,is2)log(max(is1,is2))) work

355:    The IS's do not need to be sorted.

357:    The operations are performed separately on each MPI rank

359: .seealso: [](sec_scatter), `IS`, `ISDestroy()`, `ISView()`, `ISDifference()`, `ISSum()`, `ISExpand()`, `ISConcatenate()`
360: @*/
361: PetscErrorCode ISIntersect(IS is1, IS is2, IS *isout)
362: {
363:   PetscInt        i, n1, n2, nout, *iout;
364:   const PetscInt *i1, *i2;
365:   IS              is1sorted = NULL, is2sorted = NULL;
366:   PetscBool       sorted, lsorted;
367:   MPI_Comm        comm;

369:   PetscFunctionBegin;
372:   PetscCheckSameComm(is1, 1, is2, 2);
374:   PetscCall(PetscObjectGetComm((PetscObject)is1, &comm));

376:   PetscCall(ISGetLocalSize(is1, &n1));
377:   PetscCall(ISGetLocalSize(is2, &n2));
378:   if (n1 < n2) {
379:     IS       tempis = is1;
380:     PetscInt ntemp  = n1;

382:     is1 = is2;
383:     is2 = tempis;
384:     n1  = n2;
385:     n2  = ntemp;
386:   }
387:   PetscCall(ISSorted(is1, &lsorted));
388:   PetscCall(MPIU_Allreduce(&lsorted, &sorted, 1, MPIU_BOOL, MPI_LAND, comm));
389:   if (!sorted) {
390:     PetscCall(ISDuplicate(is1, &is1sorted));
391:     PetscCall(ISSort(is1sorted));
392:     PetscCall(ISGetIndices(is1sorted, &i1));
393:   } else {
394:     is1sorted = is1;
395:     PetscCall(PetscObjectReference((PetscObject)is1));
396:     PetscCall(ISGetIndices(is1, &i1));
397:   }
398:   PetscCall(ISSorted(is2, &lsorted));
399:   PetscCall(MPIU_Allreduce(&lsorted, &sorted, 1, MPIU_BOOL, MPI_LAND, comm));
400:   if (!sorted) {
401:     PetscCall(ISDuplicate(is2, &is2sorted));
402:     PetscCall(ISSort(is2sorted));
403:     PetscCall(ISGetIndices(is2sorted, &i2));
404:   } else {
405:     is2sorted = is2;
406:     PetscCall(PetscObjectReference((PetscObject)is2));
407:     PetscCall(ISGetIndices(is2, &i2));
408:   }

410:   PetscCall(PetscMalloc1(n2, &iout));

412:   for (i = 0, nout = 0; i < n2; i++) {
413:     PetscInt key = i2[i];
414:     PetscInt loc;

416:     PetscCall(ISLocate(is1sorted, key, &loc));
417:     if (loc >= 0) {
418:       if (!nout || iout[nout - 1] < key) iout[nout++] = key;
419:     }
420:   }
421:   PetscCall(PetscRealloc(nout * sizeof(PetscInt), &iout));

423:   /* create the new IS containing the sum */
424:   PetscCall(ISCreateGeneral(comm, nout, iout, PETSC_OWN_POINTER, isout));

426:   PetscCall(ISRestoreIndices(is2sorted, &i2));
427:   PetscCall(ISDestroy(&is2sorted));
428:   PetscCall(ISRestoreIndices(is1sorted, &i1));
429:   PetscCall(ISDestroy(&is1sorted));
430:   PetscFunctionReturn(PETSC_SUCCESS);
431: }

433: PetscErrorCode ISIntersect_Caching_Internal(IS is1, IS is2, IS *isect)
434: {
435:   PetscFunctionBegin;
436:   *isect = NULL;
437:   if (is2 && is1) {
438:     char          composeStr[33] = {0};
439:     PetscObjectId is2id;

441:     PetscCall(PetscObjectGetId((PetscObject)is2, &is2id));
442:     PetscCall(PetscSNPrintf(composeStr, 32, "ISIntersect_Caching_%" PetscInt64_FMT, is2id));
443:     PetscCall(PetscObjectQuery((PetscObject)is1, composeStr, (PetscObject *)isect));
444:     if (*isect == NULL) {
445:       PetscCall(ISIntersect(is1, is2, isect));
446:       PetscCall(PetscObjectCompose((PetscObject)is1, composeStr, (PetscObject)*isect));
447:     } else {
448:       PetscCall(PetscObjectReference((PetscObject)*isect));
449:     }
450:   }
451:   PetscFunctionReturn(PETSC_SUCCESS);
452: }

454: /*@
455:    ISConcatenate - Forms a new `IS` by locally concatenating the indices from an `IS` list without reordering.

457:    Collective.

459:    Input Parameters:
460: +  comm    - communicator of the concatenated IS.
461: .  len     - size of islist array (nonnegative)
462: -  islist  - array of index sets

464:    Output Parameters:
465: .  isout   - The concatenated index set; empty, if len == 0.

467:    Level: intermediate

469:    Notes:
470:     The semantics of calling this on comm imply that the comms of the members of islist also contain this rank.

472: .seealso: [](sec_scatter), `IS`, `ISDifference()`, `ISSum()`, `ISExpand()`, `ISIntersect()`
473: @*/
474: PetscErrorCode ISConcatenate(MPI_Comm comm, PetscInt len, const IS islist[], IS *isout)
475: {
476:   PetscInt        i, n, N;
477:   const PetscInt *iidx;
478:   PetscInt       *idx;

480:   PetscFunctionBegin;
482:   if (PetscDefined(USE_DEBUG)) {
483:     for (i = 0; i < len; ++i)
485:   }
487:   if (!len) {
488:     PetscCall(ISCreateGeneral(comm, 0, NULL, PETSC_OWN_POINTER, isout));
489:     PetscFunctionReturn(PETSC_SUCCESS);
490:   }
491:   PetscCheck(len >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Negative array length: %" PetscInt_FMT, len);
492:   N = 0;
493:   for (i = 0; i < len; ++i) {
494:     if (islist[i]) {
495:       PetscCall(ISGetLocalSize(islist[i], &n));
496:       N += n;
497:     }
498:   }
499:   PetscCall(PetscMalloc1(N, &idx));
500:   N = 0;
501:   for (i = 0; i < len; ++i) {
502:     if (islist[i]) {
503:       PetscCall(ISGetLocalSize(islist[i], &n));
504:       PetscCall(ISGetIndices(islist[i], &iidx));
505:       PetscCall(PetscArraycpy(idx + N, iidx, n));
506:       PetscCall(ISRestoreIndices(islist[i], &iidx));
507:       N += n;
508:     }
509:   }
510:   PetscCall(ISCreateGeneral(comm, N, idx, PETSC_OWN_POINTER, isout));
511:   PetscFunctionReturn(PETSC_SUCCESS);
512: }

514: /*@
515:    ISListToPair  - Convert an `IS` list to a pair of `IS` of equal length defining an equivalent integer multimap.
516:                    Each `IS` on the input list is assigned an integer j so that all of the indices of that `IS` are
517:                    mapped to j.

519:    Collective.

521:    Input Parameters:
522: +  comm    -  `MPI_Comm`
523: .  listlen -  `IS` list length
524: -  islist  -  `IS` list

526:    Output Parameters:
527: +  xis -  domain `IS`
528: -  yis -  range  `I`S

530:   Level: developer

532:   Notes:
533:   The global integers assigned to the `IS` of the local input list might not correspond to the
534:   local numbers of the `IS` on that list, but the two *orderings* are the same: the global
535:   integers assigned to the `IS` on the local list form a strictly increasing sequence.

537:   The `IS` on the input list can belong to subcommunicators of comm, and the subcommunicators
538:   on the input `IS` list are assumed to be in a "deadlock-free" order.

540:   Local lists of `PetscObject`s (or their subcomms) on a comm are "deadlock-free" if subcomm1
541:   precedes subcomm2 on any local list, then it precedes subcomm2 on all ranks.
542:   Equivalently, the local numbers of the subcomms on each local list are drawn from some global
543:   numbering. This is ensured, for example, by `ISPairToList()`.

545: .seealso: `IS`, `ISPairToList()`
546: @*/
547: PetscErrorCode ISListToPair(MPI_Comm comm, PetscInt listlen, IS islist[], IS *xis, IS *yis)
548: {
549:   PetscInt        ncolors, *colors, i, leni, len, *xinds, *yinds, k, j;
550:   const PetscInt *indsi;

552:   PetscFunctionBegin;
553:   PetscCall(PetscMalloc1(listlen, &colors));
554:   PetscCall(PetscObjectsListGetGlobalNumbering(comm, listlen, (PetscObject *)islist, &ncolors, colors));
555:   len = 0;
556:   for (i = 0; i < listlen; ++i) {
557:     PetscCall(ISGetLocalSize(islist[i], &leni));
558:     len += leni;
559:   }
560:   PetscCall(PetscMalloc1(len, &xinds));
561:   PetscCall(PetscMalloc1(len, &yinds));
562:   k = 0;
563:   for (i = 0; i < listlen; ++i) {
564:     PetscCall(ISGetLocalSize(islist[i], &leni));
565:     PetscCall(ISGetIndices(islist[i], &indsi));
566:     for (j = 0; j < leni; ++j) {
567:       xinds[k] = indsi[j];
568:       yinds[k] = colors[i];
569:       ++k;
570:     }
571:   }
572:   PetscCall(PetscFree(colors));
573:   PetscCall(ISCreateGeneral(comm, len, xinds, PETSC_OWN_POINTER, xis));
574:   PetscCall(ISCreateGeneral(comm, len, yinds, PETSC_OWN_POINTER, yis));
575:   PetscFunctionReturn(PETSC_SUCCESS);
576: }

578: /*@
579:    ISPairToList - Convert an `IS` pair encoding an integer map to a list of `IS`.
580:                   Each `IS` on the output list contains the preimage for each index on the second input `IS`.
581:                   The `IS` on the output list are constructed on the subcommunicators of the input `IS` pair.
582:                   Each subcommunicator corresponds to the preimage of some index j -- this subcomm contains
583:                   exactly the ranks that assign some indices i to j.  This is essentially the inverse of
584:                   `ISListToPair()`.

586:   Collective

588:   Input Parameters:
589: + xis -  domain IS
590: - yis -  range IS

592:   Output Parameters:
593: + listlen -  length of islist
594: - islist  -  list of ISs breaking up indis by color

596:   Level: developer

598:   Notes:
599:   xis and yis must be of the same length and have congruent communicators.

601:   The resulting `IS` have subcommunicators in a "deadlock-free" order (see `ISListToPair()`).

603: .seealso: `IS`, `ISListToPair()`
604:  @*/
605: PetscErrorCode ISPairToList(IS xis, IS yis, PetscInt *listlen, IS **islist)
606: {
607:   IS              indis = xis, coloris = yis;
608:   PetscInt       *inds, *colors, llen, ilen, lstart, lend, lcount, l;
609:   PetscMPIInt     rank, size, llow, lhigh, low, high, color, subsize;
610:   const PetscInt *ccolors, *cinds;
611:   MPI_Comm        comm, subcomm;

613:   PetscFunctionBegin;
616:   PetscCheckSameComm(xis, 1, yis, 2);
619:   PetscCall(PetscObjectGetComm((PetscObject)xis, &comm));
620:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
621:   PetscCallMPI(MPI_Comm_rank(comm, &size));
622:   /* Extract, copy and sort the local indices and colors on the color. */
623:   PetscCall(ISGetLocalSize(coloris, &llen));
624:   PetscCall(ISGetLocalSize(indis, &ilen));
625:   PetscCheck(llen == ilen, comm, PETSC_ERR_ARG_SIZ, "Incompatible IS sizes: %" PetscInt_FMT " and %" PetscInt_FMT, ilen, llen);
626:   PetscCall(ISGetIndices(coloris, &ccolors));
627:   PetscCall(ISGetIndices(indis, &cinds));
628:   PetscCall(PetscMalloc2(ilen, &inds, llen, &colors));
629:   PetscCall(PetscArraycpy(inds, cinds, ilen));
630:   PetscCall(PetscArraycpy(colors, ccolors, llen));
631:   PetscCall(PetscSortIntWithArray(llen, colors, inds));
632:   /* Determine the global extent of colors. */
633:   llow   = 0;
634:   lhigh  = -1;
635:   lstart = 0;
636:   lcount = 0;
637:   while (lstart < llen) {
638:     lend = lstart + 1;
639:     while (lend < llen && colors[lend] == colors[lstart]) ++lend;
640:     llow  = PetscMin(llow, colors[lstart]);
641:     lhigh = PetscMax(lhigh, colors[lstart]);
642:     ++lcount;
643:   }
644:   PetscCall(MPIU_Allreduce(&llow, &low, 1, MPI_INT, MPI_MIN, comm));
645:   PetscCall(MPIU_Allreduce(&lhigh, &high, 1, MPI_INT, MPI_MAX, comm));
646:   *listlen = 0;
647:   if (low <= high) {
648:     if (lcount > 0) {
649:       *listlen = lcount;
650:       if (!*islist) PetscCall(PetscMalloc1(lcount, islist));
651:     }
652:     /*
653:      Traverse all possible global colors, and participate in the subcommunicators
654:      for the locally-supported colors.
655:      */
656:     lcount = 0;
657:     lstart = 0;
658:     lend   = 0;
659:     for (l = low; l <= high; ++l) {
660:       /*
661:        Find the range of indices with the same color, which is not smaller than l.
662:        Observe that, since colors is sorted, and is a subsequence of [low,high],
663:        as soon as we find a new color, it is >= l.
664:        */
665:       if (lstart < llen) {
666:         /* The start of the next locally-owned color is identified.  Now look for the end. */
667:         if (lstart == lend) {
668:           lend = lstart + 1;
669:           while (lend < llen && colors[lend] == colors[lstart]) ++lend;
670:         }
671:         /* Now check whether the identified color segment matches l. */
672:         PetscCheck(colors[lstart] >= l, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Locally owned color %" PetscInt_FMT " at location %" PetscInt_FMT " is < than the next global color %" PetscInt_FMT, colors[lstart], lcount, l);
673:       }
674:       color = (PetscMPIInt)(colors[lstart] == l);
675:       /* Check whether a proper subcommunicator exists. */
676:       PetscCall(MPIU_Allreduce(&color, &subsize, 1, MPI_INT, MPI_SUM, comm));

678:       if (subsize == 1) subcomm = PETSC_COMM_SELF;
679:       else if (subsize == size) subcomm = comm;
680:       else {
681:         /* a proper communicator is necessary, so we create it. */
682:         PetscCallMPI(MPI_Comm_split(comm, color, rank, &subcomm));
683:       }
684:       if (colors[lstart] == l) {
685:         /* If we have l among the local colors, we create an IS to hold the corresponding indices. */
686:         PetscCall(ISCreateGeneral(subcomm, lend - lstart, inds + lstart, PETSC_COPY_VALUES, *islist + lcount));
687:         /* Position lstart at the beginning of the next local color. */
688:         lstart = lend;
689:         /* Increment the counter of the local colors split off into an IS. */
690:         ++lcount;
691:       }
692:       if (subsize > 0 && subsize < size) {
693:         /*
694:          Irrespective of color, destroy the split off subcomm:
695:          a subcomm used in the IS creation above is duplicated
696:          into a proper PETSc comm.
697:          */
698:         PetscCallMPI(MPI_Comm_free(&subcomm));
699:       }
700:     } /* for (l = low; l < high; ++l) */
701:   }   /* if (low <= high) */
702:   PetscCall(PetscFree2(inds, colors));
703:   PetscFunctionReturn(PETSC_SUCCESS);
704: }

706: /*@
707:    ISEmbed - Embed `IS` a into `IS` b by finding the locations in b that have the same indices as in a.
708:              If c is the `IS` of these locations, we have a = b*c, regarded as a composition of the
709:              corresponding `ISLocalToGlobalMap`s.

711:   Not collective.

713:   Input Parameters:
714: + a    -  `IS` to embed
715: . b    -  `IS` to embed into
716: - drop -  flag indicating whether to drop a's indices that are not in b.

718:   Output Parameters:
719: . c    -  local embedding indices

721:   Level: developer

723:   Notes:
724:   If some of a's global indices are not among b's indices the embedding is impossible.  The local indices of a
725:   corresponding to these global indices are either mapped to -1 (if !drop) or are omitted (if drop).  In the former
726:   case the size of c is that same as that of a, in the latter case c's size may be smaller.

728:   The resulting `IS` is sequential, since the index substitution it encodes is purely local.

730: .seealso: `IS`, `ISLocalToGlobalMapping`
731:  @*/
732: PetscErrorCode ISEmbed(IS a, IS b, PetscBool drop, IS *c)
733: {
734:   ISLocalToGlobalMapping     ltog;
735:   ISGlobalToLocalMappingMode gtoltype = IS_GTOLM_DROP;
736:   PetscInt                   alen, clen, *cindices, *cindices2;
737:   const PetscInt            *aindices;

739:   PetscFunctionBegin;
743:   PetscCall(ISLocalToGlobalMappingCreateIS(b, &ltog));
744:   PetscCall(ISGetLocalSize(a, &alen));
745:   PetscCall(ISGetIndices(a, &aindices));
746:   PetscCall(PetscMalloc1(alen, &cindices));
747:   if (!drop) gtoltype = IS_GTOLM_MASK;
748:   PetscCall(ISGlobalToLocalMappingApply(ltog, gtoltype, alen, aindices, &clen, cindices));
749:   PetscCall(ISRestoreIndices(a, &aindices));
750:   PetscCall(ISLocalToGlobalMappingDestroy(&ltog));
751:   if (clen != alen) {
752:     cindices2 = cindices;
753:     PetscCall(PetscMalloc1(clen, &cindices));
754:     PetscCall(PetscArraycpy(cindices, cindices2, clen));
755:     PetscCall(PetscFree(cindices2));
756:   }
757:   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, clen, cindices, PETSC_OWN_POINTER, c));
758:   PetscFunctionReturn(PETSC_SUCCESS);
759: }

761: /*@
762:   ISSortPermutation  -  calculate the permutation of the indices into a nondecreasing order.

764:   Not collective.

766:   Input Parameters:
767: + f      -  `IS` to sort
768: - always -  build the permutation even when f's indices are nondecreasing.

770:   Output Parameter:
771: . h    -  permutation or NULL, if f is nondecreasing and always == `PETSC_FALSE`.

773:   Level: advanced

775:   Notes:
776:   Indices in f are unchanged. f[h[i]] is the i-th smallest f index.

778:   If always == `PETSC_FALSE`, an extra check is performed to see whether
779:   the f indices are nondecreasing. h is built on `PETSC_COMM_SELF`, since
780:   the permutation has a local meaning only.

782: .seealso: `IS`, `ISLocalToGlobalMapping`, `ISSort()`
783:  @*/
784: PetscErrorCode ISSortPermutation(IS f, PetscBool always, IS *h)
785: {
786:   const PetscInt *findices;
787:   PetscInt        fsize, *hindices, i;
788:   PetscBool       isincreasing;

790:   PetscFunctionBegin;
793:   PetscCall(ISGetLocalSize(f, &fsize));
794:   PetscCall(ISGetIndices(f, &findices));
795:   *h = NULL;
796:   if (!always) {
797:     isincreasing = PETSC_TRUE;
798:     for (i = 1; i < fsize; ++i) {
799:       if (findices[i] <= findices[i - 1]) {
800:         isincreasing = PETSC_FALSE;
801:         break;
802:       }
803:     }
804:     if (isincreasing) {
805:       PetscCall(ISRestoreIndices(f, &findices));
806:       PetscFunctionReturn(PETSC_SUCCESS);
807:     }
808:   }
809:   PetscCall(PetscMalloc1(fsize, &hindices));
810:   for (i = 0; i < fsize; ++i) hindices[i] = i;
811:   PetscCall(PetscSortIntWithPermutation(fsize, findices, hindices));
812:   PetscCall(ISRestoreIndices(f, &findices));
813:   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, fsize, hindices, PETSC_OWN_POINTER, h));
814:   PetscCall(ISSetInfo(*h, IS_PERMUTATION, IS_LOCAL, PETSC_FALSE, PETSC_TRUE));
815:   PetscFunctionReturn(PETSC_SUCCESS);
816: }