Actual source code: index.c

petsc-3.3-p7 2013-05-11
  2: /*  
  3:    Defines the abstract operations on index sets, i.e. the public interface. 
  4: */
  5: #include <petsc-private/isimpl.h>      /*I "petscis.h" I*/

  7: /* Logging support */
  8: PetscClassId  IS_CLASSID;

 12: /*@
 13:    ISIdentity - Determines whether index set is the identity mapping.

 15:    Collective on IS

 17:    Input Parmeters:
 18: .  is - the index set

 20:    Output Parameters:
 21: .  ident - PETSC_TRUE if an identity, else PETSC_FALSE

 23:    Level: intermediate

 25:    Concepts: identity mapping
 26:    Concepts: index sets^is identity

 28: .seealso: ISSetIdentity()
 29: @*/
 30: PetscErrorCode  ISIdentity(IS is,PetscBool  *ident)
 31: {

 37:   *ident = is->isidentity;
 38:   if (*ident) return(0);
 39:   if (is->ops->identity) {
 40:     (*is->ops->identity)(is,ident);
 41:   }
 42:   return(0);
 43: }

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

 50:    Logically Collective on IS

 52:    Input Parmeters:
 53: .  is - the index set

 55:    Level: intermediate

 57:    Concepts: identity mapping
 58:    Concepts: index sets^is identity

 60: .seealso: ISIdentity()
 61: @*/
 62: PetscErrorCode  ISSetIdentity(IS is)
 63: {
 66:   is->isidentity = PETSC_TRUE;
 67:   return(0);
 68: }

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

 75:    Not Collective

 77:    Input Parmeters:
 78: +  is - the index set
 79: .  gstart - global start
 80: .  gend - global end

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

 86:    Level: developer

 88:    Concepts: index sets^is contiguous

 90: .seealso: ISGetLocalSize(), VecGetOwnershipRange()
 91: @*/
 92: PetscErrorCode  ISContiguousLocal(IS is,PetscInt gstart,PetscInt gend,PetscInt *start,PetscBool *contig)
 93: {

100:   if (is->ops->contiguous) {
101:     (*is->ops->contiguous)(is,gstart,gend,start,contig);
102:   } else {
103:     *start = -1;
104:     *contig = PETSC_FALSE;
105:   }
106:   return(0);
107: }

111: /*@
112:    ISPermutation - PETSC_TRUE or PETSC_FALSE depending on whether the 
113:    index set has been declared to be a permutation.

115:    Logically Collective on IS

117:    Input Parmeters:
118: .  is - the index set

120:    Output Parameters:
121: .  perm - PETSC_TRUE if a permutation, else PETSC_FALSE

123:    Level: intermediate

125:   Concepts: permutation
126:   Concepts: index sets^is permutation

128: .seealso: ISSetPermutation()
129: @*/
130: PetscErrorCode  ISPermutation(IS is,PetscBool  *perm)
131: {
135:   *perm = (PetscBool) is->isperm;
136:   return(0);
137: }

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

144:    Logically Collective on IS

146:    Input Parmeters:
147: .  is - the index set

149:    Level: intermediate

151:   Concepts: permutation
152:   Concepts: index sets^permutation

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

157: .seealso: ISPermutation()
158: @*/
159: PetscErrorCode  ISSetPermutation(IS is)
160: {
163: #if defined(PETSC_USE_DEBUG)
164:   {
165:     PetscMPIInt    size;

168:     MPI_Comm_size(((PetscObject)is)->comm,&size);
169:     if (size == 1) {
170:       PetscInt       i,n,*idx;
171:       const PetscInt *iidx;
172: 
173:       ISGetSize(is,&n);
174:       PetscMalloc(n*sizeof(PetscInt),&idx);
175:       ISGetIndices(is,&iidx);
176:       PetscMemcpy(idx,iidx,n*sizeof(PetscInt));
177:       PetscSortInt(n,idx);
178:       for (i=0; i<n; i++) {
179:         if (idx[i] != i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Index set is not a permutation");
180:       }
181:       PetscFree(idx);
182:       ISRestoreIndices(is,&iidx);
183:     }
184:   }
185: #endif
186:   is->isperm = PETSC_TRUE;
187:   return(0);
188: }

192: /*@
193:    ISDestroy - Destroys an index set.

195:    Collective on IS

197:    Input Parameters:
198: .  is - the index set

200:    Level: beginner

202: .seealso: ISCreateGeneral(), ISCreateStride(), ISCreateBlocked()
203: @*/
204: PetscErrorCode  ISDestroy(IS *is)
205: {

209:   if (!*is) return(0);
211:   if (--((PetscObject)(*is))->refct > 0) {*is = 0; return(0);}
212:   if ((*is)->complement) {
213:     PetscInt refcnt;
214:     PetscObjectGetReference((PetscObject)((*is)->complement), &refcnt);
215:     if (refcnt > 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Nonlocal IS has not been restored");
216:     ISDestroy(&(*is)->complement);
217:   }
218:   if ((*is)->ops->destroy) {
219:     (*(*is)->ops->destroy)(*is);
220:   }
221:   /* Destroy local representations of offproc data. */
222:   PetscFree((*is)->total);
223:   PetscFree((*is)->nonlocal);
224:   PetscHeaderDestroy(is);
225:   return(0);
226: }

230: /*@
231:    ISInvertPermutation - Creates a new permutation that is the inverse of 
232:                          a given permutation.

234:    Collective on IS

236:    Input Parameter:
237: +  is - the index set
238: -  nlocal - number of indices on this processor in result (ignored for 1 proccessor) or
239:             use PETSC_DECIDE

241:    Output Parameter:
242: .  isout - the inverse permutation

244:    Level: intermediate

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

249:    Concepts: inverse permutation
250:    Concepts: permutation^inverse
251:    Concepts: index sets^inverting
252: @*/
253: PetscErrorCode  ISInvertPermutation(IS is,PetscInt nlocal,IS *isout)
254: {

260:   if (!is->isperm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Not a permutation, must call ISSetPermutation() on the IS first");
261:   (*is->ops->invertpermutation)(is,nlocal,isout);
262:   ISSetPermutation(*isout);
263:   return(0);
264: }

268: /*@
269:    ISGetSize - Returns the global length of an index set. 

271:    Not Collective

273:    Input Parameter:
274: .  is - the index set

276:    Output Parameter:
277: .  size - the global size

279:    Level: beginner

281:    Concepts: size^of index set
282:    Concepts: index sets^size

284: @*/
285: PetscErrorCode  ISGetSize(IS is,PetscInt *size)
286: {

292:   (*is->ops->getsize)(is,size);
293:   return(0);
294: }

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

301:    Not Collective

303:    Input Parameter:
304: .  is - the index set

306:    Output Parameter:
307: .  size - the local size

309:    Level: beginner

311:    Concepts: size^of index set
312:    Concepts: local size^of index set
313:    Concepts: index sets^local size
314:   
315: @*/
316: PetscErrorCode  ISGetLocalSize(IS is,PetscInt *size)
317: {

323:   (*is->ops->getlocalsize)(is,size);
324:   return(0);
325: }

329: /*@C
330:    ISGetIndices - Returns a pointer to the indices.  The user should call 
331:    ISRestoreIndices() after having looked at the indices.  The user should 
332:    NOT change the indices.

334:    Not Collective

336:    Input Parameter:
337: .  is - the index set

339:    Output Parameter:
340: .  ptr - the location to put the pointer to the indices

342:    Fortran Note:
343:    This routine is used differently from Fortran
344: $    IS          is
345: $    integer     is_array(1)
346: $    PetscOffset i_is
347: $    int         ierr
348: $       call ISGetIndices(is,is_array,i_is,ierr)
349: $
350: $   Access first local entry in list
351: $      value = is_array(i_is + 1)
352: $
353: $      ...... other code
354: $       call ISRestoreIndices(is,is_array,i_is,ierr)

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

359:    Level: intermediate

361:    Concepts: index sets^getting indices
362:    Concepts: indices of index set

364: .seealso: ISRestoreIndices(), ISGetIndicesF90()
365: @*/
366: PetscErrorCode  ISGetIndices(IS is,const PetscInt *ptr[])
367: {

373:   (*is->ops->getindices)(is,ptr);
374:   return(0);
375: }

379: /*@C
380:    ISRestoreIndices - Restores an index set to a usable state after a call 
381:                       to ISGetIndices().

383:    Not Collective

385:    Input Parameters:
386: +  is - the index set
387: -  ptr - the pointer obtained by ISGetIndices()

389:    Fortran Note:
390:    This routine is used differently from Fortran
391: $    IS          is
392: $    integer     is_array(1)
393: $    PetscOffset i_is
394: $    int         ierr
395: $       call ISGetIndices(is,is_array,i_is,ierr)
396: $
397: $   Access first local entry in list
398: $      value = is_array(i_is + 1)
399: $
400: $      ...... other code
401: $       call ISRestoreIndices(is,is_array,i_is,ierr)

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

406:    Level: intermediate

408: .seealso: ISGetIndices(), ISRestoreIndicesF90()
409: @*/
410: PetscErrorCode  ISRestoreIndices(IS is,const PetscInt *ptr[])
411: {

417:   if (is->ops->restoreindices) {
418:     (*is->ops->restoreindices)(is,ptr);
419:   }
420:   return(0);
421: }

425: static PetscErrorCode ISGatherTotal_Private(IS is)
426: {
428:   PetscInt       i,n,N;
429:   const PetscInt *lindices;
430:   MPI_Comm       comm;
431:   PetscMPIInt    rank,size,*sizes = PETSC_NULL,*offsets = PETSC_NULL,nn;


436:   PetscObjectGetComm((PetscObject)is,&comm);
437:   MPI_Comm_size(comm,&size);
438:   MPI_Comm_rank(comm,&rank);
439:   ISGetLocalSize(is,&n);
440:   PetscMalloc2(size,PetscMPIInt,&sizes,size,PetscMPIInt,&offsets);
441: 
442:   nn   = PetscMPIIntCast(n);
443:   MPI_Allgather(&nn,1,MPIU_INT,sizes,1,MPIU_INT,comm);
444:   offsets[0] = 0;
445:   for (i=1;i<size; ++i) offsets[i] = offsets[i-1] + sizes[i-1];
446:   N = offsets[size-1] + sizes[size-1];
447: 
448:   PetscMalloc(N*sizeof(PetscInt),&(is->total));
449:   ISGetIndices(is,&lindices);
450:   MPI_Allgatherv((void*)lindices,nn,MPIU_INT,is->total,sizes,offsets,MPIU_INT,comm);
451:   ISRestoreIndices(is,&lindices);
452:   is->local_offset = offsets[rank];
453:   PetscFree2(sizes,offsets);

455:   return(0);
456: }

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

463:    Collective on IS

465:    Input Parameter:
466: .  is - the index set

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

472:    Level: intermediate

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

481:    Concepts: index sets^getting nonlocal indices
482: .seealso: ISRestoreTotalIndices(), ISGetNonlocalIndices(), ISGetSize()
483: @*/
484: PetscErrorCode ISGetTotalIndices(IS is, const PetscInt *indices[])
485: {
487:   PetscMPIInt    size;

492:   MPI_Comm_size(((PetscObject)is)->comm, &size);
493:   if(size == 1) {
494:     (*is->ops->getindices)(is,indices);
495:   }
496:   else {
497:     if(!is->total) {
498:       ISGatherTotal_Private(is);
499:     }
500:     *indices = is->total;
501:   }
502:   return(0);
503: }

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

510:    Not Collective.

512:    Input Parameter:
513: +  is - the index set
514: -  indices - index array; must be the array obtained with ISGetTotalIndices()

516:    Level: intermediate

518:    Concepts: index sets^getting nonlocal indices
519:    Concepts: index sets^restoring nonlocal indices
520: .seealso: ISRestoreTotalIndices(), ISGetNonlocalIndices()
521: @*/
522: PetscErrorCode  ISRestoreTotalIndices(IS is, const PetscInt *indices[])
523: {
525:   PetscMPIInt size;
529:   MPI_Comm_size(((PetscObject)is)->comm, &size);
530:   if(size == 1) {
531:     (*is->ops->restoreindices)(is,indices);
532:   }
533:   else {
534:     if(is->total != *indices) {
535:       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Index array pointer being restored does not point to the array obtained from the IS.");
536:     }
537:   }
538:   return(0);
539: }
542: /*@C
543:    ISGetNonlocalIndices - Retrieve an array of indices from remote processors
544:                        in this communicator.

546:    Collective on IS

548:    Input Parameter:
549: .  is - the index set

551:    Output Parameter:
552: .  indices - indices with rank 0 indices first, and so on,  omitting 
553:              the current rank.  Total number of indices is the difference
554:              total and local, obtained with ISGetSize() and ISGetLocalSize(),
555:              respectively.

557:    Level: intermediate

559:    Notes: restore the indices using ISRestoreNonlocalIndices().   
560:           The same scalability considerations as those for ISGetTotalIndices
561:           apply here.

563:    Concepts: index sets^getting nonlocal indices
564: .seealso: ISGetTotalIndices(), ISRestoreNonlocalIndices(), ISGetSize(), ISGetLocalSize().
565: @*/
566: PetscErrorCode  ISGetNonlocalIndices(IS is, const PetscInt *indices[])
567: {
569:   PetscMPIInt size;
570:   PetscInt n, N;
574:   MPI_Comm_size(((PetscObject)is)->comm, &size);
575:   if(size == 1) {
576:       *indices = PETSC_NULL;
577:   }
578:   else {
579:     if(!is->total) {
580:       ISGatherTotal_Private(is);
581:     }
582:     ISGetLocalSize(is,&n);
583:     ISGetSize(is,&N);
584:     PetscMalloc(sizeof(PetscInt)*(N-n), &(is->nonlocal));
585:     PetscMemcpy(is->nonlocal, is->total, sizeof(PetscInt)*is->local_offset);
586:     PetscMemcpy(is->nonlocal+is->local_offset, is->total+is->local_offset+n, sizeof(PetscInt)*(N - is->local_offset - n));
587:     *indices = is->nonlocal;
588:   }
589:   return(0);
590: }

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

597:    Not Collective.

599:    Input Parameter:
600: +  is - the index set
601: -  indices - index array; must be the array obtained with ISGetNonlocalIndices()

603:    Level: intermediate

605:    Concepts: index sets^getting nonlocal indices
606:    Concepts: index sets^restoring nonlocal indices
607: .seealso: ISGetTotalIndices(), ISGetNonlocalIndices(), ISRestoreTotalIndices()
608: @*/
609: PetscErrorCode  ISRestoreNonlocalIndices(IS is, const PetscInt *indices[])
610: {

615:   if(is->nonlocal != *indices) {
616:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Index array pointer being restored does not point to the array obtained from the IS.");
617:   }
618:   return(0);
619: }

623: /*@
624:    ISGetNonlocalIS - Gather all nonlocal indices for this IS and present 
625:                      them as another sequential index set.  


628:    Collective on IS

630:    Input Parameter:
631: .  is - the index set

633:    Output Parameter:
634: .  complement - sequential IS with indices identical to the result of
635:                 ISGetNonlocalIndices()

637:    Level: intermediate

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

643:    Concepts: index sets^getting nonlocal indices
644: .seealso: ISGetNonlocalIndices(), ISRestoreNonlocalIndices(),  ISAllGather(), ISGetSize()
645: @*/
646: PetscErrorCode  ISGetNonlocalIS(IS is, IS *complement)
647: {

653:   /* Check if the complement exists already. */
654:   if(is->complement) {
655:     *complement = is->complement;
656:     PetscObjectReference((PetscObject)(is->complement));
657:   }
658:   else {
659:     PetscInt       N, n;
660:     const PetscInt *idx;
661:     ISGetSize(is, &N);
662:     ISGetLocalSize(is,&n);
663:     ISGetNonlocalIndices(is, &idx);
664:     ISCreateGeneral(PETSC_COMM_SELF, N-n,idx, PETSC_USE_POINTER, &(is->complement));
665:     PetscObjectReference((PetscObject)is->complement);
666:     *complement = is->complement;
667:   }
668:   return(0);
669: }


674: /*@
675:    ISRestoreNonlocalIS - Restore the IS obtained with ISGetNonlocalIS().

677:    Not collective.

679:    Input Parameter:
680: +  is         - the index set
681: -  complement - index set of is's nonlocal indices

683:    Level: intermediate


686:    Concepts: index sets^getting nonlocal indices
687:    Concepts: index sets^restoring nonlocal indices
688: .seealso: ISGetNonlocalIS(), ISGetNonlocalIndices(), ISRestoreNonlocalIndices()
689: @*/
690: PetscErrorCode  ISRestoreNonlocalIS(IS is, IS *complement)
691: {
693:   PetscInt       refcnt;

698:   if(*complement != is->complement) {
699:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Complement IS being restored was not obtained with ISGetNonlocalIS()");
700:   }
701:   PetscObjectGetReference((PetscObject)(is->complement), &refcnt);
702:   if(refcnt <= 1) {
703:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Duplicate call to ISRestoreNonlocalIS() detected");
704:   }
705:   PetscObjectDereference((PetscObject)(is->complement));
706:   return(0);
707: }

711: /*@C
712:    ISView - Displays an index set.

714:    Collective on IS

716:    Input Parameters:
717: +  is - the index set
718: -  viewer - viewer used to display the set, for example PETSC_VIEWER_STDOUT_SELF.

720:    Level: intermediate

722: .seealso: PetscViewerASCIIOpen()
723: @*/
724: PetscErrorCode  ISView(IS is,PetscViewer viewer)
725: {

730:   if (!viewer) {
731:     PetscViewerASCIIGetStdout(((PetscObject)is)->comm,&viewer);
732:   }
735: 
736:   (*is->ops->view)(is,viewer);
737:   return(0);
738: }

742: /*@
743:    ISSort - Sorts the indices of an index set.

745:    Collective on IS

747:    Input Parameters:
748: .  is - the index set

750:    Level: intermediate

752:    Concepts: index sets^sorting
753:    Concepts: sorting^index set

755: .seealso: ISSorted()
756: @*/
757: PetscErrorCode  ISSort(IS is)
758: {

763:   (*is->ops->sort)(is);
764:   return(0);
765: }

769: /*@
770:    ISToGeneral - Converts an IS object of any type to ISGENERAL type

772:    Collective on IS

774:    Input Parameters:
775: .  is - the index set

777:    Level: intermediate

779:    Concepts: index sets^sorting
780:    Concepts: sorting^index set

782: .seealso: ISSorted()
783: @*/
784: PetscErrorCode  ISToGeneral(IS is)
785: {

790:   if (is->ops->togeneral) {
791:     (*is->ops->togeneral)(is);
792:   } else SETERRQ1(((PetscObject)is)->comm,PETSC_ERR_SUP,"Not written for this type %s",((PetscObject)is)->type_name);
793:   return(0);
794: }

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

801:    Collective on IS

803:    Input Parameter:
804: .  is - the index set

806:    Output Parameter:
807: .  flg - output flag, either PETSC_TRUE if the index set is sorted, 
808:          or PETSC_FALSE otherwise.

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

814:    Level: intermediate

816: .seealso: ISSort()
817: @*/
818: PetscErrorCode  ISSorted(IS is,PetscBool  *flg)
819: {

825:   (*is->ops->sorted)(is,flg);
826:   return(0);
827: }

831: /*@
832:    ISDuplicate - Creates a duplicate copy of an index set.

834:    Collective on IS

836:    Input Parmeters:
837: .  is - the index set

839:    Output Parameters:
840: .  isnew - the copy of the index set

842:    Notes:
843:    ISDuplicate() does not copy the index set, but rather allocates storage
844:    for the new one.  Use ISCopy() to copy an index set.

846:    Level: beginner

848:    Concepts: index sets^duplicating

850: .seealso: ISCreateGeneral(), ISCopy()
851: @*/
852: PetscErrorCode  ISDuplicate(IS is,IS *newIS)
853: {

859:   (*is->ops->duplicate)(is,newIS);
860:   return(0);
861: }

865: /*@
866:    ISCopy - Copies an index set.

868:    Collective on IS

870:    Input Parmeters:
871: .  is - the index set

873:    Output Parameters:
874: .  isy - the copy of the index set

876:    Level: beginner

878:    Concepts: index sets^copying

880: .seealso: ISDuplicate()
881: @*/
882: PetscErrorCode  ISCopy(IS is,IS isy)
883: {

890:   if (is == isy) return(0);
891:   (*is->ops->copy)(is,isy);
892:   isy->isperm     = is->isperm;
893:   isy->max        = is->max;
894:   isy->min        = is->min;
895:   isy->isidentity = is->isidentity;
896:   return(0);
897: }

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

904:    Collective on IS and comm

906:    Input Arguments:
907: + is - index set
908: . comm - communicator for new index set
909: - mode - copy semantics, PETSC_USE_POINTER for no-copy if possible, otherwise PETSC_COPY_VALUES

911:    Output Arguments:
912: . newis - new IS on comm

914:    Level: advanced

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

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

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

924: .seealso: ISSplit()
925: @*/
926: PetscErrorCode  ISOnComm(IS is,MPI_Comm comm,PetscCopyMode mode,IS *newis)
927: {
929:   PetscMPIInt match;

934:   MPI_Comm_compare(((PetscObject)is)->comm,comm,&match);
935:   if (mode != PETSC_COPY_VALUES && (match == MPI_IDENT || match == MPI_CONGRUENT)) {
936:     PetscObjectReference((PetscObject)is);
937:     *newis = is;
938:   } else {
939:     (*is->ops->oncomm)(is,comm,mode,newis);
940:   }
941:   return(0);
942: }

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

949:    Logicall Collective on IS

951:    Input Arguments:
952: + is - index set
953: - bs - block size

955:    Level: intermediate

957: .seealso: ISGetBlockSize(), ISCreateBlock()
958: @*/
959: PetscErrorCode  ISSetBlockSize(IS is,PetscInt bs)
960: {

966:   if (bs < 1) SETERRQ1(((PetscObject)is)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Block size %D, must be positive",bs);
967:   (*is->ops->setblocksize)(is,bs);
968:   return(0);
969: }

973: /*@
974:    ISGetBlockSize - Returns the number of elements in a block.

976:    Not Collective

978:    Input Parameter:
979: .  is - the index set

981:    Output Parameter:
982: .  size - the number of elements in a block

984:    Level: intermediate

986:    Concepts: IS^block size
987:    Concepts: index sets^block size

989: .seealso: ISBlockGetSize(), ISGetSize(), ISCreateBlock(), ISSetBlockSize()
990: @*/
991: PetscErrorCode  ISGetBlockSize(IS is,PetscInt *size)
992: {
994:   *size = is->bs;
995:   return(0);
996: }

998: EXTERN_C_BEGIN
1001: PetscErrorCode ISGetIndicesCopy(IS is, PetscInt idx[])
1002: {
1004:   PetscInt       len,i;
1005:   const PetscInt *ptr;

1008:   ISGetSize(is,&len);
1009:   ISGetIndices(is,&ptr);
1010:   for(i=0;i<len;i++) idx[i] = ptr[i];
1011:   ISRestoreIndices(is,&ptr);
1012:   return(0);
1013: }
1014: EXTERN_C_END

1016: /*MC
1017:     ISGetIndicesF90 - Accesses the elements of an index set from Fortran90.
1018:     The users should call ISRestoreIndicesF90() after having looked at the
1019:     indices.  The user should NOT change the indices.

1021:     Synopsis:
1022:     ISGetIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)

1024:     Not collective

1026:     Input Parameter:
1027: .   x - index set

1029:     Output Parameters:
1030: +   xx_v - the Fortran90 pointer to the array
1031: -   ierr - error code

1033:     Example of Usage: 
1034: .vb
1035:     PetscScalar, pointer xx_v(:)
1036:     ....
1037:     call ISGetIndicesF90(x,xx_v,ierr)
1038:     a = xx_v(3)
1039:     call ISRestoreIndicesF90(x,xx_v,ierr)
1040: .ve

1042:     Notes:
1043:     Not yet supported for all F90 compilers.

1045:     Level: intermediate

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

1049:   Concepts: index sets^getting indices in f90
1050:   Concepts: indices of index set in f90

1052: M*/

1054: /*MC
1055:     ISRestoreIndicesF90 - Restores an index set to a usable state after
1056:     a call to ISGetIndicesF90().

1058:     Synopsis:
1059:     ISRestoreIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)

1061:     Not collective

1063:     Input Parameters:
1064: .   x - index set
1065: .   xx_v - the Fortran90 pointer to the array

1067:     Output Parameter:
1068: .   ierr - error code


1071:     Example of Usage: 
1072: .vb
1073:     PetscScalar, pointer xx_v(:)
1074:     ....
1075:     call ISGetIndicesF90(x,xx_v,ierr)
1076:     a = xx_v(3)
1077:     call ISRestoreIndicesF90(x,xx_v,ierr)
1078: .ve
1079:    
1080:     Notes:
1081:     Not yet supported for all F90 compilers.

1083:     Level: intermediate

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

1087: M*/

1089: /*MC
1090:     ISBlockGetIndicesF90 - Accesses the elements of an index set from Fortran90.
1091:     The users should call ISBlockRestoreIndicesF90() after having looked at the
1092:     indices.  The user should NOT change the indices.

1094:     Synopsis:
1095:     ISBlockGetIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)

1097:     Not collective

1099:     Input Parameter:
1100: .   x - index set

1102:     Output Parameters:
1103: +   xx_v - the Fortran90 pointer to the array
1104: -   ierr - error code
1105:     Example of Usage: 
1106: .vb
1107:     PetscScalar, pointer xx_v(:)
1108:     ....
1109:     call ISBlockGetIndicesF90(x,xx_v,ierr)
1110:     a = xx_v(3)
1111:     call ISBlockRestoreIndicesF90(x,xx_v,ierr)
1112: .ve

1114:     Notes:
1115:     Not yet supported for all F90 compilers

1117:     Level: intermediate

1119: .seealso:  ISBlockRestoreIndicesF90(), ISGetIndices(), ISRestoreIndices(),
1120:            ISRestoreIndices()

1122:   Concepts: index sets^getting block indices in f90
1123:   Concepts: indices of index set in f90
1124:   Concepts: block^ indices of index set in f90

1126: M*/

1128: /*MC
1129:     ISBlockRestoreIndicesF90 - Restores an index set to a usable state after
1130:     a call to ISBlockGetIndicesF90().

1132:     Synopsis:
1133:     ISBlockRestoreIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)

1135:     Not Collective

1137:     Input Parameters:
1138: +   x - index set
1139: -   xx_v - the Fortran90 pointer to the array

1141:     Output Parameter:
1142: .   ierr - error code

1144:     Example of Usage: 
1145: .vb
1146:     PetscScalar, pointer xx_v(:)
1147:     ....
1148:     call ISBlockGetIndicesF90(x,xx_v,ierr)
1149:     a = xx_v(3)
1150:     call ISBlockRestoreIndicesF90(x,xx_v,ierr)
1151: .ve
1152:    
1153:     Notes:
1154:     Not yet supported for all F90 compilers

1156:     Level: intermediate

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

1160: M*/