Actual source code: ex1f90.F

  1: !
  2: !  Description: Creates an index set based on a set of integers. Views that index set
  3: !  and then destroys it.
  4: !
  5: !/*T
  6: !    Concepts: index sets^manipulating a general index set;
  7: !    Concepts: Fortran90^accessing indices of index set;
  8: !T*/
  9: !
 10: !  The following include statements are required for Fortran programs
 11: !  that use PETSc index sets:
 12: !     petsc.h  - base PETSc routines
 13: !     petscis.h     - index sets (IS objects)
 14: !     petscis.h90   - to allow access to Fortran90 features of index sets
 15: !
 16:       program main
 17:       implicit none

 19:  #include finclude/petsc.h
 20:  #include finclude/petscis.h
 21: #include "finclude/petscis.h90"

 23:       integer ierr,indices(5),rank,n
 24:       integer, pointer :: idx(:)
 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 ISGetIndicesF90(is,idx,ierr)

 52:       if (associated(idx)) then
 53:          write (*,*) 'Association check passed'
 54:       else
 55:          write (*,*) 'Association check failed'
 56:       endif

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

 61:       write(6,50) idx
 62:  50   format(5I3)

 64:       write(6,100) rank,idx(1),idx(5)
 65:  100  format('[',i5,'] First index = ',i5,' fifth index = ',i5)
 66: 
 67: !   Once we no longer need access to the indices they should
 68: !   returned to the system

 70:       call ISRestoreIndicesF90(is,idx,ierr)
 71: 
 72: !   All PETSc objects should be destroyed once they are
 73: !   no longer needed

 75:       call ISDestroy(is,ierr)
 76:       call PetscFinalize(ierr)
 77:       end

 79: