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 IS

 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:    Notes:
 19:    Negative values are removed from the lists. is2 may have values
 20:    that are not in is1. This requires O(imax-imin) memory and O(imax-imin)
 21:    work, where imin and imax are the bounds on the indices in is1.

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

 25:    Level: intermediate

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

 38:   if (!is2) {
 39:     ISDuplicate(is1, isout);
 40:     return 0;
 41:   }

 44:   ISGetIndices(is1, &i1);
 45:   ISGetLocalSize(is1, &n1);

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

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

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

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

 89:   PetscBTDestroy(&mask);
 90:   return 0;
 91: }

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

 96:    Only sequential version (at the moment)

 98:    Input Parameters:
 99: +  is1 - index set to be extended
100: -  is2 - index values to be added

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

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

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

110:    Level: intermediate

112: .seealso: `ISDestroy()`, `ISView()`, `ISDifference()`, `ISExpand()`

114: @*/
115: PetscErrorCode ISSum(IS is1, IS is2, IS *is3)
116: {
117:   MPI_Comm        comm;
118:   PetscBool       f;
119:   PetscMPIInt     size;
120:   const PetscInt *i1, *i2;
121:   PetscInt        n1, n2, n3, p1, p2, *iout;

125:   PetscObjectGetComm((PetscObject)(is1), &comm);
126:   MPI_Comm_size(comm, &size);

129:   ISSorted(is1, &f);
131:   ISSorted(is2, &f);

134:   ISGetLocalSize(is1, &n1);
135:   ISGetLocalSize(is2, &n2);
136:   if (!n2) {
137:     ISDuplicate(is1, is3);
138:     return 0;
139:   }
140:   ISGetIndices(is1, &i1);
141:   ISGetIndices(is2, &i2);

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

190:   if (n3 == n1) { /* no new elements to be added */
191:     ISRestoreIndices(is1, &i1);
192:     ISRestoreIndices(is2, &i2);
193:     ISDuplicate(is1, is3);
194:     return 0;
195:   }
196:   PetscMalloc1(n3, &iout);

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

234:   ISRestoreIndices(is1, &i1);
235:   ISRestoreIndices(is2, &i2);
236:   ISCreateGeneral(comm, n3, iout, PETSC_OWN_POINTER, is3);
237:   return 0;
238: }

240: /*@
241:    ISExpand - Computes the union of two index sets, by concatenating 2 lists and
242:    removing duplicates.

244:    Collective on IS

246:    Input Parameters:
247: +  is1 - first index set
248: -  is2 - index values to be added

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

253:    Notes:
254:    Negative values are removed from the lists. This requires O(imax-imin)
255:    memory and O(imax-imin) work, where imin and imax are the bounds on the
256:    indices in is1 and is2.

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

260:    Level: intermediate

262: .seealso: `ISDestroy()`, `ISView()`, `ISDifference()`, `ISSum()`

264: @*/
265: PetscErrorCode ISExpand(IS is1, IS is2, IS *isout)
266: {
267:   PetscInt        i, n1, n2, imin, imax, nout, *iout;
268:   const PetscInt *i1, *i2;
269:   PetscBT         mask;
270:   MPI_Comm        comm;


277:   if (!is1) {
278:     ISDuplicate(is2, isout);
279:     return 0;
280:   }
281:   if (!is2) {
282:     ISDuplicate(is1, isout);
283:     return 0;
284:   }
285:   ISGetIndices(is1, &i1);
286:   ISGetLocalSize(is1, &n1);
287:   ISGetIndices(is2, &i2);
288:   ISGetLocalSize(is2, &n2);

290:   /* Create a bit mask array to contain required values */
291:   if (n1 || n2) {
292:     imin = PETSC_MAX_INT;
293:     imax = 0;
294:     for (i = 0; i < n1; i++) {
295:       if (i1[i] < 0) continue;
296:       imin = PetscMin(imin, i1[i]);
297:       imax = PetscMax(imax, i1[i]);
298:     }
299:     for (i = 0; i < n2; i++) {
300:       if (i2[i] < 0) continue;
301:       imin = PetscMin(imin, i2[i]);
302:       imax = PetscMax(imax, i2[i]);
303:     }
304:   } else imin = imax = 0;

306:   PetscMalloc1(n1 + n2, &iout);
307:   nout = 0;
308:   PetscBTCreate(imax - imin, &mask);
309:   /* Put the values from is1 */
310:   for (i = 0; i < n1; i++) {
311:     if (i1[i] < 0) continue;
312:     if (!PetscBTLookupSet(mask, i1[i] - imin)) iout[nout++] = i1[i];
313:   }
314:   ISRestoreIndices(is1, &i1);
315:   /* Put the values from is2 */
316:   for (i = 0; i < n2; i++) {
317:     if (i2[i] < 0) continue;
318:     if (!PetscBTLookupSet(mask, i2[i] - imin)) iout[nout++] = i2[i];
319:   }
320:   ISRestoreIndices(is2, &i2);

322:   /* create the new IS containing the sum */
323:   PetscObjectGetComm((PetscObject)is1, &comm);
324:   ISCreateGeneral(comm, nout, iout, PETSC_OWN_POINTER, isout);

326:   PetscBTDestroy(&mask);
327:   return 0;
328: }

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

333:    Collective on IS

335:    Input Parameters:
336: +  is1 - first index set
337: -  is2 - second index set

339:    Output Parameters:
340: .  isout - the sorted intersection of is1 and is2

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

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

348:    Level: intermediate

350: .seealso: `ISDestroy()`, `ISView()`, `ISDifference()`, `ISSum()`, `ISExpand()`
351: @*/
352: PetscErrorCode ISIntersect(IS is1, IS is2, IS *isout)
353: {
354:   PetscInt        i, n1, n2, nout, *iout;
355:   const PetscInt *i1, *i2;
356:   IS              is1sorted = NULL, is2sorted = NULL;
357:   PetscBool       sorted, lsorted;
358:   MPI_Comm        comm;

364:   PetscObjectGetComm((PetscObject)is1, &comm);

366:   ISGetLocalSize(is1, &n1);
367:   ISGetLocalSize(is2, &n2);
368:   if (n1 < n2) {
369:     IS       tempis = is1;
370:     PetscInt ntemp  = n1;

372:     is1 = is2;
373:     is2 = tempis;
374:     n1  = n2;
375:     n2  = ntemp;
376:   }
377:   ISSorted(is1, &lsorted);
378:   MPIU_Allreduce(&lsorted, &sorted, 1, MPIU_BOOL, MPI_LAND, comm);
379:   if (!sorted) {
380:     ISDuplicate(is1, &is1sorted);
381:     ISSort(is1sorted);
382:     ISGetIndices(is1sorted, &i1);
383:   } else {
384:     is1sorted = is1;
385:     PetscObjectReference((PetscObject)is1);
386:     ISGetIndices(is1, &i1);
387:   }
388:   ISSorted(is2, &lsorted);
389:   MPIU_Allreduce(&lsorted, &sorted, 1, MPIU_BOOL, MPI_LAND, comm);
390:   if (!sorted) {
391:     ISDuplicate(is2, &is2sorted);
392:     ISSort(is2sorted);
393:     ISGetIndices(is2sorted, &i2);
394:   } else {
395:     is2sorted = is2;
396:     PetscObjectReference((PetscObject)is2);
397:     ISGetIndices(is2, &i2);
398:   }

400:   PetscMalloc1(n2, &iout);

402:   for (i = 0, nout = 0; i < n2; i++) {
403:     PetscInt key = i2[i];
404:     PetscInt loc;

406:     ISLocate(is1sorted, key, &loc);
407:     if (loc >= 0) {
408:       if (!nout || iout[nout - 1] < key) iout[nout++] = key;
409:     }
410:   }
411:   PetscRealloc(nout * sizeof(PetscInt), &iout);

413:   /* create the new IS containing the sum */
414:   ISCreateGeneral(comm, nout, iout, PETSC_OWN_POINTER, isout);

416:   ISRestoreIndices(is2sorted, &i2);
417:   ISDestroy(&is2sorted);
418:   ISRestoreIndices(is1sorted, &i1);
419:   ISDestroy(&is1sorted);
420:   return 0;
421: }

423: PetscErrorCode ISIntersect_Caching_Internal(IS is1, IS is2, IS *isect)
424: {
425:   *isect = NULL;
426:   if (is2 && is1) {
427:     char          composeStr[33] = {0};
428:     PetscObjectId is2id;

430:     PetscObjectGetId((PetscObject)is2, &is2id);
431:     PetscSNPrintf(composeStr, 32, "ISIntersect_Caching_%" PetscInt64_FMT, is2id);
432:     PetscObjectQuery((PetscObject)is1, composeStr, (PetscObject *)isect);
433:     if (*isect == NULL) {
434:       ISIntersect(is1, is2, isect);
435:       PetscObjectCompose((PetscObject)is1, composeStr, (PetscObject)*isect);
436:     } else {
437:       PetscObjectReference((PetscObject)*isect);
438:     }
439:   }
440:   return 0;
441: }

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

446:    Collective.

448:    Input Parameters:
449: +  comm    - communicator of the concatenated IS.
450: .  len     - size of islist array (nonnegative)
451: -  islist  - array of index sets

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

456:    Notes:
457:     The semantics of calling this on comm imply that the comms of the members if islist also contain this rank.

459:    Level: intermediate

461: .seealso: `ISDifference()`, `ISSum()`, `ISExpand()`

463: @*/
464: PetscErrorCode ISConcatenate(MPI_Comm comm, PetscInt len, const IS islist[], IS *isout)
465: {
466:   PetscInt        i, n, N;
467:   const PetscInt *iidx;
468:   PetscInt       *idx;

471:   if (PetscDefined(USE_DEBUG)) {
472:     for (i = 0; i < len; ++i)
474:   }
476:   if (!len) {
477:     ISCreateGeneral(comm, 0, NULL, PETSC_OWN_POINTER, isout);
478:     return 0;
479:   }
481:   N = 0;
482:   for (i = 0; i < len; ++i) {
483:     if (islist[i]) {
484:       ISGetLocalSize(islist[i], &n);
485:       N += n;
486:     }
487:   }
488:   PetscMalloc1(N, &idx);
489:   N = 0;
490:   for (i = 0; i < len; ++i) {
491:     if (islist[i]) {
492:       ISGetLocalSize(islist[i], &n);
493:       ISGetIndices(islist[i], &iidx);
494:       PetscArraycpy(idx + N, iidx, n);
495:       ISRestoreIndices(islist[i], &iidx);
496:       N += n;
497:     }
498:   }
499:   ISCreateGeneral(comm, N, idx, PETSC_OWN_POINTER, isout);
500:   return 0;
501: }

503: /*@
504:    ISListToPair     -    convert an IS list to a pair of ISs of equal length defining an equivalent integer multimap.
505:                         Each IS on the input list is assigned an integer j so that all of the indices of that IS are
506:                         mapped to j.

508:   Collective.

510:   Input arguments:
511: + comm    -  MPI_Comm
512: . listlen -  IS list length
513: - islist  -  IS list

515:   Output arguments:
516: + xis -  domain IS
517: - yis -  range  IS

519:   Level: advanced

521:   Notes:
522:   The global integers assigned to the ISs of the local input list might not correspond to the
523:   local numbers of the ISs on that list, but the two *orderings* are the same: the global
524:   integers assigned to the ISs on the local list form a strictly increasing sequence.

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

529:   Local lists of PetscObjects (or their subcommes) on a comm are "deadlock-free" if subcomm1
530:   preceeds subcomm2 on any local list, then it preceeds subcomm2 on all ranks.
531:   Equivalently, the local numbers of the subcomms on each local list are drawn from some global
532:   numbering. This is ensured, for example, by ISPairToList().

534: .seealso `ISPairToList()`
535: @*/
536: PetscErrorCode ISListToPair(MPI_Comm comm, PetscInt listlen, IS islist[], IS *xis, IS *yis)
537: {
538:   PetscInt        ncolors, *colors, i, leni, len, *xinds, *yinds, k, j;
539:   const PetscInt *indsi;

541:   PetscMalloc1(listlen, &colors);
542:   PetscObjectsListGetGlobalNumbering(comm, listlen, (PetscObject *)islist, &ncolors, colors);
543:   len = 0;
544:   for (i = 0; i < listlen; ++i) {
545:     ISGetLocalSize(islist[i], &leni);
546:     len += leni;
547:   }
548:   PetscMalloc1(len, &xinds);
549:   PetscMalloc1(len, &yinds);
550:   k = 0;
551:   for (i = 0; i < listlen; ++i) {
552:     ISGetLocalSize(islist[i], &leni);
553:     ISGetIndices(islist[i], &indsi);
554:     for (j = 0; j < leni; ++j) {
555:       xinds[k] = indsi[j];
556:       yinds[k] = colors[i];
557:       ++k;
558:     }
559:   }
560:   PetscFree(colors);
561:   ISCreateGeneral(comm, len, xinds, PETSC_OWN_POINTER, xis);
562:   ISCreateGeneral(comm, len, yinds, PETSC_OWN_POINTER, yis);
563:   return 0;
564: }

566: /*@
567:    ISPairToList   -   convert an IS pair encoding an integer map to a list of ISs.
568:                      Each IS on the output list contains the preimage for each index on the second input IS.
569:                      The ISs on the output list are constructed on the subcommunicators of the input IS pair.
570:                      Each subcommunicator corresponds to the preimage of some index j -- this subcomm contains
571:                      exactly the ranks that assign some indices i to j.  This is essentially the inverse of
572:                      ISListToPair().

574:   Collective on indis.

576:   Input arguments:
577: + xis -  domain IS
578: - yis -  range IS

580:   Output arguments:
581: + listlen -  length of islist
582: - islist  -  list of ISs breaking up indis by color

584:   Note:
585:     xis and yis must be of the same length and have congruent communicators.

587:     The resulting ISs have subcommunicators in a "deadlock-free" order (see ISListToPair()).

589:   Level: advanced

591: .seealso `ISListToPair()`
592:  @*/
593: PetscErrorCode ISPairToList(IS xis, IS yis, PetscInt *listlen, IS **islist)
594: {
595:   IS              indis = xis, coloris = yis;
596:   PetscInt       *inds, *colors, llen, ilen, lstart, lend, lcount, l;
597:   PetscMPIInt     rank, size, llow, lhigh, low, high, color, subsize;
598:   const PetscInt *ccolors, *cinds;
599:   MPI_Comm        comm, subcomm;

606:   PetscObjectGetComm((PetscObject)xis, &comm);
607:   MPI_Comm_rank(comm, &rank);
608:   MPI_Comm_rank(comm, &size);
609:   /* Extract, copy and sort the local indices and colors on the color. */
610:   ISGetLocalSize(coloris, &llen);
611:   ISGetLocalSize(indis, &ilen);
613:   ISGetIndices(coloris, &ccolors);
614:   ISGetIndices(indis, &cinds);
615:   PetscMalloc2(ilen, &inds, llen, &colors);
616:   PetscArraycpy(inds, cinds, ilen);
617:   PetscArraycpy(colors, ccolors, llen);
618:   PetscSortIntWithArray(llen, colors, inds);
619:   /* Determine the global extent of colors. */
620:   llow   = 0;
621:   lhigh  = -1;
622:   lstart = 0;
623:   lcount = 0;
624:   while (lstart < llen) {
625:     lend = lstart + 1;
626:     while (lend < llen && colors[lend] == colors[lstart]) ++lend;
627:     llow  = PetscMin(llow, colors[lstart]);
628:     lhigh = PetscMax(lhigh, colors[lstart]);
629:     ++lcount;
630:   }
631:   MPIU_Allreduce(&llow, &low, 1, MPI_INT, MPI_MIN, comm);
632:   MPIU_Allreduce(&lhigh, &high, 1, MPI_INT, MPI_MAX, comm);
633:   *listlen = 0;
634:   if (low <= high) {
635:     if (lcount > 0) {
636:       *listlen = lcount;
637:       if (!*islist) PetscMalloc1(lcount, islist);
638:     }
639:     /*
640:      Traverse all possible global colors, and participate in the subcommunicators
641:      for the locally-supported colors.
642:      */
643:     lcount = 0;
644:     lstart = 0;
645:     lend   = 0;
646:     for (l = low; l <= high; ++l) {
647:       /*
648:        Find the range of indices with the same color, which is not smaller than l.
649:        Observe that, since colors is sorted, and is a subsequence of [low,high],
650:        as soon as we find a new color, it is >= l.
651:        */
652:       if (lstart < llen) {
653:         /* The start of the next locally-owned color is identified.  Now look for the end. */
654:         if (lstart == lend) {
655:           lend = lstart + 1;
656:           while (lend < llen && colors[lend] == colors[lstart]) ++lend;
657:         }
658:         /* Now check whether the identified color segment matches l. */
660:       }
661:       color = (PetscMPIInt)(colors[lstart] == l);
662:       /* Check whether a proper subcommunicator exists. */
663:       MPIU_Allreduce(&color, &subsize, 1, MPI_INT, MPI_SUM, comm);

665:       if (subsize == 1) subcomm = PETSC_COMM_SELF;
666:       else if (subsize == size) subcomm = comm;
667:       else {
668:         /* a proper communicator is necessary, so we create it. */
669:         MPI_Comm_split(comm, color, rank, &subcomm);
670:       }
671:       if (colors[lstart] == l) {
672:         /* If we have l among the local colors, we create an IS to hold the corresponding indices. */
673:         ISCreateGeneral(subcomm, lend - lstart, inds + lstart, PETSC_COPY_VALUES, *islist + lcount);
674:         /* Position lstart at the beginning of the next local color. */
675:         lstart = lend;
676:         /* Increment the counter of the local colors split off into an IS. */
677:         ++lcount;
678:       }
679:       if (subsize > 0 && subsize < size) {
680:         /*
681:          Irrespective of color, destroy the split off subcomm:
682:          a subcomm used in the IS creation above is duplicated
683:          into a proper PETSc comm.
684:          */
685:         MPI_Comm_free(&subcomm);
686:       }
687:     } /* for (l = low; l < high; ++l) */
688:   }   /* if (low <= high) */
689:   PetscFree2(inds, colors);
690:   return 0;
691: }

693: /*@
694:    ISEmbed   -   embed IS a into IS b by finding the locations in b that have the same indices as in a.
695:                  If c is the IS of these locations, we have a = b*c, regarded as a composition of the
696:                  corresponding ISLocalToGlobalMaps.

698:   Not collective.

700:   Input arguments:
701: + a    -  IS to embed
702: . b    -  IS to embed into
703: - drop -  flag indicating whether to drop a's indices that are not in b.

705:   Output arguments:
706: . c    -  local embedding indices

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

713:   The resulting IS is sequential, since the index substition it encodes is purely local.

715:   Level: advanced

717: .seealso `ISLocalToGlobalMapping`
718:  @*/
719: PetscErrorCode ISEmbed(IS a, IS b, PetscBool drop, IS *c)
720: {
721:   ISLocalToGlobalMapping     ltog;
722:   ISGlobalToLocalMappingMode gtoltype = IS_GTOLM_DROP;
723:   PetscInt                   alen, clen, *cindices, *cindices2;
724:   const PetscInt            *aindices;

729:   ISLocalToGlobalMappingCreateIS(b, &ltog);
730:   ISGetLocalSize(a, &alen);
731:   ISGetIndices(a, &aindices);
732:   PetscMalloc1(alen, &cindices);
733:   if (!drop) gtoltype = IS_GTOLM_MASK;
734:   ISGlobalToLocalMappingApply(ltog, gtoltype, alen, aindices, &clen, cindices);
735:   ISRestoreIndices(a, &aindices);
736:   ISLocalToGlobalMappingDestroy(&ltog);
737:   if (clen != alen) {
738:     cindices2 = cindices;
739:     PetscMalloc1(clen, &cindices);
740:     PetscArraycpy(cindices, cindices2, clen);
741:     PetscFree(cindices2);
742:   }
743:   ISCreateGeneral(PETSC_COMM_SELF, clen, cindices, PETSC_OWN_POINTER, c);
744:   return 0;
745: }

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

750:   Not collective.

752:   Input arguments:
753: + f      -  IS to sort
754: - always -  build the permutation even when f's indices are nondecreasing.

756:   Output argument:
757: . h    -  permutation or NULL, if f is nondecreasing and always == PETSC_FALSE.

759:   Note: Indices in f are unchanged. f[h[i]] is the i-th smallest f index.
760:         If always == PETSC_FALSE, an extra check is peformed to see whether
761:         the f indices are nondecreasing. h is built on PETSC_COMM_SELF, since
762:         the permutation has a local meaning only.

764:   Level: advanced

766: .seealso `ISLocalToGlobalMapping`, `ISSort()`
767:  @*/
768: PetscErrorCode ISSortPermutation(IS f, PetscBool always, IS *h)
769: {
770:   const PetscInt *findices;
771:   PetscInt        fsize, *hindices, i;
772:   PetscBool       isincreasing;

776:   ISGetLocalSize(f, &fsize);
777:   ISGetIndices(f, &findices);
778:   *h = NULL;
779:   if (!always) {
780:     isincreasing = PETSC_TRUE;
781:     for (i = 1; i < fsize; ++i) {
782:       if (findices[i] <= findices[i - 1]) {
783:         isincreasing = PETSC_FALSE;
784:         break;
785:       }
786:     }
787:     if (isincreasing) {
788:       ISRestoreIndices(f, &findices);
789:       return 0;
790:     }
791:   }
792:   PetscMalloc1(fsize, &hindices);
793:   for (i = 0; i < fsize; ++i) hindices[i] = i;
794:   PetscSortIntWithPermutation(fsize, findices, hindices);
795:   ISRestoreIndices(f, &findices);
796:   ISCreateGeneral(PETSC_COMM_SELF, fsize, hindices, PETSC_OWN_POINTER, h);
797:   ISSetInfo(*h, IS_PERMUTATION, IS_LOCAL, PETSC_FALSE, PETSC_TRUE);
798:   return 0;
799: }