Actual source code: ex1f.F90

  1: !
  2: !  Description: Creates an index set based on a set of integers. Views that index set
  3: !  and then destroys it.
  4: !
  5: !
  6: !
  7:       program main
  8: #include <petsc/finclude/petscis.h>
  9:       use petscis
 10:       implicit none

 12:       PetscErrorCode ierr
 13:       PetscInt n,indices(5),index1,index5
 14:       PetscMPIInt rank
 15:       IS          is
 16:       PetscInt, pointer :: indices2(:)

 18:       PetscCallA(PetscInitialize(ierr))
 19:       PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr))

 21: !  Create an index set with 5 entries. Each processor creates
 22: !  its own index set with its own list of integers.

 24:       indices(1) = rank + 1
 25:       indices(2) = rank + 2
 26:       indices(3) = rank + 3
 27:       indices(4) = rank + 4
 28:       indices(5) = rank + 5

 30: !     if using 64-bit integers cannot pass 5 into routine expecting an integer*8
 31:       n = 5
 32:       PetscCallA(ISCreateGeneral(PETSC_COMM_SELF,n,indices,PETSC_COPY_VALUES,is,ierr))

 34: !  Print the index set to stdout

 36:       PetscCallA(ISView(is,PETSC_VIEWER_STDOUT_SELF,ierr))

 38: !  Get the number of indices in the set

 40:       PetscCallA(ISGetLocalSize(is,n,ierr))

 42: !   Get the indices in the index set

 44:       PetscCallA(ISGetIndicesF90(is,indices2,ierr))

 46: !   Now any code that needs access to the list of integers
 47: !   has access to it here

 49: !
 50: !      Bug in IRIX64-F90 libraries - write/format cannot handle integer(integer*8 + integer)
 51: !

 53:       index1 = indices(1)
 54:       index5 = indices(5)
 55:       write(6,100) rank,index1,index5
 56:  100  format('[',i5,'] First index = ',i5,' fifth index = ',i5)

 58: !   Once we no longer need access to the indices they should
 59: !   returned to the system

 61:       PetscCallA(ISRestoreIndicesF90(is,indices2,ierr))

 63: !   All PETSc objects should be destroyed once they are
 64: !   no longer needed

 66:       PetscCallA(ISDestroy(is,ierr))
 67:       PetscCallA(PetscFinalize(ierr))
 68:       end

 70: !/*TEST
 71: !
 72: !   test:
 73: !      filter: sort -b
 74: !      filter_output: sort -b
 75: !
 76: !TEST*/