Actual source code: index.c

petsc-master 2017-04-26
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: For parallel index sets this does the complete parallel permutation, but the
388:     code is not efficient for huge index sets (10,000,000 indices).

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

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

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

414:    Not Collective

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

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

422:    Level: beginner

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

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

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

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

442:    Not Collective

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

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

450:    Level: beginner

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

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

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

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

473:    Not Collective

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

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

481:    Fortran Note:
482:    This routine is used differently from Fortran
483: $    IS          is
484: $    integer     is_array(1)
485: $    PetscOffset i_is
486: $    int         ierr
487: $       call ISGetIndices(is,is_array,i_is,ierr)
488: $
489: $   Access first local entry in list
490: $      value = is_array(i_is + 1)
491: $
492: $      ...... other code
493: $       call ISRestoreIndices(is,is_array,i_is,ierr)

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

498:    Level: intermediate

500:    Concepts: index sets^getting indices
501:    Concepts: indices of index set

503: .seealso: ISRestoreIndices(), ISGetIndicesF90()
504: @*/
505: PetscErrorCode  ISGetIndices(IS is,const PetscInt *ptr[])
506: {

512:   (*is->ops->getindices)(is,ptr);
513:   return(0);
514: }

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

519:    Not Collective

521:    Input Parameter:
522: .  is - the index set

524:    Output Parameter:
525: +   min - the minimum value
526: -   max - the maximum value

528:    Level: intermediate

530:    Concepts: index sets^getting indices
531:    Concepts: indices of index set

533: .seealso: ISGetIndices(), ISRestoreIndices(), ISGetIndicesF90()
534: @*/
535: PetscErrorCode  ISGetMinMax(IS is,PetscInt *min,PetscInt *max)
536: {
539:   if (min) *min = is->min;
540:   if (max) *max = is->max;
541:   return(0);
542: }

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

547:   Not Collective

549:   Input Parameter:
550: + is - the index set
551: - key - the search key

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

556:   Level: intermediate
557: @*/
558: PetscErrorCode ISLocate(IS is, PetscInt key, PetscInt *location)
559: {

563:   if (is->ops->locate) {
564:     (*is->ops->locate)(is,key,location);
565:   } else {
566:     PetscInt       numIdx;
567:     PetscBool      sorted;
568:     const PetscInt *idx;

570:     ISGetLocalSize(is,&numIdx);
571:     ISGetIndices(is,&idx);
572:     ISSorted(is,&sorted);
573:     if (sorted) {
574:       PetscFindInt(key,numIdx,idx,location);
575:     } else {
576:       PetscInt i;

578:       *location = -1;
579:       for (i = 0; i < numIdx; i++) {
580:         if (idx[i] == key) {
581:           *location = i;
582:           break;
583:         }
584:       }
585:     }
586:     ISRestoreIndices(is,&idx);
587:   }
588:   return(0);
589: }

591: /*@C
592:    ISRestoreIndices - Restores an index set to a usable state after a call
593:                       to ISGetIndices().

595:    Not Collective

597:    Input Parameters:
598: +  is - the index set
599: -  ptr - the pointer obtained by ISGetIndices()

601:    Fortran Note:
602:    This routine is used differently from Fortran
603: $    IS          is
604: $    integer     is_array(1)
605: $    PetscOffset i_is
606: $    int         ierr
607: $       call ISGetIndices(is,is_array,i_is,ierr)
608: $
609: $   Access first local entry in list
610: $      value = is_array(i_is + 1)
611: $
612: $      ...... other code
613: $       call ISRestoreIndices(is,is_array,i_is,ierr)

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

618:    Level: intermediate

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

623: .seealso: ISGetIndices(), ISRestoreIndicesF90()
624: @*/
625: PetscErrorCode  ISRestoreIndices(IS is,const PetscInt *ptr[])
626: {

632:   if (is->ops->restoreindices) {
633:     (*is->ops->restoreindices)(is,ptr);
634:   }
635:   return(0);
636: }

638: static PetscErrorCode ISGatherTotal_Private(IS is)
639: {
641:   PetscInt       i,n,N;
642:   const PetscInt *lindices;
643:   MPI_Comm       comm;
644:   PetscMPIInt    rank,size,*sizes = NULL,*offsets = NULL,nn;


649:   PetscObjectGetComm((PetscObject)is,&comm);
650:   MPI_Comm_size(comm,&size);
651:   MPI_Comm_rank(comm,&rank);
652:   ISGetLocalSize(is,&n);
653:   PetscMalloc2(size,&sizes,size,&offsets);

655:   PetscMPIIntCast(n,&nn);
656:   MPI_Allgather(&nn,1,MPI_INT,sizes,1,MPI_INT,comm);
657:   offsets[0] = 0;
658:   for (i=1; i<size; ++i) offsets[i] = offsets[i-1] + sizes[i-1];
659:   N = offsets[size-1] + sizes[size-1];

661:   PetscMalloc1(N,&(is->total));
662:   ISGetIndices(is,&lindices);
663:   MPI_Allgatherv((void*)lindices,nn,MPIU_INT,is->total,sizes,offsets,MPIU_INT,comm);
664:   ISRestoreIndices(is,&lindices);
665:   is->local_offset = offsets[rank];
666:   PetscFree2(sizes,offsets);
667:   return(0);
668: }

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

673:    Collective on IS

675:    Input Parameter:
676: .  is - the index set

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

682:    Level: intermediate

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

691:    Concepts: index sets^getting nonlocal indices
692: .seealso: ISRestoreTotalIndices(), ISGetNonlocalIndices(), ISGetSize()
693: @*/
694: PetscErrorCode ISGetTotalIndices(IS is, const PetscInt *indices[])
695: {
697:   PetscMPIInt    size;

702:   MPI_Comm_size(PetscObjectComm((PetscObject)is), &size);
703:   if (size == 1) {
704:     (*is->ops->getindices)(is,indices);
705:   } else {
706:     if (!is->total) {
707:       ISGatherTotal_Private(is);
708:     }
709:     *indices = is->total;
710:   }
711:   return(0);
712: }

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

717:    Not Collective.

719:    Input Parameter:
720: +  is - the index set
721: -  indices - index array; must be the array obtained with ISGetTotalIndices()

723:    Level: intermediate

725:    Concepts: index sets^getting nonlocal indices
726:    Concepts: index sets^restoring nonlocal indices
727: .seealso: ISRestoreTotalIndices(), ISGetNonlocalIndices()
728: @*/
729: PetscErrorCode  ISRestoreTotalIndices(IS is, const PetscInt *indices[])
730: {
732:   PetscMPIInt    size;

737:   MPI_Comm_size(PetscObjectComm((PetscObject)is), &size);
738:   if (size == 1) {
739:     (*is->ops->restoreindices)(is,indices);
740:   } else {
741:     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.");
742:   }
743:   return(0);
744: }
745: /*@C
746:    ISGetNonlocalIndices - Retrieve an array of indices from remote processors
747:                        in this communicator.

749:    Collective on IS

751:    Input Parameter:
752: .  is - the index set

754:    Output Parameter:
755: .  indices - indices with rank 0 indices first, and so on,  omitting
756:              the current rank.  Total number of indices is the difference
757:              total and local, obtained with ISGetSize() and ISGetLocalSize(),
758:              respectively.

760:    Level: intermediate

762:    Notes: restore the indices using ISRestoreNonlocalIndices().
763:           The same scalability considerations as those for ISGetTotalIndices
764:           apply here.

766:    Concepts: index sets^getting nonlocal indices
767: .seealso: ISGetTotalIndices(), ISRestoreNonlocalIndices(), ISGetSize(), ISGetLocalSize().
768: @*/
769: PetscErrorCode  ISGetNonlocalIndices(IS is, const PetscInt *indices[])
770: {
772:   PetscMPIInt    size;
773:   PetscInt       n, N;

778:   MPI_Comm_size(PetscObjectComm((PetscObject)is), &size);
779:   if (size == 1) *indices = NULL;
780:   else {
781:     if (!is->total) {
782:       ISGatherTotal_Private(is);
783:     }
784:     ISGetLocalSize(is,&n);
785:     ISGetSize(is,&N);
786:     PetscMalloc1(N-n, &(is->nonlocal));
787:     PetscMemcpy(is->nonlocal, is->total, sizeof(PetscInt)*is->local_offset);
788:     PetscMemcpy(is->nonlocal+is->local_offset, is->total+is->local_offset+n, sizeof(PetscInt)*(N - is->local_offset - n));
789:     *indices = is->nonlocal;
790:   }
791:   return(0);
792: }

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

797:    Not Collective.

799:    Input Parameter:
800: +  is - the index set
801: -  indices - index array; must be the array obtained with ISGetNonlocalIndices()

803:    Level: intermediate

805:    Concepts: index sets^getting nonlocal indices
806:    Concepts: index sets^restoring nonlocal indices
807: .seealso: ISGetTotalIndices(), ISGetNonlocalIndices(), ISRestoreTotalIndices()
808: @*/
809: PetscErrorCode  ISRestoreNonlocalIndices(IS is, const PetscInt *indices[])
810: {
814:   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.");
815:   return(0);
816: }

818: /*@
819:    ISGetNonlocalIS - Gather all nonlocal indices for this IS and present
820:                      them as another sequential index set.


823:    Collective on IS

825:    Input Parameter:
826: .  is - the index set

828:    Output Parameter:
829: .  complement - sequential IS with indices identical to the result of
830:                 ISGetNonlocalIndices()

832:    Level: intermediate

834:    Notes: complement represents the result of ISGetNonlocalIndices as an IS.
835:           Therefore scalability issues similar to ISGetNonlocalIndices apply.
836:           The resulting IS must be restored using ISRestoreNonlocalIS().

838:    Concepts: index sets^getting nonlocal indices
839: .seealso: ISGetNonlocalIndices(), ISRestoreNonlocalIndices(),  ISAllGather(), ISGetSize()
840: @*/
841: PetscErrorCode  ISGetNonlocalIS(IS is, IS *complement)
842: {

848:   /* Check if the complement exists already. */
849:   if (is->complement) {
850:     *complement = is->complement;
851:     PetscObjectReference((PetscObject)(is->complement));
852:   } else {
853:     PetscInt       N, n;
854:     const PetscInt *idx;
855:     ISGetSize(is, &N);
856:     ISGetLocalSize(is,&n);
857:     ISGetNonlocalIndices(is, &idx);
858:     ISCreateGeneral(PETSC_COMM_SELF, N-n,idx, PETSC_USE_POINTER, &(is->complement));
859:     PetscObjectReference((PetscObject)is->complement);
860:     *complement = is->complement;
861:   }
862:   return(0);
863: }


866: /*@
867:    ISRestoreNonlocalIS - Restore the IS obtained with ISGetNonlocalIS().

869:    Not collective.

871:    Input Parameter:
872: +  is         - the index set
873: -  complement - index set of is's nonlocal indices

875:    Level: intermediate


878:    Concepts: index sets^getting nonlocal indices
879:    Concepts: index sets^restoring nonlocal indices
880: .seealso: ISGetNonlocalIS(), ISGetNonlocalIndices(), ISRestoreNonlocalIndices()
881: @*/
882: PetscErrorCode  ISRestoreNonlocalIS(IS is, IS *complement)
883: {
885:   PetscInt       refcnt;

890:   if (*complement != is->complement) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Complement IS being restored was not obtained with ISGetNonlocalIS()");
891:   PetscObjectGetReference((PetscObject)(is->complement), &refcnt);
892:   if (refcnt <= 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Duplicate call to ISRestoreNonlocalIS() detected");
893:   PetscObjectDereference((PetscObject)(is->complement));
894:   return(0);
895: }

897: /*@C
898:    ISView - Displays an index set.

900:    Collective on IS

902:    Input Parameters:
903: +  is - the index set
904: -  viewer - viewer used to display the set, for example PETSC_VIEWER_STDOUT_SELF.

906:    Level: intermediate

908: .seealso: PetscViewerASCIIOpen()
909: @*/
910: PetscErrorCode  ISView(IS is,PetscViewer viewer)
911: {

916:   if (!viewer) {
917:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)is),&viewer);
918:   }

922:   PetscObjectPrintClassNamePrefixType((PetscObject)is,viewer);
923:   (*is->ops->view)(is,viewer);
924:   return(0);
925: }

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

930:   Collective on PetscViewer

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

936:   Level: intermediate

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

943:   Concepts: index set^loading from file

945: .seealso: PetscViewerBinaryOpen(), ISView(), MatLoad(), VecLoad()
946: @*/
947: PetscErrorCode ISLoad(IS is, PetscViewer viewer)
948: {
949:   PetscBool      isbinary, ishdf5;

955:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERBINARY, &isbinary);
956:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
957:   if (!isbinary && !ishdf5) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
958:   if (!((PetscObject)is)->type_name) {ISSetType(is, ISGENERAL);}
959:   (*is->ops->load)(is, viewer);
960:   return(0);
961: }

963: /*@
964:    ISSort - Sorts the indices of an index set.

966:    Collective on IS

968:    Input Parameters:
969: .  is - the index set

971:    Level: intermediate

973:    Concepts: index sets^sorting
974:    Concepts: sorting^index set

976: .seealso: ISSortRemoveDups(), ISSorted()
977: @*/
978: PetscErrorCode  ISSort(IS is)
979: {

984:   (*is->ops->sort)(is);
985:   return(0);
986: }

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

991:   Collective on IS

993:   Input Parameters:
994: . is - the index set

996:   Level: intermediate

998:   Concepts: index sets^sorting
999:   Concepts: sorting^index set

1001: .seealso: ISSort(), ISSorted()
1002: @*/
1003: PetscErrorCode ISSortRemoveDups(IS is)
1004: {

1009:   (*is->ops->sortremovedups)(is);
1010:   return(0);
1011: }

1013: /*@
1014:    ISToGeneral - Converts an IS object of any type to ISGENERAL type

1016:    Collective on IS

1018:    Input Parameters:
1019: .  is - the index set

1021:    Level: intermediate

1023:    Concepts: index sets^sorting
1024:    Concepts: sorting^index set

1026: .seealso: ISSorted()
1027: @*/
1028: PetscErrorCode  ISToGeneral(IS is)
1029: {

1034:   if (is->ops->togeneral) {
1035:     (*is->ops->togeneral)(is);
1036:   } else SETERRQ1(PetscObjectComm((PetscObject)is),PETSC_ERR_SUP,"Not written for this type %s",((PetscObject)is)->type_name);
1037:   return(0);
1038: }

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

1043:    Collective on IS

1045:    Input Parameter:
1046: .  is - the index set

1048:    Output Parameter:
1049: .  flg - output flag, either PETSC_TRUE if the index set is sorted,
1050:          or PETSC_FALSE otherwise.

1052:    Notes: For parallel IS objects this only indicates if the local part of the IS
1053:           is sorted. So some processors may return PETSC_TRUE while others may
1054:           return PETSC_FALSE.

1056:    Level: intermediate

1058: .seealso: ISSort(), ISSortRemoveDups()
1059: @*/
1060: PetscErrorCode  ISSorted(IS is,PetscBool  *flg)
1061: {

1067:   (*is->ops->sorted)(is,flg);
1068:   return(0);
1069: }

1071: /*@
1072:    ISDuplicate - Creates a duplicate copy of an index set.

1074:    Collective on IS

1076:    Input Parmeters:
1077: .  is - the index set

1079:    Output Parameters:
1080: .  isnew - the copy of the index set

1082:    Level: beginner

1084:    Concepts: index sets^duplicating

1086: .seealso: ISCreateGeneral(), ISCopy()
1087: @*/
1088: PetscErrorCode  ISDuplicate(IS is,IS *newIS)
1089: {

1095:   (*is->ops->duplicate)(is,newIS);
1096:   (*newIS)->isidentity = is->isidentity;
1097:   (*newIS)->isperm     = is->isperm;
1098:   return(0);
1099: }

1101: /*@
1102:    ISCopy - Copies an index set.

1104:    Collective on IS

1106:    Input Parmeters:
1107: .  is - the index set

1109:    Output Parameters:
1110: .  isy - the copy of the index set

1112:    Level: beginner

1114:    Concepts: index sets^copying

1116: .seealso: ISDuplicate()
1117: @*/
1118: PetscErrorCode  ISCopy(IS is,IS isy)
1119: {

1126:   if (is == isy) return(0);
1127:   (*is->ops->copy)(is,isy);
1128:   isy->isperm     = is->isperm;
1129:   isy->max        = is->max;
1130:   isy->min        = is->min;
1131:   isy->isidentity = is->isidentity;
1132:   return(0);
1133: }

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

1138:    Collective on IS and comm

1140:    Input Arguments:
1141: + is - index set
1142: . comm - communicator for new index set
1143: - mode - copy semantics, PETSC_USE_POINTER for no-copy if possible, otherwise PETSC_COPY_VALUES

1145:    Output Arguments:
1146: . newis - new IS on comm

1148:    Level: advanced

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

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

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

1158: .seealso: ISSplit()
1159: @*/
1160: PetscErrorCode  ISOnComm(IS is,MPI_Comm comm,PetscCopyMode mode,IS *newis)
1161: {
1163:   PetscMPIInt    match;

1168:   MPI_Comm_compare(PetscObjectComm((PetscObject)is),comm,&match);
1169:   if (mode != PETSC_COPY_VALUES && (match == MPI_IDENT || match == MPI_CONGRUENT)) {
1170:     PetscObjectReference((PetscObject)is);
1171:     *newis = is;
1172:   } else {
1173:     (*is->ops->oncomm)(is,comm,mode,newis);
1174:   }
1175:   return(0);
1176: }

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

1181:    Logicall Collective on IS

1183:    Input Arguments:
1184: + is - index set
1185: - bs - block size

1187:    Level: intermediate

1189: .seealso: ISGetBlockSize(), ISCreateBlock()
1190: @*/
1191: PetscErrorCode  ISSetBlockSize(IS is,PetscInt bs)
1192: {

1198:   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)is),PETSC_ERR_ARG_OUTOFRANGE,"Block size %D, must be positive",bs);
1199:   (*is->ops->setblocksize)(is,bs);
1200:   return(0);
1201: }

1203: /*@
1204:    ISGetBlockSize - Returns the number of elements in a block.

1206:    Not Collective

1208:    Input Parameter:
1209: .  is - the index set

1211:    Output Parameter:
1212: .  size - the number of elements in a block

1214:    Level: intermediate

1216:    Concepts: IS^block size
1217:    Concepts: index sets^block size

1219: .seealso: ISBlockGetSize(), ISGetSize(), ISCreateBlock(), ISSetBlockSize()
1220: @*/
1221: PetscErrorCode  ISGetBlockSize(IS is,PetscInt *size)
1222: {

1226:   PetscLayoutGetBlockSize(is->map, size);
1227:   return(0);
1228: }

1230: PetscErrorCode ISGetIndicesCopy(IS is, PetscInt idx[])
1231: {
1233:   PetscInt       len,i;
1234:   const PetscInt *ptr;

1237:   ISGetSize(is,&len);
1238:   ISGetIndices(is,&ptr);
1239:   for (i=0; i<len; i++) idx[i] = ptr[i];
1240:   ISRestoreIndices(is,&ptr);
1241:   return(0);
1242: }

1244: /*MC
1245:     ISGetIndicesF90 - Accesses the elements of an index set from Fortran90.
1246:     The users should call ISRestoreIndicesF90() after having looked at the
1247:     indices.  The user should NOT change the indices.

1249:     Synopsis:
1250:     ISGetIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)

1252:     Not collective

1254:     Input Parameter:
1255: .   x - index set

1257:     Output Parameters:
1258: +   xx_v - the Fortran90 pointer to the array
1259: -   ierr - error code

1261:     Example of Usage:
1262: .vb
1263:     PetscScalar, pointer xx_v(:)
1264:     ....
1265:     call ISGetIndicesF90(x,xx_v,ierr)
1266:     a = xx_v(3)
1267:     call ISRestoreIndicesF90(x,xx_v,ierr)
1268: .ve

1270:     Notes:
1271:     Not yet supported for all F90 compilers.

1273:     Level: intermediate

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

1277:   Concepts: index sets^getting indices in f90
1278:   Concepts: indices of index set in f90

1280: M*/

1282: /*MC
1283:     ISRestoreIndicesF90 - Restores an index set to a usable state after
1284:     a call to ISGetIndicesF90().

1286:     Synopsis:
1287:     ISRestoreIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)

1289:     Not collective

1291:     Input Parameters:
1292: .   x - index set
1293: .   xx_v - the Fortran90 pointer to the array

1295:     Output Parameter:
1296: .   ierr - error code


1299:     Example of Usage:
1300: .vb
1301:     PetscScalar, pointer xx_v(:)
1302:     ....
1303:     call ISGetIndicesF90(x,xx_v,ierr)
1304:     a = xx_v(3)
1305:     call ISRestoreIndicesF90(x,xx_v,ierr)
1306: .ve

1308:     Notes:
1309:     Not yet supported for all F90 compilers.

1311:     Level: intermediate

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

1315: M*/

1317: /*MC
1318:     ISBlockGetIndicesF90 - Accesses the elements of an index set from Fortran90.
1319:     The users should call ISBlockRestoreIndicesF90() after having looked at the
1320:     indices.  The user should NOT change the indices.

1322:     Synopsis:
1323:     ISBlockGetIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)

1325:     Not collective

1327:     Input Parameter:
1328: .   x - index set

1330:     Output Parameters:
1331: +   xx_v - the Fortran90 pointer to the array
1332: -   ierr - error code
1333:     Example of Usage:
1334: .vb
1335:     PetscScalar, pointer xx_v(:)
1336:     ....
1337:     call ISBlockGetIndicesF90(x,xx_v,ierr)
1338:     a = xx_v(3)
1339:     call ISBlockRestoreIndicesF90(x,xx_v,ierr)
1340: .ve

1342:     Notes:
1343:     Not yet supported for all F90 compilers

1345:     Level: intermediate

1347: .seealso:  ISBlockRestoreIndicesF90(), ISGetIndices(), ISRestoreIndices(),
1348:            ISRestoreIndices()

1350:   Concepts: index sets^getting block indices in f90
1351:   Concepts: indices of index set in f90
1352:   Concepts: block^ indices of index set in f90

1354: M*/

1356: /*MC
1357:     ISBlockRestoreIndicesF90 - Restores an index set to a usable state after
1358:     a call to ISBlockGetIndicesF90().

1360:     Synopsis:
1361:     ISBlockRestoreIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)

1363:     Not Collective

1365:     Input Parameters:
1366: +   x - index set
1367: -   xx_v - the Fortran90 pointer to the array

1369:     Output Parameter:
1370: .   ierr - error code

1372:     Example of Usage:
1373: .vb
1374:     PetscScalar, pointer xx_v(:)
1375:     ....
1376:     call ISBlockGetIndicesF90(x,xx_v,ierr)
1377:     a = xx_v(3)
1378:     call ISBlockRestoreIndicesF90(x,xx_v,ierr)
1379: .ve

1381:     Notes:
1382:     Not yet supported for all F90 compilers

1384:     Level: intermediate

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

1388: M*/