Actual source code: ex1f.F

  1: !
  2: !  "$Id: ex1f.F,v 1.24 2001/01/15 21:44:33 bsmith Exp $";
  3: !
  4: !  Description: Creates an index set based on a set of integers. Views that index set
  5: !  and then destroys it.
  6: !
  7: !/*T
  8: !    Concepts: index sets^manipulating a general index set;
  9: !T*/
 10: !
 11: !
 12: !  The following include statements are required for Fortran programs
 13: !  that use PETSc index sets:
 14: !     petsc.h  - base PETSc routines
 15: !     petscis.h     - index sets (IS objects)
 16: !
 17:       program main
 18:       implicit none

 20:  #include include/finclude/petsc.h
 21:  #include include/finclude/petscis.h

 23:       integer     ierr,indices(5),rank,n,index1,index5
 24:       PetscOffset ix
 25:       IS          is

 27:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 28:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)

 30: !  Create an index set with 5 entries. Each processor creates
 31: !  its own index set with its own list of integers.
 32: 
 33:       indices(1) = rank + 1
 34:       indices(2) = rank + 2
 35:       indices(3) = rank + 3
 36:       indices(4) = rank + 4
 37:       indices(5) = rank + 5
 38:       call ISCreateGeneral(PETSC_COMM_SELF,5,indices,is,ierr)

 40: !  Print the index set to stdout

 42:       call ISView(is,PETSC_VIEWER_STDOUT_SELF,ierr)

 44: !  Get the number of indices in the set

 46:       call ISGetLocalSize(is,n,ierr)

 48: !   Get the indices in the index set

 50:       call ISGetIndices(is,indices,ix,ierr)

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

 55: !
 56: !      Bug in IRIX64-F90 libraries - write/format cannot handle integer(integer*8 + integer)
 57: !

 59:       index1 = indices(ix+1)
 60:       index5 = indices(ix+5)
 61:       write(6,100) rank,index1,index5
 62:  100  format('[',i5,'] First index = ',i5,' fifth index = ',i5)
 63: 
 64: !   Once we no longer need access to the indices they should
 65: !   returned to the system

 67:       call ISRestoreIndices(is,indices,ix,ierr)
 68: 
 69: !   All PETSc objects should be destroyed once they are
 70: !   no longer needed

 72:       call ISDestroy(is,ierr)
 73:       call PetscFinalize(ierr)
 74:       end

 76: