Actual source code: isdiff.c

petsc-3.12.0 2019-09-29
Report Typos and Errors

  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 Parameter:
 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:    Level: intermediate


 26: .seealso: ISDestroy(), ISView(), ISSum(), ISExpand()

 28: @*/
 29: PetscErrorCode  ISDifference(IS is1,IS is2,IS *isout)
 30: {
 32:   PetscInt       i,n1,n2,imin,imax,nout,*iout;
 33:   const PetscInt *i1,*i2;
 34:   PetscBT        mask;
 35:   MPI_Comm       comm;


 42:   ISGetIndices(is1,&i1);
 43:   ISGetLocalSize(is1,&n1);

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

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

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

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

 87:   PetscBTDestroy(&mask);
 88:   return(0);
 89: }

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

 94:    Only sequential version (at the moment)

 96:    Input Parameter:
 97: +  is1 - index set to be extended
 98: -  is2 - index values to be added

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

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

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

108:    Level: intermediate

110: .seealso: ISDestroy(), ISView(), ISDifference(), ISExpand()


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

126:   PetscObjectGetComm((PetscObject)(is1),&comm);
127:   MPI_Comm_size(comm,&size);
128:   if (size>1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently only for uni-processor IS");

130:   ISSorted(is1,&f);
131:   if (!f) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Arg 1 is not sorted");
132:   ISSorted(is2,&f);
133:   if (!f) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Arg 2 is not sorted");

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

144:   p1 = 0; p2 = 0; n3 = 0;
145:   do {
146:     if (p1==n1) { /* cleanup for is2 */ n3 += n2-p2; break;
147:     } else {
148:       while (p2<n2 && i2[p2]<i1[p1]) {
149:         n3++; p2++;
150:       }
151:       if (p2==n2) {
152:         /* cleanup for is1 */
153:         n3 += n1-p1; break;
154:       } else {
155:         if (i2[p2]==i1[p1]) { n3++; p1++; p2++; }
156:       }
157:     }
158:     if (p2==n2) {
159:       /* cleanup for is1 */
160:       n3 += n1-p1; break;
161:     } else {
162:       while (p1<n1 && i1[p1]<i2[p2]) {
163:         n3++; p1++;
164:       }
165:       if (p1==n1) {
166:         /* clean up for is2 */
167:         n3 += n2-p2; break;
168:       } else {
169:         if (i1[p1]==i2[p2]) { n3++; p1++; p2++; }
170:       }
171:     }
172:   } while (p1<n1 || p2<n2);

174:   if (n3==n1) { /* no new elements to be added */
175:     ISRestoreIndices(is1,&i1);
176:     ISRestoreIndices(is2,&i2);
177:     ISDuplicate(is1,is3);
178:     return(0);
179:   }
180:   PetscMalloc1(n3,&iout);

182:   p1 = 0; p2 = 0; n3 = 0;
183:   do {
184:     if (p1==n1) { /* cleanup for is2 */
185:       while (p2<n2) iout[n3++] = i2[p2++];
186:       break;
187:     } else {
188:       while (p2<n2 && i2[p2]<i1[p1]) iout[n3++] = i2[p2++];
189:       if (p2==n2) { /* cleanup for is1 */
190:         while (p1<n1) iout[n3++] = i1[p1++];
191:         break;
192:       } else {
193:         if (i2[p2]==i1[p1]) { iout[n3++] = i1[p1++]; p2++; }
194:       }
195:     }
196:     if (p2==n2) { /* cleanup for is1 */
197:       while (p1<n1) iout[n3++] = i1[p1++];
198:       break;
199:     } else {
200:       while (p1<n1 && i1[p1]<i2[p2]) iout[n3++] = i1[p1++];
201:       if (p1==n1) { /* clean up for is2 */
202:         while (p2<n2) iout[n3++] = i2[p2++];
203:         break;
204:       } else {
205:         if (i1[p1]==i2[p2]) { iout[n3++] = i1[p1++]; p2++; }
206:       }
207:     }
208:   } while (p1<n1 || p2<n2);

210:   ISRestoreIndices(is1,&i1);
211:   ISRestoreIndices(is2,&i2);
212:   ISCreateGeneral(comm,n3,iout,PETSC_OWN_POINTER,is3);
213:   return(0);
214: }

216: /*@
217:    ISExpand - Computes the union of two index sets, by concatenating 2 lists and
218:    removing duplicates.

220:    Collective on IS

222:    Input Parameter:
223: +  is1 - first index set
224: -  is2 - index values to be added

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

229:    Notes:
230:    Negative values are removed from the lists. This requires O(imax-imin)
231:    memory and O(imax-imin) work, where imin and imax are the bounds on the
232:    indices in is1 and is2.

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

236:    Level: intermediate

238: .seealso: ISDestroy(), ISView(), ISDifference(), ISSum()


241: @*/
242: PetscErrorCode ISExpand(IS is1,IS is2,IS *isout)
243: {
245:   PetscInt       i,n1,n2,imin,imax,nout,*iout;
246:   const PetscInt *i1,*i2;
247:   PetscBT        mask;
248:   MPI_Comm       comm;


255:   if (!is1 && !is2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Both arguments cannot be NULL");
256:   if (!is1) {ISDuplicate(is2, isout);return(0);}
257:   if (!is2) {ISDuplicate(is1, isout);return(0);}
258:   ISGetIndices(is1,&i1);
259:   ISGetLocalSize(is1,&n1);
260:   ISGetIndices(is2,&i2);
261:   ISGetLocalSize(is2,&n2);

263:   /* Create a bit mask array to contain required values */
264:   if (n1 || n2) {
265:     imin = PETSC_MAX_INT;
266:     imax = 0;
267:     for (i=0; i<n1; i++) {
268:       if (i1[i] < 0) continue;
269:       imin = PetscMin(imin,i1[i]);
270:       imax = PetscMax(imax,i1[i]);
271:     }
272:     for (i=0; i<n2; i++) {
273:       if (i2[i] < 0) continue;
274:       imin = PetscMin(imin,i2[i]);
275:       imax = PetscMax(imax,i2[i]);
276:     }
277:   } else imin = imax = 0;

279:   PetscMalloc1(n1+n2,&iout);
280:   nout = 0;
281:   PetscBTCreate(imax-imin,&mask);
282:   /* Put the values from is1 */
283:   for (i=0; i<n1; i++) {
284:     if (i1[i] < 0) continue;
285:     if (!PetscBTLookupSet(mask,i1[i] - imin)) iout[nout++] = i1[i];
286:   }
287:   ISRestoreIndices(is1,&i1);
288:   /* Put the values from is2 */
289:   for (i=0; i<n2; i++) {
290:     if (i2[i] < 0) continue;
291:     if (!PetscBTLookupSet(mask,i2[i] - imin)) iout[nout++] = i2[i];
292:   }
293:   ISRestoreIndices(is2,&i2);

295:   /* create the new IS containing the sum */
296:   PetscObjectGetComm((PetscObject)is1,&comm);
297:   ISCreateGeneral(comm,nout,iout,PETSC_OWN_POINTER,isout);

299:   PetscBTDestroy(&mask);
300:   return(0);
301: }

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

306:    Collective on IS

308:    Input Parameter:
309: +  is1 - first index set
310: -  is2 - second index set

312:    Output Parameters:
313: .  isout - the sorted intersection of is1 and is2

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

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

321:    Level: intermediate

323: .seealso: ISDestroy(), ISView(), ISDifference(), ISSum(), ISExpand()


326: @*/
327: PetscErrorCode ISIntersect(IS is1,IS is2,IS *isout)
328: {
330:   PetscInt       i,n1,n2,nout,*iout;
331:   const PetscInt *i1,*i2;
332:   IS             is1sorted = NULL, is2sorted = NULL;
333:   PetscBool      sorted, lsorted;
334:   MPI_Comm       comm;

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

343:   ISGetLocalSize(is1,&n1);
344:   ISGetLocalSize(is2,&n2);
345:   if (n1 < n2) {
346:     IS       tempis = is1;
347:     PetscInt ntemp = n1;

349:     is1 = is2;
350:     is2 = tempis;
351:     n1  = n2;
352:     n2  = ntemp;
353:   }
354:   ISSorted(is1,&lsorted);
355:   MPIU_Allreduce(&lsorted,&sorted,1,MPIU_BOOL,MPI_LAND,comm);
356:   if (!sorted) {
357:     ISDuplicate(is1,&is1sorted);
358:     ISSort(is1sorted);
359:     ISGetIndices(is1sorted,&i1);
360:   } else {
361:     is1sorted = is1;
362:     PetscObjectReference((PetscObject)is1);
363:     ISGetIndices(is1,&i1);
364:   }
365:   ISSorted(is2,&lsorted);
366:   MPIU_Allreduce(&lsorted,&sorted,1,MPIU_BOOL,MPI_LAND,comm);
367:   if (!sorted) {
368:     ISDuplicate(is2,&is2sorted);
369:     ISSort(is2sorted);
370:     ISGetIndices(is2sorted,&i2);
371:   } else {
372:     is2sorted = is2;
373:     PetscObjectReference((PetscObject)is2);
374:     ISGetIndices(is2,&i2);
375:   }

377:   PetscMalloc1(n2,&iout);

379:   for (i = 0, nout = 0; i < n2; i++) {
380:     PetscInt key = i2[i];
381:     PetscInt loc;

383:     ISLocate(is1sorted,key,&loc);
384:     if (loc >= 0) {
385:       if (!nout || iout[nout-1] < key) {
386:         iout[nout++] = key;
387:       }
388:     }
389:   }
390:   PetscRealloc(nout*sizeof(PetscInt),&iout);

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

395:   ISRestoreIndices(is2sorted,&i2);
396:   ISDestroy(&is2sorted);
397:   ISRestoreIndices(is1sorted,&i1);
398:   ISDestroy(&is1sorted);
399:   return(0);
400: }

402: PetscErrorCode ISIntersect_Caching_Internal(IS is1, IS is2, IS *isect)
403: {

407:   *isect = NULL;
408:   if (is2 && is1) {
409:     char           composeStr[33] = {0};
410:     PetscObjectId  is2id;

412:     PetscObjectGetId((PetscObject)is2,&is2id);
413:     PetscSNPrintf(composeStr,32,"ISIntersect_Caching_%x",is2id);
414:     PetscObjectQuery((PetscObject) is1, composeStr, (PetscObject *) isect);
415:     if (*isect == NULL) {
416:       ISIntersect(is1, is2, isect);
417:       PetscObjectCompose((PetscObject) is1, composeStr, (PetscObject) *isect);
418:     } else {
419:       PetscObjectReference((PetscObject) *isect);
420:     }
421:   }
422:   return(0);
423: }

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


429:    Collective.

431:    Input Parameter:
432: +  comm    - communicator of the concatenated IS.
433: .  len     - size of islist array (nonnegative)
434: -  islist  - array of index sets

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

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

442:    Level: intermediate

444: .seealso: ISDifference(), ISSum(), ISExpand()


447: @*/
448: PetscErrorCode ISConcatenate(MPI_Comm comm, PetscInt len, const IS islist[], IS *isout)
449: {
451:   PetscInt i,n,N;
452:   const PetscInt *iidx;
453:   PetscInt *idx;

457: #if defined(PETSC_USE_DEBUG)
459: #endif
461:   if (!len) {
462:     ISCreateStride(comm, 0,0,0, isout);
463:     return(0);
464:   }
465:   if (len < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Negative array length: %D", len);
466:   N = 0;
467:   for (i = 0; i < len; ++i) {
468:     if (islist[i]) {
469:       ISGetLocalSize(islist[i], &n);
470:       N   += n;
471:     }
472:   }
473:   PetscMalloc1(N, &idx);
474:   N = 0;
475:   for (i = 0; i < len; ++i) {
476:     if (islist[i]) {
477:       ISGetLocalSize(islist[i], &n);
478:       ISGetIndices(islist[i], &iidx);
479:       PetscArraycpy(idx+N,iidx, n);
480:       ISRestoreIndices(islist[i], &iidx);
481:       N   += n;
482:     }
483:   }
484:   ISCreateGeneral(comm, N, idx, PETSC_OWN_POINTER, isout);
485:   return(0);
486: }

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


494:   Collective.

496:   Input arguments:
497: + comm    -  MPI_Comm
498: . listlen -  IS list length
499: - islist  -  IS list

501:   Output arguments:
502: + xis -  domain IS
503: - yis -  range  IS

505:   Level: advanced

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

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

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

520: .seealso ISPairToList()
521: @*/
522: PetscErrorCode ISListToPair(MPI_Comm comm, PetscInt listlen, IS islist[], IS *xis, IS *yis)
523: {
525:   PetscInt       ncolors, *colors,i, leni,len,*xinds, *yinds,k,j;
526:   const PetscInt *indsi;

529:   PetscMalloc1(listlen, &colors);
530:   PetscObjectsListGetGlobalNumbering(comm, listlen, (PetscObject*)islist,&ncolors, colors);
531:   len  = 0;
532:   for (i = 0; i < listlen; ++i) {
533:     ISGetLocalSize(islist[i], &leni);
534:     len += leni;
535:   }
536:   PetscMalloc1(len, &xinds);
537:   PetscMalloc1(len, &yinds);
538:   k    = 0;
539:   for (i = 0; i < listlen; ++i) {
540:     ISGetLocalSize(islist[i], &leni);
541:     ISGetIndices(islist[i],&indsi);
542:     for (j = 0; j < leni; ++j) {
543:       xinds[k] = indsi[j];
544:       yinds[k] = colors[i];
545:       ++k;
546:     }
547:   }
548:   PetscFree(colors);
549:   ISCreateGeneral(comm,len,xinds,PETSC_OWN_POINTER,xis);
550:   ISCreateGeneral(comm,len,yinds,PETSC_OWN_POINTER,yis);
551:   return(0);
552: }


555: /*@
556:    ISPairToList   -   convert an IS pair encoding an integer map to a list of ISs.
557:                      Each IS on the output list contains the preimage for each index on the second input IS.
558:                      The ISs on the output list are constructed on the subcommunicators of the input IS pair.
559:                      Each subcommunicator corresponds to the preimage of some index j -- this subcomm contains
560:                      exactly the ranks that assign some indices i to j.  This is essentially the inverse of
561:                      ISListToPair().

563:   Collective on indis.

565:   Input arguments:
566: + xis -  domain IS
567: - yis -  range IS

569:   Output arguments:
570: + listlen -  length of islist
571: - islist  -  list of ISs breaking up indis by color

573:   Note:
574: + xis and yis must be of the same length and have congruent communicators.
575: - The resulting ISs have subcommunicators in a "deadlock-free" order (see ISListToPair()).

577:   Level: advanced

579: .seealso ISListToPair()
580:  @*/
581: PetscErrorCode ISPairToList(IS xis, IS yis, PetscInt *listlen, IS **islist)
582: {
584:   IS             indis = xis, coloris = yis;
585:   PetscInt       *inds, *colors, llen, ilen, lstart, lend, lcount,l;
586:   PetscMPIInt    rank, size, llow, lhigh, low, high,color,subsize;
587:   const PetscInt *ccolors, *cinds;
588:   MPI_Comm       comm, subcomm;

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

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


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

688:   Not collective.

690:   Input arguments:
691: + a    -  IS to embed
692: . b    -  IS to embed into
693: - drop -  flag indicating whether to drop a's indices that are not in b.

695:   Output arguments:
696: . c    -  local embedding indices

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

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

705:   Level: advanced

707: .seealso ISLocalToGlobalMapping
708:  @*/
709: PetscErrorCode ISEmbed(IS a, IS b, PetscBool drop, IS *c)
710: {
711:   PetscErrorCode             ierr;
712:   ISLocalToGlobalMapping     ltog;
713:   ISGlobalToLocalMappingMode gtoltype = IS_GTOLM_DROP;
714:   PetscInt                   alen, clen, *cindices, *cindices2;
715:   const PetscInt             *aindices;

721:   ISLocalToGlobalMappingCreateIS(b, &ltog);
722:   ISGetLocalSize(a, &alen);
723:   ISGetIndices(a, &aindices);
724:   PetscMalloc1(alen, &cindices);
725:   if (!drop) gtoltype = IS_GTOLM_MASK;
726:   ISGlobalToLocalMappingApply(ltog,gtoltype,alen,aindices,&clen,cindices);
727:   ISLocalToGlobalMappingDestroy(&ltog);
728:   if (clen != alen) {
729:     cindices2 = cindices;
730:     PetscMalloc1(clen, &cindices);
731:     PetscArraycpy(cindices,cindices2,clen);
732:     PetscFree(cindices2);
733:   }
734:   ISCreateGeneral(PETSC_COMM_SELF,clen,cindices,PETSC_OWN_POINTER,c);
735:   return(0);
736: }


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

742:   Not collective.

744:   Input arguments:
745: + f      -  IS to sort
746: - always -  build the permutation even when f's indices are nondecreasing.

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


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

757:   Level: advanced

759: .seealso ISLocalToGlobalMapping, ISSort()
760:  @*/
761: PetscErrorCode ISSortPermutation(IS f,PetscBool always,IS *h)
762: {
763:   PetscErrorCode  ierr;
764:   const PetscInt  *findices;
765:   PetscInt        fsize,*hindices,i;
766:   PetscBool       isincreasing;

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