Actual source code: index.c

  1: #define PETSCVEC_DLL

  3: /*  
  4:    Defines the abstract operations on index sets, i.e. the public interface. 
  5: */
 6:  #include src/vec/is/isimpl.h

  8: /* Logging support */
  9: PetscCookie  IS_COOKIE = 0;

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

 16:    Collective on IS

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

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

 24:    Level: intermediate

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

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

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

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

 51:    Collective on IS

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

 56:    Level: intermediate

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

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

 73: /*@
 74:    ISPermutation - PETSC_TRUE or PETSC_FALSE depending on whether the 
 75:    index set has been declared to be a permutation.

 77:    Collective on IS

 79:    Input Parmeters:
 80: .  is - the index set

 82:    Output Parameters:
 83: .  perm - PETSC_TRUE if a permutation, else PETSC_FALSE

 85:    Level: intermediate

 87:   Concepts: permutation
 88:   Concepts: index sets^is permutation

 90: .seealso: ISSetPermutation()
 91: @*/
 92: PetscErrorCode  ISPermutation(IS is,PetscTruth *perm)
 93: {
 97:   *perm = (PetscTruth) is->isperm;
 98:   return(0);
 99: }

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

106:    Collective on IS

108:    Input Parmeters:
109: .  is - the index set

111:    Level: intermediate

113:   Concepts: permutation
114:   Concepts: index sets^permutation

116: .seealso: ISPermutation()
117: @*/
118: PetscErrorCode  ISSetPermutation(IS is)
119: {
122:   is->isperm = PETSC_TRUE;
123:   return(0);
124: }

128: /*@
129:    ISDestroy - Destroys an index set.

131:    Collective on IS

133:    Input Parameters:
134: .  is - the index set

136:    Level: beginner

138: .seealso: ISCreateGeneral(), ISCreateStride(), ISCreateBlocked()
139: @*/
140: PetscErrorCode  ISDestroy(IS is)
141: {

146:   if (--is->refct > 0) return(0);
147:   (*is->ops->destroy)(is);
148:   return(0);
149: }

153: /*@
154:    ISInvertPermutation - Creates a new permutation that is the inverse of 
155:                          a given permutation.

157:    Collective on IS

159:    Input Parameter:
160: +  is - the index set
161: -  nlocal - number of indices on this processor in result (ignored for 1 proccessor) or
162:             use PETSC_DECIDE

164:    Output Parameter:
165: .  isout - the inverse permutation

167:    Level: intermediate

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

172:    Concepts: inverse permutation
173:    Concepts: permutation^inverse
174:    Concepts: index sets^inverting
175: @*/
176: PetscErrorCode  ISInvertPermutation(IS is,PetscInt nlocal,IS *isout)
177: {

183:   if (!is->isperm) SETERRQ(PETSC_ERR_ARG_WRONG,"not a permutation");
184:   (*is->ops->invertpermutation)(is,nlocal,isout);
185:   ISSetPermutation(*isout);
186:   return(0);
187: }

191: /*@
192:    ISGetSize - Returns the global length of an index set. 

194:    Not Collective

196:    Input Parameter:
197: .  is - the index set

199:    Output Parameter:
200: .  size - the global size

202:    Level: beginner

204:    Concepts: size^of index set
205:    Concepts: index sets^size

207: @*/
208: PetscErrorCode  ISGetSize(IS is,PetscInt *size)
209: {

215:   (*is->ops->getsize)(is,size);
216:   return(0);
217: }

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

224:    Not Collective

226:    Input Parameter:
227: .  is - the index set

229:    Output Parameter:
230: .  size - the local size

232:    Level: beginner

234:    Concepts: size^of index set
235:    Concepts: local size^of index set
236:    Concepts: index sets^local size
237:   
238: @*/
239: PetscErrorCode  ISGetLocalSize(IS is,PetscInt *size)
240: {

246:   (*is->ops->getlocalsize)(is,size);
247:   return(0);
248: }

252: /*@C
253:    ISGetIndices - Returns a pointer to the indices.  The user should call 
254:    ISRestoreIndices() after having looked at the indices.  The user should 
255:    NOT change the indices.

257:    Not Collective

259:    Input Parameter:
260: .  is - the index set

262:    Output Parameter:
263: .  ptr - the location to put the pointer to the indices

265:    Fortran Note:
266:    This routine is used differently from Fortran
267: $    IS          is
268: $    integer     is_array(1)
269: $    PetscOffset i_is
270: $    int         ierr
271: $       call ISGetIndices(is,is_array,i_is,ierr)
272: $
273: $   Access first local entry in list
274: $      value = is_array(i_is + 1)
275: $
276: $      ...... other code
277: $       call ISRestoreIndices(is,is_array,i_is,ierr)

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

282:    Level: intermediate

284:    Concepts: index sets^getting indices
285:    Concepts: indices of index set

287: .seealso: ISRestoreIndices(), ISGetIndicesF90()
288: @*/
289: PetscErrorCode  ISGetIndices(IS is,PetscInt *ptr[])
290: {

296:   (*is->ops->getindices)(is,ptr);
297:   return(0);
298: }

302: /*@C
303:    ISRestoreIndices - Restores an index set to a usable state after a call 
304:                       to ISGetIndices().

306:    Not Collective

308:    Input Parameters:
309: +  is - the index set
310: -  ptr - the pointer obtained by ISGetIndices()

312:    Fortran Note:
313:    This routine is used differently from Fortran
314: $    IS          is
315: $    integer     is_array(1)
316: $    PetscOffset i_is
317: $    int         ierr
318: $       call ISGetIndices(is,is_array,i_is,ierr)
319: $
320: $   Access first local entry in list
321: $      value = is_array(i_is + 1)
322: $
323: $      ...... other code
324: $       call ISRestoreIndices(is,is_array,i_is,ierr)

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

329:    Level: intermediate

331: .seealso: ISGetIndices(), ISRestoreIndicesF90()
332: @*/
333: PetscErrorCode  ISRestoreIndices(IS is,PetscInt *ptr[])
334: {

340:   if (is->ops->restoreindices) {
341:     (*is->ops->restoreindices)(is,ptr);
342:   }
343:   return(0);
344: }

348: /*@C
349:    ISView - Displays an index set.

351:    Collective on IS

353:    Input Parameters:
354: +  is - the index set
355: -  viewer - viewer used to display the set, for example PETSC_VIEWER_STDOUT_SELF.

357:    Level: intermediate

359: .seealso: PetscViewerASCIIOpen()
360: @*/
361: PetscErrorCode  ISView(IS is,PetscViewer viewer)
362: {

367:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(is->comm);
370: 
371:   (*is->ops->view)(is,viewer);
372:   return(0);
373: }

377: /*@
378:    ISSort - Sorts the indices of an index set.

380:    Collective on IS

382:    Input Parameters:
383: .  is - the index set

385:    Level: intermediate

387:    Concepts: index sets^sorting
388:    Concepts: sorting^index set

390: .seealso: ISSorted()
391: @*/
392: PetscErrorCode  ISSort(IS is)
393: {

398:   (*is->ops->sortindices)(is);
399:   return(0);
400: }

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

407:    Collective on IS

409:    Input Parameter:
410: .  is - the index set

412:    Output Parameter:
413: .  flg - output flag, either PETSC_TRUE if the index set is sorted, 
414:          or PETSC_FALSE otherwise.

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

420:    Level: intermediate

422: .seealso: ISSort()
423: @*/
424: PetscErrorCode  ISSorted(IS is,PetscTruth *flg)
425: {

431:   (*is->ops->sorted)(is,flg);
432:   return(0);
433: }

437: /*@
438:    ISDuplicate - Creates a duplicate copy of an index set.

440:    Collective on IS

442:    Input Parmeters:
443: .  is - the index set

445:    Output Parameters:
446: .  isnew - the copy of the index set

448:    Level: beginner

450:    Concepts: index sets^duplicating

452: .seealso: ISCreateGeneral()
453: @*/
454: PetscErrorCode  ISDuplicate(IS is,IS *newIS)
455: {

461:   (*is->ops->duplicate)(is,newIS);
462:   return(0);
463: }


466: /*MC
467:     ISGetIndicesF90 - Accesses the elements of an index set from Fortran90.
468:     The users should call ISRestoreIndicesF90() after having looked at the
469:     indices.  The user should NOT change the indices.

471:     Synopsis:
472:     ISGetIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)

474:     Not collective

476:     Input Parameter:
477: .   x - index set

479:     Output Parameters:
480: +   xx_v - the Fortran90 pointer to the array
481: -   ierr - error code

483:     Example of Usage: 
484: .vb
485:     PetscScalar, pointer xx_v(:)
486:     ....
487:     call ISGetIndicesF90(x,xx_v,ierr)
488:     a = xx_v(3)
489:     call ISRestoreIndicesF90(x,xx_v,ierr)
490: .ve

492:     Notes:
493:     Not yet supported for all F90 compilers.

495:     Level: intermediate

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

499:   Concepts: index sets^getting indices in f90
500:   Concepts: indices of index set in f90

502: M*/

504: /*MC
505:     ISRestoreIndicesF90 - Restores an index set to a usable state after
506:     a call to ISGetIndicesF90().

508:     Synopsis:
509:     ISRestoreIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)

511:     Not collective

513:     Input Parameters:
514: .   x - index set
515: .   xx_v - the Fortran90 pointer to the array

517:     Output Parameter:
518: .   ierr - error code


521:     Example of Usage: 
522: .vb
523:     PetscScalar, pointer xx_v(:)
524:     ....
525:     call ISGetIndicesF90(x,xx_v,ierr)
526:     a = xx_v(3)
527:     call ISRestoreIndicesF90(x,xx_v,ierr)
528: .ve
529:    
530:     Notes:
531:     Not yet supported for all F90 compilers.

533:     Level: intermediate

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

537: M*/

539: /*MC
540:     ISBlockGetIndicesF90 - Accesses the elements of an index set from Fortran90.
541:     The users should call ISBlockRestoreIndicesF90() after having looked at the
542:     indices.  The user should NOT change the indices.

544:     Synopsis:
545:     ISBlockGetIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)

547:     Not collective

549:     Input Parameter:
550: .   x - index set

552:     Output Parameters:
553: +   xx_v - the Fortran90 pointer to the array
554: -   ierr - error code
555:     Example of Usage: 
556: .vb
557:     PetscScalar, pointer xx_v(:)
558:     ....
559:     call ISBlockGetIndicesF90(x,xx_v,ierr)
560:     a = xx_v(3)
561:     call ISBlockRestoreIndicesF90(x,xx_v,ierr)
562: .ve

564:     Notes:
565:     Not yet supported for all F90 compilers

567:     Level: intermediate

569: .seealso:  ISBlockRestoreIndicesF90(), ISGetIndices(), ISRestoreIndices(),
570:            ISRestoreIndices()

572:   Concepts: index sets^getting block indices in f90
573:   Concepts: indices of index set in f90
574:   Concepts: block^ indices of index set in f90

576: M*/

578: /*MC
579:     ISBlockRestoreIndicesF90 - Restores an index set to a usable state after
580:     a call to ISBlockGetIndicesF90().

582:     Synopsis:
583:     ISBlockRestoreIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)

585:     Input Parameters:
586: +   x - index set
587: -   xx_v - the Fortran90 pointer to the array

589:     Output Parameter:
590: .   ierr - error code

592:     Example of Usage: 
593: .vb
594:     PetscScalar, pointer xx_v(:)
595:     ....
596:     call ISBlockGetIndicesF90(x,xx_v,ierr)
597:     a = xx_v(3)
598:     call ISBlockRestoreIndicesF90(x,xx_v,ierr)
599: .ve
600:    
601:     Notes:
602:     Not yet supported for all F90 compilers

604:     Level: intermediate

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

608: M*/