Actual source code: index.c

petsc-3.10.5 2019-03-28
Report Typos and Errors

  2: /*
  3:    Defines the abstract operations on index sets, i.e. the public interface.
  4: */
  5:  #include <petsc/private/isimpl.h>
  6:  #include <petscviewer.h>
  7:  #include <petscsf.h>

  9: /* Logging support */
 10: PetscClassId IS_CLASSID;

 12: /*@
 13:    ISRenumber - Renumbers an index set (with multiplicities) in a contiguous way.

 15:    Collective on IS

 17:    Input Parmeters:
 18: .  subset - the index set
 19: .  subset_mult - the multiplcity of each entry in subset (optional, can be NULL)

 21:    Output Parameters:
 22: .  N - the maximum entry of the new IS
 23: .  subset_n - the new IS

 25:    Level: intermediate

 27: .seealso:
 28: @*/
 29: PetscErrorCode ISRenumber(IS subset, IS subset_mult, PetscInt *N, IS *subset_n)
 30: {
 31:   PetscSF        sf;
 32:   PetscLayout    map;
 33:   const PetscInt *idxs;
 34:   PetscInt       *leaf_data,*root_data,*gidxs;
 35:   PetscInt       N_n,n,i,lbounds[2],gbounds[2],Nl;
 36:   PetscInt       n_n,nlocals,start,first_index;
 37:   PetscMPIInt    commsize;
 38:   PetscBool      first_found;

 42:   ISGetLocalSize(subset,&n);
 43:   if (subset_mult) {
 45:     ISGetLocalSize(subset,&i);
 46:     if (i != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Local subset and multiplicity sizes don't match! %d != %d",n,i);
 47:   }
 48:   /* create workspace layout for computing global indices of subset */
 49:   ISGetIndices(subset,&idxs);
 50:   lbounds[0] = lbounds[1] = 0;
 51:   for (i=0;i<n;i++) {
 52:     if (idxs[i] < lbounds[0]) lbounds[0] = idxs[i];
 53:     else if (idxs[i] > lbounds[1]) lbounds[1] = idxs[i];
 54:   }
 55:   lbounds[0] = -lbounds[0];
 56:   MPIU_Allreduce(lbounds,gbounds,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)subset));
 57:   gbounds[0] = -gbounds[0];
 58:   N_n= gbounds[1] - gbounds[0] + 1;
 59:   PetscLayoutCreate(PetscObjectComm((PetscObject)subset),&map);
 60:   PetscLayoutSetBlockSize(map,1);
 61:   PetscLayoutSetSize(map,N_n);
 62:   PetscLayoutSetUp(map);
 63:   PetscLayoutGetLocalSize(map,&Nl);

 65:   /* create sf : leaf_data == multiplicity of indexes, root data == global index in layout */
 66:   PetscMalloc2(n,&leaf_data,Nl,&root_data);
 67:   if (subset_mult) {
 68:     const PetscInt* idxs_mult;

 70:     ISGetIndices(subset_mult,&idxs_mult);
 71:     PetscMemcpy(leaf_data,idxs_mult,n*sizeof(PetscInt));
 72:     ISRestoreIndices(subset_mult,&idxs_mult);
 73:   } else {
 74:     for (i=0;i<n;i++) leaf_data[i] = 1;
 75:   }
 76:   /* local size of new subset */
 77:   n_n = 0;
 78:   for (i=0;i<n;i++) n_n += leaf_data[i];

 80:   /* global indexes in layout */
 81:   PetscMalloc1(n_n,&gidxs); /* allocating possibly extra space in gidxs which will be used later */
 82:   for (i=0;i<n;i++) gidxs[i] = idxs[i] - gbounds[0];
 83:   ISRestoreIndices(subset,&idxs);
 84:   PetscSFCreate(PetscObjectComm((PetscObject)subset),&sf);
 85:   PetscSFSetGraphLayout(sf,map,n,NULL,PETSC_COPY_VALUES,gidxs);
 86:   PetscLayoutDestroy(&map);

 88:   /* reduce from leaves to roots */
 89:   PetscMemzero(root_data,Nl*sizeof(PetscInt));
 90:   PetscSFReduceBegin(sf,MPIU_INT,leaf_data,root_data,MPI_MAX);
 91:   PetscSFReduceEnd(sf,MPIU_INT,leaf_data,root_data,MPI_MAX);

 93:   /* count indexes in local part of layout */
 94:   nlocals = 0;
 95:   first_index = -1;
 96:   first_found = PETSC_FALSE;
 97:   for (i=0;i<Nl;i++) {
 98:     if (!first_found && root_data[i]) {
 99:       first_found = PETSC_TRUE;
100:       first_index = i;
101:     }
102:     nlocals += root_data[i];
103:   }

105:   /* cumulative of number of indexes and size of subset without holes */
106: #if defined(PETSC_HAVE_MPI_EXSCAN)
107:   start = 0;
108:   MPI_Exscan(&nlocals,&start,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)subset));
109: #else
110:   MPI_Scan(&nlocals,&start,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)subset));
111:   start = start-nlocals;
112: #endif

114:   if (N) { /* compute total size of new subset if requested */
115:     *N = start + nlocals;
116:     MPI_Comm_size(PetscObjectComm((PetscObject)subset),&commsize);
117:     MPI_Bcast(N,1,MPIU_INT,commsize-1,PetscObjectComm((PetscObject)subset));
118:   }

120:   /* adapt root data with cumulative */
121:   if (first_found) {
122:     PetscInt old_index;

124:     root_data[first_index] += start;
125:     old_index = first_index;
126:     for (i=first_index+1;i<Nl;i++) {
127:       if (root_data[i]) {
128:         root_data[i] += root_data[old_index];
129:         old_index = i;
130:       }
131:     }
132:   }

134:   /* from roots to leaves */
135:   PetscSFBcastBegin(sf,MPIU_INT,root_data,leaf_data);
136:   PetscSFBcastEnd(sf,MPIU_INT,root_data,leaf_data);
137:   PetscSFDestroy(&sf);

139:   /* create new IS with global indexes without holes */
140:   if (subset_mult) {
141:     const PetscInt* idxs_mult;
142:     PetscInt        cum;

144:     cum = 0;
145:     ISGetIndices(subset_mult,&idxs_mult);
146:     for (i=0;i<n;i++) {
147:       PetscInt j;
148:       for (j=0;j<idxs_mult[i];j++) gidxs[cum++] = leaf_data[i] - idxs_mult[i] + j;
149:     }
150:     ISRestoreIndices(subset_mult,&idxs_mult);
151:   } else {
152:     for (i=0;i<n;i++) {
153:       gidxs[i] = leaf_data[i]-1;
154:     }
155:   }
156:   ISCreateGeneral(PetscObjectComm((PetscObject)subset),n_n,gidxs,PETSC_OWN_POINTER,subset_n);
157:   PetscFree2(leaf_data,root_data);
158:   return(0);
159: }

161: /*@
162:    ISIdentity - Determines whether index set is the identity mapping.

164:    Collective on IS

166:    Input Parmeters:
167: .  is - the index set

169:    Output Parameters:
170: .  ident - PETSC_TRUE if an identity, else PETSC_FALSE

172:    Level: intermediate

174:    Concepts: identity mapping
175:    Concepts: index sets^is identity

177: .seealso: ISSetIdentity()
178: @*/
179: PetscErrorCode  ISIdentity(IS is,PetscBool  *ident)
180: {

186:   *ident = is->isidentity;
187:   if (*ident) return(0);
188:   if (is->ops->identity) {
189:     (*is->ops->identity)(is,ident);
190:   }
191:   return(0);
192: }

194: /*@
195:    ISSetIdentity - Informs the index set that it is an identity.

197:    Logically Collective on IS

199:    Input Parmeters:
200: .  is - the index set

202:    Level: intermediate

204:    Concepts: identity mapping
205:    Concepts: index sets^is identity

207: .seealso: ISIdentity()
208: @*/
209: PetscErrorCode  ISSetIdentity(IS is)
210: {

215:   is->isidentity = PETSC_TRUE;
216:   ISSetPermutation(is);
217:   return(0);
218: }

220: /*@
221:    ISContiguousLocal - Locates an index set with contiguous range within a global range, if possible

223:    Not Collective

225:    Input Parmeters:
226: +  is - the index set
227: .  gstart - global start
228: .  gend - global end

230:    Output Parameters:
231: +  start - start of contiguous block, as an offset from gstart
232: -  contig - PETSC_TRUE if the index set refers to contiguous entries on this process, else PETSC_FALSE

234:    Level: developer

236:    Concepts: index sets^is contiguous

238: .seealso: ISGetLocalSize(), VecGetOwnershipRange()
239: @*/
240: PetscErrorCode  ISContiguousLocal(IS is,PetscInt gstart,PetscInt gend,PetscInt *start,PetscBool *contig)
241: {

248:   if (is->ops->contiguous) {
249:     (*is->ops->contiguous)(is,gstart,gend,start,contig);
250:   } else {
251:     *start  = -1;
252:     *contig = PETSC_FALSE;
253:   }
254:   return(0);
255: }

257: /*@
258:    ISPermutation - PETSC_TRUE or PETSC_FALSE depending on whether the
259:    index set has been declared to be a permutation.

261:    Logically Collective on IS

263:    Input Parmeters:
264: .  is - the index set

266:    Output Parameters:
267: .  perm - PETSC_TRUE if a permutation, else PETSC_FALSE

269:    Level: intermediate

271:   Concepts: permutation
272:   Concepts: index sets^is permutation

274: .seealso: ISSetPermutation()
275: @*/
276: PetscErrorCode  ISPermutation(IS is,PetscBool  *perm)
277: {
281:   *perm = (PetscBool) is->isperm;
282:   return(0);
283: }

285: /*@
286:    ISSetPermutation - Informs the index set that it is a permutation.

288:    Logically Collective on IS

290:    Input Parmeters:
291: .  is - the index set

293:    Level: intermediate

295:   Concepts: permutation
296:   Concepts: index sets^permutation

298:    The debug version of the libraries (./configure --with-debugging=1) checks if the
299:   index set is actually a permutation. The optimized version just believes you.

301: .seealso: ISPermutation()
302: @*/
303: PetscErrorCode  ISSetPermutation(IS is)
304: {
307: #if defined(PETSC_USE_DEBUG)
308:   {
309:     PetscMPIInt    size;

312:     MPI_Comm_size(PetscObjectComm((PetscObject)is),&size);
313:     if (size == 1) {
314:       PetscInt       i,n,*idx;
315:       const PetscInt *iidx;

317:       ISGetSize(is,&n);
318:       PetscMalloc1(n,&idx);
319:       ISGetIndices(is,&iidx);
320:       PetscMemcpy(idx,iidx,n*sizeof(PetscInt));
321:       PetscSortInt(n,idx);
322:       for (i=0; i<n; i++) {
323:         if (idx[i] != i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Index set is not a permutation");
324:       }
325:       PetscFree(idx);
326:       ISRestoreIndices(is,&iidx);
327:     }
328:   }
329: #endif
330:   is->isperm = PETSC_TRUE;
331:   return(0);
332: }

334: /*@
335:    ISDestroy - Destroys an index set.

337:    Collective on IS

339:    Input Parameters:
340: .  is - the index set

342:    Level: beginner

344: .seealso: ISCreateGeneral(), ISCreateStride(), ISCreateBlocked()
345: @*/
346: PetscErrorCode  ISDestroy(IS *is)
347: {

351:   if (!*is) return(0);
353:   if (--((PetscObject)(*is))->refct > 0) {*is = 0; return(0);}
354:   if ((*is)->complement) {
355:     PetscInt refcnt;
356:     PetscObjectGetReference((PetscObject)((*is)->complement), &refcnt);
357:     if (refcnt > 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Nonlocal IS has not been restored");
358:     ISDestroy(&(*is)->complement);
359:   }
360:   if ((*is)->ops->destroy) {
361:     (*(*is)->ops->destroy)(*is);
362:   }
363:   PetscLayoutDestroy(&(*is)->map);
364:   /* Destroy local representations of offproc data. */
365:   PetscFree((*is)->total);
366:   PetscFree((*is)->nonlocal);
367:   PetscHeaderDestroy(is);
368:   return(0);
369: }

371: /*@
372:    ISInvertPermutation - Creates a new permutation that is the inverse of
373:                          a given permutation.

375:    Collective on IS

377:    Input Parameter:
378: +  is - the index set
379: -  nlocal - number of indices on this processor in result (ignored for 1 proccessor) or
380:             use PETSC_DECIDE

382:    Output Parameter:
383: .  isout - the inverse permutation

385:    Level: intermediate

387:    Notes:
388:     For parallel index sets this does the complete parallel permutation, but the
389:     code is not efficient for huge index sets (10,000,000 indices).

391:    Concepts: inverse permutation
392:    Concepts: permutation^inverse
393:    Concepts: index sets^inverting
394: @*/
395: PetscErrorCode  ISInvertPermutation(IS is,PetscInt nlocal,IS *isout)
396: {

402:   if (!is->isperm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Not a permutation, must call ISSetPermutation() on the IS first");
403:   if (is->isidentity) {
404:     ISDuplicate(is,isout);
405:   } else {
406:     (*is->ops->invertpermutation)(is,nlocal,isout);
407:     ISSetPermutation(*isout);
408:   }
409:   return(0);
410: }

412: /*@
413:    ISGetSize - Returns the global length of an index set.

415:    Not Collective

417:    Input Parameter:
418: .  is - the index set

420:    Output Parameter:
421: .  size - the global size

423:    Level: beginner

425:    Concepts: size^of index set
426:    Concepts: index sets^size

428: @*/
429: PetscErrorCode  ISGetSize(IS is,PetscInt *size)
430: {

436:   (*is->ops->getsize)(is,size);
437:   return(0);
438: }

440: /*@
441:    ISGetLocalSize - Returns the local (processor) length of an index set.

443:    Not Collective

445:    Input Parameter:
446: .  is - the index set

448:    Output Parameter:
449: .  size - the local size

451:    Level: beginner

453:    Concepts: size^of index set
454:    Concepts: local size^of index set
455:    Concepts: index sets^local size

457: @*/
458: PetscErrorCode  ISGetLocalSize(IS is,PetscInt *size)
459: {

465:   (*is->ops->getlocalsize)(is,size);
466:   return(0);
467: }

469: /*@C
470:    ISGetIndices - Returns a pointer to the indices.  The user should call
471:    ISRestoreIndices() after having looked at the indices.  The user should
472:    NOT change the indices.

474:    Not Collective

476:    Input Parameter:
477: .  is - the index set

479:    Output Parameter:
480: .  ptr - the location to put the pointer to the indices

482:    Fortran Note:
483:    This routine has two different interfaces from Fortran; the first is not recommend, it does not require Fortran 90
484: $    IS          is
485: $    integer     is_array(1)
486: $    PetscOffset i_is
487: $    int         ierr
488: $       call ISGetIndices(is,is_array,i_is,ierr)
489: $
490: $   Access first local entry in list
491: $      value = is_array(i_is + 1)
492: $
493: $      ...... other code
494: $       call ISRestoreIndices(is,is_array,i_is,ierr)
495:    The second Fortran interface is recommended.
496: $          use petscisdef
497: $          PetscInt, pointer :: array(:)
498: $          PetscErrorCode  ierr
499: $          IS       i
500: $          call ISGetIndicesF90(i,array,ierr)



504:    See the Fortran chapter of the users manual and
505:    petsc/src/is/examples/[tutorials,tests] for details.

507:    Level: intermediate

509:    Concepts: index sets^getting indices
510:    Concepts: indices of index set

512: .seealso: ISRestoreIndices(), ISGetIndicesF90()
513: @*/
514: PetscErrorCode  ISGetIndices(IS is,const PetscInt *ptr[])
515: {

521:   (*is->ops->getindices)(is,ptr);
522:   return(0);
523: }

525: /*@C
526:    ISGetMinMax - Gets the minimum and maximum values in an IS

528:    Not Collective

530:    Input Parameter:
531: .  is - the index set

533:    Output Parameter:
534: +   min - the minimum value
535: -   max - the maximum value

537:    Level: intermediate

539:    Notes:
540:     Empty index sets return min=PETSC_MAX_INT and max=PETSC_MIN_INT.
541:     In parallel, it returns the min and max of the local portion of the IS

543:    Concepts: index sets^getting indices
544:    Concepts: indices of index set

546: .seealso: ISGetIndices(), ISRestoreIndices(), ISGetIndicesF90()
547: @*/
548: PetscErrorCode  ISGetMinMax(IS is,PetscInt *min,PetscInt *max)
549: {
552:   if (min) *min = is->min;
553:   if (max) *max = is->max;
554:   return(0);
555: }

557: /*@
558:   ISLocate - determine the location of an index within the local component of an index set

560:   Not Collective

562:   Input Parameter:
563: + is - the index set
564: - key - the search key

566:   Output Parameter:
567: . location - if >= 0, a location within the index set that is equal to the key, otherwise the key is not in the index set

569:   Level: intermediate
570: @*/
571: PetscErrorCode ISLocate(IS is, PetscInt key, PetscInt *location)
572: {

576:   if (is->ops->locate) {
577:     (*is->ops->locate)(is,key,location);
578:   } else {
579:     PetscInt       numIdx;
580:     PetscBool      sorted;
581:     const PetscInt *idx;

583:     ISGetLocalSize(is,&numIdx);
584:     ISGetIndices(is,&idx);
585:     ISSorted(is,&sorted);
586:     if (sorted) {
587:       PetscFindInt(key,numIdx,idx,location);
588:     } else {
589:       PetscInt i;

591:       *location = -1;
592:       for (i = 0; i < numIdx; i++) {
593:         if (idx[i] == key) {
594:           *location = i;
595:           break;
596:         }
597:       }
598:     }
599:     ISRestoreIndices(is,&idx);
600:   }
601:   return(0);
602: }

604: /*@C
605:    ISRestoreIndices - Restores an index set to a usable state after a call
606:                       to ISGetIndices().

608:    Not Collective

610:    Input Parameters:
611: +  is - the index set
612: -  ptr - the pointer obtained by ISGetIndices()

614:    Fortran Note:
615:    This routine is used differently from Fortran
616: $    IS          is
617: $    integer     is_array(1)
618: $    PetscOffset i_is
619: $    int         ierr
620: $       call ISGetIndices(is,is_array,i_is,ierr)
621: $
622: $   Access first local entry in list
623: $      value = is_array(i_is + 1)
624: $
625: $      ...... other code
626: $       call ISRestoreIndices(is,is_array,i_is,ierr)

628:    See the Fortran chapter of the users manual and
629:    petsc/src/is/examples/[tutorials,tests] for details.

631:    Level: intermediate

633:    Note:
634:    This routine zeros out ptr. This is to prevent accidental us of the array after it has been restored.

636: .seealso: ISGetIndices(), ISRestoreIndicesF90()
637: @*/
638: PetscErrorCode  ISRestoreIndices(IS is,const PetscInt *ptr[])
639: {

645:   if (is->ops->restoreindices) {
646:     (*is->ops->restoreindices)(is,ptr);
647:   }
648:   return(0);
649: }

651: static PetscErrorCode ISGatherTotal_Private(IS is)
652: {
654:   PetscInt       i,n,N;
655:   const PetscInt *lindices;
656:   MPI_Comm       comm;
657:   PetscMPIInt    rank,size,*sizes = NULL,*offsets = NULL,nn;


662:   PetscObjectGetComm((PetscObject)is,&comm);
663:   MPI_Comm_size(comm,&size);
664:   MPI_Comm_rank(comm,&rank);
665:   ISGetLocalSize(is,&n);
666:   PetscMalloc2(size,&sizes,size,&offsets);

668:   PetscMPIIntCast(n,&nn);
669:   MPI_Allgather(&nn,1,MPI_INT,sizes,1,MPI_INT,comm);
670:   offsets[0] = 0;
671:   for (i=1; i<size; ++i) offsets[i] = offsets[i-1] + sizes[i-1];
672:   N = offsets[size-1] + sizes[size-1];

674:   PetscMalloc1(N,&(is->total));
675:   ISGetIndices(is,&lindices);
676:   MPI_Allgatherv((void*)lindices,nn,MPIU_INT,is->total,sizes,offsets,MPIU_INT,comm);
677:   ISRestoreIndices(is,&lindices);
678:   is->local_offset = offsets[rank];
679:   PetscFree2(sizes,offsets);
680:   return(0);
681: }

683: /*@C
684:    ISGetTotalIndices - Retrieve an array containing all indices across the communicator.

686:    Collective on IS

688:    Input Parameter:
689: .  is - the index set

691:    Output Parameter:
692: .  indices - total indices with rank 0 indices first, and so on; total array size is
693:              the same as returned with ISGetSize().

695:    Level: intermediate

697:    Notes:
698:     this is potentially nonscalable, but depends on the size of the total index set
699:      and the size of the communicator. This may be feasible for index sets defined on
700:      subcommunicators, such that the set size does not grow with PETSC_WORLD_COMM.
701:      Note also that there is no way to tell where the local part of the indices starts
702:      (use ISGetIndices() and ISGetNonlocalIndices() to retrieve just the local and just
703:       the nonlocal part (complement), respectively).

705:    Concepts: index sets^getting nonlocal indices
706: .seealso: ISRestoreTotalIndices(), ISGetNonlocalIndices(), ISGetSize()
707: @*/
708: PetscErrorCode ISGetTotalIndices(IS is, const PetscInt *indices[])
709: {
711:   PetscMPIInt    size;

716:   MPI_Comm_size(PetscObjectComm((PetscObject)is), &size);
717:   if (size == 1) {
718:     (*is->ops->getindices)(is,indices);
719:   } else {
720:     if (!is->total) {
721:       ISGatherTotal_Private(is);
722:     }
723:     *indices = is->total;
724:   }
725:   return(0);
726: }

728: /*@C
729:    ISRestoreTotalIndices - Restore the index array obtained with ISGetTotalIndices().

731:    Not Collective.

733:    Input Parameter:
734: +  is - the index set
735: -  indices - index array; must be the array obtained with ISGetTotalIndices()

737:    Level: intermediate

739:    Concepts: index sets^getting nonlocal indices
740:    Concepts: index sets^restoring nonlocal indices
741: .seealso: ISRestoreTotalIndices(), ISGetNonlocalIndices()
742: @*/
743: PetscErrorCode  ISRestoreTotalIndices(IS is, const PetscInt *indices[])
744: {
746:   PetscMPIInt    size;

751:   MPI_Comm_size(PetscObjectComm((PetscObject)is), &size);
752:   if (size == 1) {
753:     (*is->ops->restoreindices)(is,indices);
754:   } else {
755:     if (is->total != *indices) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Index array pointer being restored does not point to the array obtained from the IS.");
756:   }
757:   return(0);
758: }
759: /*@C
760:    ISGetNonlocalIndices - Retrieve an array of indices from remote processors
761:                        in this communicator.

763:    Collective on IS

765:    Input Parameter:
766: .  is - the index set

768:    Output Parameter:
769: .  indices - indices with rank 0 indices first, and so on,  omitting
770:              the current rank.  Total number of indices is the difference
771:              total and local, obtained with ISGetSize() and ISGetLocalSize(),
772:              respectively.

774:    Level: intermediate

776:    Notes:
777:     restore the indices using ISRestoreNonlocalIndices().
778:           The same scalability considerations as those for ISGetTotalIndices
779:           apply here.

781:    Concepts: index sets^getting nonlocal indices
782: .seealso: ISGetTotalIndices(), ISRestoreNonlocalIndices(), ISGetSize(), ISGetLocalSize().
783: @*/
784: PetscErrorCode  ISGetNonlocalIndices(IS is, const PetscInt *indices[])
785: {
787:   PetscMPIInt    size;
788:   PetscInt       n, N;

793:   MPI_Comm_size(PetscObjectComm((PetscObject)is), &size);
794:   if (size == 1) *indices = NULL;
795:   else {
796:     if (!is->total) {
797:       ISGatherTotal_Private(is);
798:     }
799:     ISGetLocalSize(is,&n);
800:     ISGetSize(is,&N);
801:     PetscMalloc1(N-n, &(is->nonlocal));
802:     PetscMemcpy(is->nonlocal, is->total, sizeof(PetscInt)*is->local_offset);
803:     PetscMemcpy(is->nonlocal+is->local_offset, is->total+is->local_offset+n, sizeof(PetscInt)*(N - is->local_offset - n));
804:     *indices = is->nonlocal;
805:   }
806:   return(0);
807: }

809: /*@C
810:    ISRestoreTotalIndices - Restore the index array obtained with ISGetNonlocalIndices().

812:    Not Collective.

814:    Input Parameter:
815: +  is - the index set
816: -  indices - index array; must be the array obtained with ISGetNonlocalIndices()

818:    Level: intermediate

820:    Concepts: index sets^getting nonlocal indices
821:    Concepts: index sets^restoring nonlocal indices
822: .seealso: ISGetTotalIndices(), ISGetNonlocalIndices(), ISRestoreTotalIndices()
823: @*/
824: PetscErrorCode  ISRestoreNonlocalIndices(IS is, const PetscInt *indices[])
825: {
829:   if (is->nonlocal != *indices) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Index array pointer being restored does not point to the array obtained from the IS.");
830:   return(0);
831: }

833: /*@
834:    ISGetNonlocalIS - Gather all nonlocal indices for this IS and present
835:                      them as another sequential index set.


838:    Collective on IS

840:    Input Parameter:
841: .  is - the index set

843:    Output Parameter:
844: .  complement - sequential IS with indices identical to the result of
845:                 ISGetNonlocalIndices()

847:    Level: intermediate

849:    Notes:
850:     complement represents the result of ISGetNonlocalIndices as an IS.
851:           Therefore scalability issues similar to ISGetNonlocalIndices apply.
852:           The resulting IS must be restored using ISRestoreNonlocalIS().

854:    Concepts: index sets^getting nonlocal indices
855: .seealso: ISGetNonlocalIndices(), ISRestoreNonlocalIndices(),  ISAllGather(), ISGetSize()
856: @*/
857: PetscErrorCode  ISGetNonlocalIS(IS is, IS *complement)
858: {

864:   /* Check if the complement exists already. */
865:   if (is->complement) {
866:     *complement = is->complement;
867:     PetscObjectReference((PetscObject)(is->complement));
868:   } else {
869:     PetscInt       N, n;
870:     const PetscInt *idx;
871:     ISGetSize(is, &N);
872:     ISGetLocalSize(is,&n);
873:     ISGetNonlocalIndices(is, &idx);
874:     ISCreateGeneral(PETSC_COMM_SELF, N-n,idx, PETSC_USE_POINTER, &(is->complement));
875:     PetscObjectReference((PetscObject)is->complement);
876:     *complement = is->complement;
877:   }
878:   return(0);
879: }


882: /*@
883:    ISRestoreNonlocalIS - Restore the IS obtained with ISGetNonlocalIS().

885:    Not collective.

887:    Input Parameter:
888: +  is         - the index set
889: -  complement - index set of is's nonlocal indices

891:    Level: intermediate


894:    Concepts: index sets^getting nonlocal indices
895:    Concepts: index sets^restoring nonlocal indices
896: .seealso: ISGetNonlocalIS(), ISGetNonlocalIndices(), ISRestoreNonlocalIndices()
897: @*/
898: PetscErrorCode  ISRestoreNonlocalIS(IS is, IS *complement)
899: {
901:   PetscInt       refcnt;

906:   if (*complement != is->complement) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Complement IS being restored was not obtained with ISGetNonlocalIS()");
907:   PetscObjectGetReference((PetscObject)(is->complement), &refcnt);
908:   if (refcnt <= 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Duplicate call to ISRestoreNonlocalIS() detected");
909:   PetscObjectDereference((PetscObject)(is->complement));
910:   return(0);
911: }

913: /*@C
914:    ISView - Displays an index set.

916:    Collective on IS

918:    Input Parameters:
919: +  is - the index set
920: -  viewer - viewer used to display the set, for example PETSC_VIEWER_STDOUT_SELF.

922:    Level: intermediate

924: .seealso: PetscViewerASCIIOpen()
925: @*/
926: PetscErrorCode  ISView(IS is,PetscViewer viewer)
927: {

932:   if (!viewer) {
933:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)is),&viewer);
934:   }

938:   PetscObjectPrintClassNamePrefixType((PetscObject)is,viewer);
939:   (*is->ops->view)(is,viewer);
940:   return(0);
941: }

943: /*@
944:   ISLoad - Loads a vector that has been stored in binary or HDF5 format with ISView().

946:   Collective on PetscViewer

948:   Input Parameters:
949: + is - the newly loaded vector, this needs to have been created with ISCreate() or some related function before a call to ISLoad().
950: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or HDF5 file viewer, obtained from PetscViewerHDF5Open()

952:   Level: intermediate

954:   Notes:
955:   IF using HDF5, you must assign the IS the same name as was used in the IS
956:   that was stored in the file using PetscObjectSetName(). Otherwise you will
957:   get the error message: "Cannot H5DOpen2() with Vec name NAMEOFOBJECT"

959:   Concepts: index set^loading from file

961: .seealso: PetscViewerBinaryOpen(), ISView(), MatLoad(), VecLoad()
962: @*/
963: PetscErrorCode ISLoad(IS is, PetscViewer viewer)
964: {
965:   PetscBool      isbinary, ishdf5;

971:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERBINARY, &isbinary);
972:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
973:   if (!isbinary && !ishdf5) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
974:   if (!((PetscObject)is)->type_name) {ISSetType(is, ISGENERAL);}
975:   (*is->ops->load)(is, viewer);
976:   return(0);
977: }

979: /*@
980:    ISSort - Sorts the indices of an index set.

982:    Collective on IS

984:    Input Parameters:
985: .  is - the index set

987:    Level: intermediate

989:    Concepts: index sets^sorting
990:    Concepts: sorting^index set

992: .seealso: ISSortRemoveDups(), ISSorted()
993: @*/
994: PetscErrorCode  ISSort(IS is)
995: {

1000:   (*is->ops->sort)(is);
1001:   return(0);
1002: }

1004: /*@
1005:   ISSortRemoveDups - Sorts the indices of an index set, removing duplicates.

1007:   Collective on IS

1009:   Input Parameters:
1010: . is - the index set

1012:   Level: intermediate

1014:   Concepts: index sets^sorting
1015:   Concepts: sorting^index set

1017: .seealso: ISSort(), ISSorted()
1018: @*/
1019: PetscErrorCode ISSortRemoveDups(IS is)
1020: {

1025:   (*is->ops->sortremovedups)(is);
1026:   return(0);
1027: }

1029: /*@
1030:    ISToGeneral - Converts an IS object of any type to ISGENERAL type

1032:    Collective on IS

1034:    Input Parameters:
1035: .  is - the index set

1037:    Level: intermediate

1039:    Concepts: index sets^sorting
1040:    Concepts: sorting^index set

1042: .seealso: ISSorted()
1043: @*/
1044: PetscErrorCode  ISToGeneral(IS is)
1045: {

1050:   if (is->ops->togeneral) {
1051:     (*is->ops->togeneral)(is);
1052:   } else SETERRQ1(PetscObjectComm((PetscObject)is),PETSC_ERR_SUP,"Not written for this type %s",((PetscObject)is)->type_name);
1053:   return(0);
1054: }

1056: /*@
1057:    ISSorted - Checks the indices to determine whether they have been sorted.

1059:    Collective on IS

1061:    Input Parameter:
1062: .  is - the index set

1064:    Output Parameter:
1065: .  flg - output flag, either PETSC_TRUE if the index set is sorted,
1066:          or PETSC_FALSE otherwise.

1068:    Notes:
1069:     For parallel IS objects this only indicates if the local part of the IS
1070:           is sorted. So some processors may return PETSC_TRUE while others may
1071:           return PETSC_FALSE.

1073:    Level: intermediate

1075: .seealso: ISSort(), ISSortRemoveDups()
1076: @*/
1077: PetscErrorCode  ISSorted(IS is,PetscBool  *flg)
1078: {

1084:   (*is->ops->sorted)(is,flg);
1085:   return(0);
1086: }

1088: /*@
1089:    ISDuplicate - Creates a duplicate copy of an index set.

1091:    Collective on IS

1093:    Input Parmeters:
1094: .  is - the index set

1096:    Output Parameters:
1097: .  isnew - the copy of the index set

1099:    Level: beginner

1101:    Concepts: index sets^duplicating

1103: .seealso: ISCreateGeneral(), ISCopy()
1104: @*/
1105: PetscErrorCode  ISDuplicate(IS is,IS *newIS)
1106: {

1112:   (*is->ops->duplicate)(is,newIS);
1113:   (*newIS)->isidentity = is->isidentity;
1114:   (*newIS)->isperm     = is->isperm;
1115:   return(0);
1116: }

1118: /*@
1119:    ISCopy - Copies an index set.

1121:    Collective on IS

1123:    Input Parmeters:
1124: .  is - the index set

1126:    Output Parameters:
1127: .  isy - the copy of the index set

1129:    Level: beginner

1131:    Concepts: index sets^copying

1133: .seealso: ISDuplicate()
1134: @*/
1135: PetscErrorCode  ISCopy(IS is,IS isy)
1136: {

1143:   if (is == isy) return(0);
1144:   (*is->ops->copy)(is,isy);
1145:   isy->isperm     = is->isperm;
1146:   isy->max        = is->max;
1147:   isy->min        = is->min;
1148:   isy->isidentity = is->isidentity;
1149:   return(0);
1150: }

1152: /*@
1153:    ISOnComm - Split a parallel IS on subcomms (usually self) or concatenate index sets on subcomms into a parallel index set

1155:    Collective on IS and comm

1157:    Input Arguments:
1158: + is - index set
1159: . comm - communicator for new index set
1160: - mode - copy semantics, PETSC_USE_POINTER for no-copy if possible, otherwise PETSC_COPY_VALUES

1162:    Output Arguments:
1163: . newis - new IS on comm

1165:    Level: advanced

1167:    Notes:
1168:    It is usually desirable to create a parallel IS and look at the local part when necessary.

1170:    This function is useful if serial ISs must be created independently, or to view many
1171:    logically independent serial ISs.

1173:    The input IS must have the same type on every process.

1175: .seealso: ISSplit()
1176: @*/
1177: PetscErrorCode  ISOnComm(IS is,MPI_Comm comm,PetscCopyMode mode,IS *newis)
1178: {
1180:   PetscMPIInt    match;

1185:   MPI_Comm_compare(PetscObjectComm((PetscObject)is),comm,&match);
1186:   if (mode != PETSC_COPY_VALUES && (match == MPI_IDENT || match == MPI_CONGRUENT)) {
1187:     PetscObjectReference((PetscObject)is);
1188:     *newis = is;
1189:   } else {
1190:     (*is->ops->oncomm)(is,comm,mode,newis);
1191:   }
1192:   return(0);
1193: }

1195: /*@
1196:    ISSetBlockSize - informs an index set that it has a given block size

1198:    Logicall Collective on IS

1200:    Input Arguments:
1201: + is - index set
1202: - bs - block size

1204:    Level: intermediate

1206: .seealso: ISGetBlockSize(), ISCreateBlock()
1207: @*/
1208: PetscErrorCode  ISSetBlockSize(IS is,PetscInt bs)
1209: {

1215:   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)is),PETSC_ERR_ARG_OUTOFRANGE,"Block size %D, must be positive",bs);
1216:   (*is->ops->setblocksize)(is,bs);
1217:   return(0);
1218: }

1220: /*@
1221:    ISGetBlockSize - Returns the number of elements in a block.

1223:    Not Collective

1225:    Input Parameter:
1226: .  is - the index set

1228:    Output Parameter:
1229: .  size - the number of elements in a block

1231:    Level: intermediate

1233:    Concepts: IS^block size
1234:    Concepts: index sets^block size

1236: .seealso: ISBlockGetSize(), ISGetSize(), ISCreateBlock(), ISSetBlockSize()
1237: @*/
1238: PetscErrorCode  ISGetBlockSize(IS is,PetscInt *size)
1239: {

1243:   PetscLayoutGetBlockSize(is->map, size);
1244:   return(0);
1245: }

1247: PetscErrorCode ISGetIndicesCopy(IS is, PetscInt idx[])
1248: {
1250:   PetscInt       len,i;
1251:   const PetscInt *ptr;

1254:   ISGetSize(is,&len);
1255:   ISGetIndices(is,&ptr);
1256:   for (i=0; i<len; i++) idx[i] = ptr[i];
1257:   ISRestoreIndices(is,&ptr);
1258:   return(0);
1259: }

1261: /*MC
1262:     ISGetIndicesF90 - Accesses the elements of an index set from Fortran90.
1263:     The users should call ISRestoreIndicesF90() after having looked at the
1264:     indices.  The user should NOT change the indices.

1266:     Synopsis:
1267:     ISGetIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)

1269:     Not collective

1271:     Input Parameter:
1272: .   x - index set

1274:     Output Parameters:
1275: +   xx_v - the Fortran90 pointer to the array
1276: -   ierr - error code

1278:     Example of Usage:
1279: .vb
1280:     PetscInt, pointer xx_v(:)
1281:     ....
1282:     call ISGetIndicesF90(x,xx_v,ierr)
1283:     a = xx_v(3)
1284:     call ISRestoreIndicesF90(x,xx_v,ierr)
1285: .ve

1287:     Level: intermediate

1289: .seealso:  ISRestoreIndicesF90(), ISGetIndices(), ISRestoreIndices()

1291:   Concepts: index sets^getting indices in f90
1292:   Concepts: indices of index set in f90

1294: M*/

1296: /*MC
1297:     ISRestoreIndicesF90 - Restores an index set to a usable state after
1298:     a call to ISGetIndicesF90().

1300:     Synopsis:
1301:     ISRestoreIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)

1303:     Not collective

1305:     Input Parameters:
1306: .   x - index set
1307: .   xx_v - the Fortran90 pointer to the array

1309:     Output Parameter:
1310: .   ierr - error code


1313:     Example of Usage:
1314: .vb
1315:     PetscInt, pointer xx_v(:)
1316:     ....
1317:     call ISGetIndicesF90(x,xx_v,ierr)
1318:     a = xx_v(3)
1319:     call ISRestoreIndicesF90(x,xx_v,ierr)
1320: .ve

1322:     Level: intermediate

1324: .seealso:  ISGetIndicesF90(), ISGetIndices(), ISRestoreIndices()

1326: M*/

1328: /*MC
1329:     ISBlockGetIndicesF90 - Accesses the elements of an index set from Fortran90.
1330:     The users should call ISBlockRestoreIndicesF90() after having looked at the
1331:     indices.  The user should NOT change the indices.

1333:     Synopsis:
1334:     ISBlockGetIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)

1336:     Not collective

1338:     Input Parameter:
1339: .   x - index set

1341:     Output Parameters:
1342: +   xx_v - the Fortran90 pointer to the array
1343: -   ierr - error code
1344:     Example of Usage:
1345: .vb
1346:     PetscInt, pointer xx_v(:)
1347:     ....
1348:     call ISBlockGetIndicesF90(x,xx_v,ierr)
1349:     a = xx_v(3)
1350:     call ISBlockRestoreIndicesF90(x,xx_v,ierr)
1351: .ve

1353:     Level: intermediate

1355: .seealso:  ISBlockRestoreIndicesF90(), ISGetIndices(), ISRestoreIndices(),
1356:            ISRestoreIndices()

1358:   Concepts: index sets^getting block indices in f90
1359:   Concepts: indices of index set in f90
1360:   Concepts: block^ indices of index set in f90

1362: M*/

1364: /*MC
1365:     ISBlockRestoreIndicesF90 - Restores an index set to a usable state after
1366:     a call to ISBlockGetIndicesF90().

1368:     Synopsis:
1369:     ISBlockRestoreIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)

1371:     Not Collective

1373:     Input Parameters:
1374: +   x - index set
1375: -   xx_v - the Fortran90 pointer to the array

1377:     Output Parameter:
1378: .   ierr - error code

1380:     Example of Usage:
1381: .vb
1382:     PetscInt, pointer xx_v(:)
1383:     ....
1384:     call ISBlockGetIndicesF90(x,xx_v,ierr)
1385:     a = xx_v(3)
1386:     call ISBlockRestoreIndicesF90(x,xx_v,ierr)
1387: .ve

1389:     Notes:
1390:     Not yet supported for all F90 compilers

1392:     Level: intermediate

1394: .seealso:  ISBlockGetIndicesF90(), ISGetIndices(), ISRestoreIndices(), ISRestoreIndicesF90()

1396: M*/